OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Class representing a general abstraction for surface finite element interpolation class. More...
#include <feinterpol2d.h>
Public Member Functions | |
FEInterpolation2d (int o, int ind1, int ind2) | |
virtual int | giveNsd () |
Returns number of spatial dimensions. More... | |
virtual double | giveArea (const FEICellGeometry &cellgeo) const |
Computes the exact area. More... | |
virtual double | giveCharacteristicLength (const FEICellGeometry &cellgeo) const |
Returns a characteristic length of the geometry, typically a diagonal or edge length. More... | |
virtual int | global2local (FloatArray &answer, const FloatArray &gcoords, const FEICellGeometry &cellgeo) |
Default implementation using Newton's method to find the local coordinates. More... | |
virtual void | giveJacobianMatrixAt (FloatMatrix &jacobianMatrix, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Gives the jacobian matrix at the local coordinates. More... | |
virtual bool | inside (const FloatArray &lcoords) const |
Boundary interpolation services. | |
Boundary is defined as entity of one dimension lower than the interpolation represents | |
virtual void | boundaryEdgeGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
virtual void | boundaryEdgeEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the basis functions on the requested boundary. More... | |
virtual double | boundaryEdgeGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation Jacobian on the requested boundary. More... | |
virtual void | boundaryEdgeLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Maps the local boundary coordinates to global. More... | |
virtual void | boundaryGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
virtual void | boundaryEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the basis functions on the requested boundary. More... | |
virtual double | boundaryEvalNormal (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the normal on the requested boundary. More... | |
virtual double | boundaryGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation Jacobian on the requested boundary. More... | |
virtual void | boundaryLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Maps the local boundary coordinates to global. More... | |
Surface interpolation services | |
virtual void | boundarySurfaceEvalN (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the array of edge interpolation functions (shape functions) at given point. More... | |
virtual void | boundarySurfaceEvaldNdx (FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point. More... | |
virtual double | boundarySurfaceEvalNormal (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the normal out of the surface at given point. More... | |
virtual void | boundarySurfaceLocal2global (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates edge global coordinates from given local ones. More... | |
virtual double | boundarySurfaceGiveTransformationJacobian (int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the edge jacobian of transformation between local and global coordinates. More... | |
virtual void | boundarySurfaceGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
Edge interpolation services. | |
virtual void | computeLocalEdgeMapping (IntArray &edgeNodes, int iedge)=0 |
void | computeEdgeMapping (IntArray &edgeNodes, IntArray &elemNodes, int iedge) |
virtual void | edgeEvalN (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates the array of edge interpolation functions (shape functions) at given point. More... | |
virtual double | edgeEvalNormal (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates the normal on the given edge. More... | |
virtual void | edgeEvaldNds (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point. More... | |
virtual void | edgeLocal2global (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates edge global coordinates from given local ones. More... | |
virtual double | edgeGiveTransformationJacobian (int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the edge Jacobian of transformation between local and global coordinates. More... | |
Public Member Functions inherited from oofem::FEInterpolation | |
FEInterpolation (int o) | |
virtual | ~FEInterpolation () |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual int | giveNumberOfNodes () const |
Returns the number of geometric nodes of the receiver. More... | |
std::string | errorInfo (const char *func) const |
virtual integrationDomain | giveIntegrationDomain () const =0 |
Returns the integration domain of the interpolator. More... | |
virtual Element_Geometry_Type | giveGeometryType () const =0 |
Returns the geometry type fo the interpolator. More... | |
int | giveInterpolationOrder () |
Returns the interpolation order. More... | |
virtual void | evalN (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates the array of interpolation functions (shape functions) at given point. More... | |
virtual double | evaldNdx (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | evald2Ndx2 (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of second derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | evaldNdxi (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | giveLocalNodeCoords (FloatMatrix &answer) |
Returns a matrix containing the local coordinates for each node corresponding to the interpolation. More... | |
virtual void | local2global (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0 |
Evaluates global coordinates from given local ones. More... | |
virtual double | giveTransformationJacobian (const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation. More... | |
virtual IntegrationRule * | giveIntegrationRule (int order) |
Sets up a suitable integration rule for numerical integrating over volume. More... | |
virtual integrationDomain | giveBoundaryEdgeIntegrationDomain (int boundary) const =0 |
Returns boundary integration domain. More... | |
virtual IntegrationRule * | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
virtual integrationDomain | giveBoundarySurfaceIntegrationDomain (int boundary) const =0 |
Returns boundary integration domain. More... | |
virtual IntegrationRule * | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
virtual double | evalNXIntegral (int boundary, const FEICellGeometry &cellgeo) |
Computes the integral . More... | |
virtual integrationDomain | giveBoundaryIntegrationDomain (int boundary) const =0 |
Returns boundary integration domain. More... | |
virtual IntegrationRule * | giveBoundaryIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
virtual int | giveKnotSpanBasisFuncMask (const IntArray &knotSpan, IntArray &mask) |
Returns indices (zero based) of nonzero basis functions for given knot span. More... | |
virtual int | giveNumberOfKnotSpanBasisFunctions (const IntArray &knotSpan) |
Returns the number of nonzero basis functions at individual knot span,. More... | |
virtual bool | hasSubPatchFormulation () |
Returns true, if receiver is formulated on sub-patch basis. More... | |
virtual const double *const * | giveKnotVector () |
Returns the subdivision of patch parametric space. More... | |
virtual int | giveNumberOfKnotSpans (int dim) |
Returns the number of knot spans of the receiver. More... | |
virtual const FloatArray * | giveKnotValues (int dim) |
Returns the knot values of the receiver. More... | |
virtual const IntArray * | giveKnotMultiplicity (int dim) |
Returns the knot multiplicity of the receiver. More... | |
virtual int | giveNumberOfEdges () const |
Returns number of edges. More... | |
Protected Attributes | |
int | xind |
int | yind |
Protected Attributes inherited from oofem::FEInterpolation | |
int | order |
Class representing a general abstraction for surface finite element interpolation class.
Definition at line 45 of file feinterpol2d.h.
|
inline |
Definition at line 51 of file feinterpol2d.h.
|
virtual |
Evaluates the basis functions on the requested boundary.
Only basis functions that are nonzero anywhere on the boundary are given. Ordering can be obtained from giveBoundaryNodes. Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
answer | Basis functions Array to be filled with the boundary nodes. |
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 45 of file feinterpol2d.C.
References edgeEvalN().
|
virtual |
Gives the boundary nodes for requested boundary number.
answer | Array to be filled with the boundary nodes. |
boundary | Boundary number. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::FEI2dLineLin.
Definition at line 40 of file feinterpol2d.C.
References computeLocalEdgeMapping().
|
virtual |
Evaluates the determinant of the transformation Jacobian on the requested boundary.
Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::FEI2dQuadQuadAxi, oofem::FEI2dQuadLinAxi, and oofem::FEI2dTrLinAxi.
Definition at line 50 of file feinterpol2d.C.
References edgeGiveTransformationJacobian().
|
virtual |
Maps the local boundary coordinates to global.
Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
answer | Global coordinates. |
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 55 of file feinterpol2d.C.
References edgeLocal2global().
|
virtual |
Evaluates the basis functions on the requested boundary.
Only basis functions that are nonzero anywhere on the boundary are given. Ordering can be obtained from giveBoundaryNodes. Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
answer | Basis functions Array to be filled with the boundary nodes. |
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 144 of file feinterpol2d.C.
References edgeEvalN().
|
virtual |
Evaluates the normal on the requested boundary.
answer | The evaluated normal. |
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 149 of file feinterpol2d.C.
References edgeEvalNormal().
|
virtual |
Gives the boundary nodes for requested boundary number.
Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
answer | Array to be filled with the boundary nodes. |
boundary | Boundary number. |
Implements oofem::FEInterpolation.
Definition at line 139 of file feinterpol2d.C.
References computeLocalEdgeMapping().
Referenced by oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement().
|
virtual |
Evaluates the determinant of the transformation Jacobian on the requested boundary.
Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::FEI2dQuadQuadAxi, oofem::FEI2dQuadLinAxi, and oofem::FEI2dTrLinAxi.
Definition at line 154 of file feinterpol2d.C.
References edgeGiveTransformationJacobian().
Referenced by oofem::Tr21Stokes::computeBoundarySurfaceLoadVector(), and oofem::Tr1BubbleStokes::computeBoundarySurfaceLoadVector().
|
virtual |
Maps the local boundary coordinates to global.
Boundaries are defined as the corner nodes for 1D geometries, edges for 2D geometries and surfaces for 3D geometries.
answer | Global coordinates. |
boundary | Boundary number. |
lcoords | The local coordinates (on the boundary local coordinate system). |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 159 of file feinterpol2d.C.
References edgeLocal2global().
Referenced by oofem::Tr21Stokes::computeBoundarySurfaceLoadVector(), and oofem::Tr1BubbleStokes::computeBoundarySurfaceLoadVector().
|
virtual |
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point.
These derivatives are in global coordinate system (where the nodal coordinates are defined).
answer | Contains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi. |
isurf | Determines the surface number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 186 of file feinterpol2d.C.
References oofem::FEInterpolation::evaldNdx().
|
virtual |
Evaluates the array of edge interpolation functions (shape functions) at given point.
answer | Contains resulting array of evaluated interpolation functions. |
isurf | Surface number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 181 of file feinterpol2d.C.
References oofem::FEInterpolation::evalN().
|
virtual |
Evaluates the normal out of the surface at given point.
answer | Contains resulting normal vector. |
isurf | Determines the surface number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 192 of file feinterpol2d.C.
References oofem::FEInterpolation::giveTransformationJacobian().
|
virtual |
Gives the boundary nodes for requested boundary number.
answer | Array to be filled with the boundary nodes. |
boundary | Boundary number. |
Implements oofem::FEInterpolation.
Definition at line 211 of file feinterpol2d.C.
References oofem::IntArray::at(), oofem::FEInterpolation::giveNumberOfNodes(), and oofem::IntArray::resize().
|
virtual |
Evaluates the edge jacobian of transformation between local and global coordinates.
isurf | Determines the surface number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 205 of file feinterpol2d.C.
References oofem::FEInterpolation::giveTransformationJacobian().
|
virtual |
Evaluates edge global coordinates from given local ones.
These derivatives are in global coordinate system (where the nodal coordinates are defined).
answer | Contains resulting global coordinates. |
isurf | Determines the surface number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 199 of file feinterpol2d.C.
References oofem::FEInterpolation::local2global().
void oofem::FEInterpolation2d::computeEdgeMapping | ( | IntArray & | edgeNodes, |
IntArray & | elemNodes, | ||
int | iedge | ||
) |
Definition at line 164 of file feinterpol2d.C.
References oofem::IntArray::at(), computeLocalEdgeMapping(), oofem::IntArray::giveSize(), and oofem::IntArray::resize().
|
pure virtual |
Implemented in oofem::FEI2dQuadQuad, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dLineHermite, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, oofem::FEI2dQuadConst, and oofem::FEI2dTrConst.
Referenced by boundaryEdgeGiveNodes(), boundaryGiveNodes(), computeEdgeMapping(), oofem::PrescribedGradientBCNeumann::integrateTangent(), and oofem::XfemElementInterface::partitionEdgeSegment().
|
pure virtual |
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point.
These derivatives are in global coordinate system (where the nodal coordinates are defined).
answer | Contains resulting array of derivatives, the member at i position contains value of . |
iedge | Determines the edge number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implemented in oofem::FEI2dQuadQuad, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dLineHermite, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, oofem::FEI2dQuadConst, and oofem::FEI2dTrConst.
Referenced by oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement(), and oofem::SurfaceTensionBoundaryCondition::computeTangentFromElement().
|
pure virtual |
Evaluates the array of edge interpolation functions (shape functions) at given point.
answer | Contains resulting array of evaluated interpolation functions. |
iedge | Edge number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implemented in oofem::FEI2dQuadQuad, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dLineHermite, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, oofem::FEI2dQuadConst, and oofem::FEI2dTrConst.
Referenced by boundaryEdgeEvalN(), and boundaryEvalN().
|
pure virtual |
Evaluates the normal on the given edge.
answer | Contains the evaluated normal. |
iedge | Determines the edge number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implemented in oofem::FEI2dQuadQuad, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dLineHermite, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, oofem::FEI2dQuadConst, and oofem::FEI2dTrConst.
Referenced by boundaryEvalNormal(), and edgeGiveTransformationJacobian().
|
virtual |
Evaluates the edge Jacobian of transformation between local and global coordinates.
iedge | Determines edge number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Reimplemented in oofem::FEI2dQuadQuadAxi, oofem::FEI2dQuadLinAxi, and oofem::FEI2dTrLinAxi.
Definition at line 175 of file feinterpol2d.C.
References edgeEvalNormal().
Referenced by boundaryEdgeGiveTransformationJacobian(), oofem::FEI2dTrLinAxi::boundaryEdgeGiveTransformationJacobian(), oofem::FEI2dQuadLinAxi::boundaryEdgeGiveTransformationJacobian(), oofem::FEI2dQuadQuadAxi::boundaryEdgeGiveTransformationJacobian(), boundaryGiveTransformationJacobian(), oofem::FEI2dTrLinAxi::boundaryGiveTransformationJacobian(), oofem::FEI2dQuadLinAxi::boundaryGiveTransformationJacobian(), oofem::FEI2dQuadQuadAxi::boundaryGiveTransformationJacobian(), oofem::TR1_2D_PFEM::computeEdgeBCSubVectorAt(), oofem::TrAxisym1_ht::computeEdgeVolumeAround(), oofem::QuadAxisym1_ht::computeEdgeVolumeAround(), oofem::Tr1Darcy::computeEdgeVolumeAround(), oofem::Quad1_ht::computeEdgeVolumeAround(), oofem::Tr1_ht::computeEdgeVolumeAround(), oofem::DKTPlate3d::computeEdgeVolumeAround(), oofem::CCTPlate::computeEdgeVolumeAround(), oofem::Tr_Warp::computeEdgeVolumeAround(), oofem::Quad1Mindlin::computeEdgeVolumeAround(), oofem::QDKTPlate::computeEdgeVolumeAround(), oofem::DKTPlate::computeEdgeVolumeAround(), oofem::Quad1MindlinShell3D::computeEdgeVolumeAround(), oofem::MITC4Shell::computeEdgeVolumeAround(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dTrLinAxi::edgeGiveTransformationJacobian(), oofem::FEI2dQuadLinAxi::edgeGiveTransformationJacobian(), and oofem::FEI2dQuadQuadAxi::edgeGiveTransformationJacobian().
|
pure virtual |
Evaluates edge global coordinates from given local ones.
These derivatives are in global coordinate system (where the nodal coordinates are defined).
answer | Contains resulting global coordinates. |
iedge | Determines edge number. |
lcoords | Array containing (local) coordinates. |
cellgeo | Underlying cell geometry. |
Implemented in oofem::FEI2dQuadQuad, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, oofem::FEI2dQuadConst, and oofem::FEI2dTrConst.
Referenced by boundaryEdgeLocal2Global(), and boundaryLocal2Global().
|
virtual |
Computes the exact area.
cellgeo | Cell geometry for the element. |
Reimplemented in oofem::FEI2dQuadQuad, oofem::FEI2dTrLin, oofem::FEI2dTrQuad, oofem::FEI2dLineHermite, oofem::FEI2dLineQuad, oofem::FEI2dLineLin, and oofem::FEI2dQuadLin.
Definition at line 60 of file feinterpol2d.C.
References OOFEM_ERROR.
Referenced by oofem::Element::computeArea().
|
inlinevirtual |
Returns a characteristic length of the geometry, typically a diagonal or edge length.
cellgeo | Underlying cell geometry. |
Reimplemented in oofem::FEI2dQuadQuad.
Definition at line 67 of file feinterpol2d.h.
Referenced by global2local().
|
virtual |
Gives the jacobian matrix at the local coordinates.
jacobianMatrix | The requested matrix. |
lcoords | Local coordinates. |
cellgeo | Element geometry. |
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::FEI2dLineQuad.
Definition at line 112 of file feinterpol2d.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FEInterpolation::evaldNdxi(), oofem::FloatMatrix::giveNumberOfRows(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FloatMatrix::resize(), xind, yind, and oofem::FloatMatrix::zero().
Referenced by oofem::Quad1MindlinShell3D::computeBmatrixAt(), and global2local().
|
inlinevirtual |
Returns number of spatial dimensions.
Implements oofem::FEInterpolation.
Definition at line 53 of file feinterpol2d.h.
|
virtual |
Default implementation using Newton's method to find the local coordinates.
Can be overloaded if desired.
Implements oofem::FEInterpolation.
Reimplemented in oofem::FEI2dLineHermite, oofem::FEI2dLineQuad, oofem::FEI2dLineLin, oofem::FEI2dQuadLin, oofem::FEI2dTrLin, and oofem::FEI2dTrConst.
Definition at line 68 of file feinterpol2d.C.
References oofem::FloatArray::add(), oofem::FloatArray::computeNorm(), giveCharacteristicLength(), giveJacobianMatrixAt(), inside(), oofem::FEInterpolation::local2global(), OOFEM_WARNING, oofem::FloatArray::resize(), oofem::FloatMatrix::solveForRhs(), and oofem::FloatArray::zero().
Referenced by oofem::TR21_2D_SUPG::LS_PCS_computeVOFFractions().
|
virtual |
Reimplemented in oofem::FEI2dQuadQuad, oofem::FEI2dTrLin, oofem::FEI2dQuadLin, oofem::FEI2dTrQuad, and oofem::FEI2dQuadConst.
Definition at line 134 of file feinterpol2d.C.
References OOFEM_ERROR.
Referenced by oofem::FEI2dTrConst::global2local(), and global2local().
|
protected |
Definition at line 48 of file feinterpol2d.h.
Referenced by oofem::FEI2dTrConst::edgeComputeLength(), oofem::FEI2dQuadConst::edgeComputeLength(), oofem::FEI2dTrLin::edgeComputeLength(), oofem::FEI2dQuadLin::edgeComputeLength(), oofem::FEI2dLineLin::edgeComputeLength(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dLineQuad::edgeEvaldNds(), oofem::FEI2dLineLin::edgeEvaldNds(), oofem::FEI2dQuadLin::edgeEvalNormal(), oofem::FEI2dTrQuad::edgeEvalNormal(), oofem::FEI2dTrLin::edgeEvalNormal(), oofem::FEI2dLineHermite::edgeEvalNormal(), oofem::FEI2dLineQuad::edgeEvalNormal(), oofem::FEI2dLineLin::edgeEvalNormal(), oofem::FEI2dQuadQuad::edgeEvalNormal(), oofem::FEI2dTrConst::edgeLocal2global(), oofem::FEI2dTrQuad::edgeLocal2global(), oofem::FEI2dQuadLin::edgeLocal2global(), oofem::FEI2dTrLin::edgeLocal2global(), oofem::FEI2dQuadQuad::edgeLocal2global(), oofem::FEI2dTrQuad::evald2Ndx2(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dQuadLin::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dLineHermite::evaldNdx(), oofem::FEI2dTrQuad::evalNXIntegral(), oofem::FEI2dQuadLin::evalNXIntegral(), oofem::FEI2dTrLin::evalNXIntegral(), oofem::FEI2dLineQuad::evalNXIntegral(), oofem::FEI2dLineLin::evalNXIntegral(), oofem::FEI2dQuadQuad::evalNXIntegral(), oofem::FEI2dQuadLin::giveArea(), oofem::FEI2dQuadQuad::giveArea(), giveJacobianMatrixAt(), oofem::FEI2dLineQuad::giveJacobianMatrixAt(), oofem::FEI2dLineHermite::giveLength(), oofem::FEI2dTrConst::giveTransformationJacobian(), oofem::FEI2dTrLin::giveTransformationJacobian(), oofem::FEI2dLineHermite::giveTransformationJacobian(), oofem::FEI2dLineLin::giveTransformationJacobian(), oofem::FEI2dLineQuad::giveTransformationJacobian(), oofem::FEI2dTrConst::global2local(), oofem::FEI2dTrLin::global2local(), oofem::FEI2dLineLin::global2local(), oofem::FEI2dQuadLin::global2local(), oofem::FEI2dLineQuad::global2local(), oofem::FEI2dLineHermite::global2local(), oofem::FEI2dTrConst::local2global(), oofem::FEI2dQuadConst::local2global(), oofem::FEI2dTrLin::local2global(), oofem::FEI2dTrQuad::local2global(), oofem::FEI2dQuadLin::local2global(), oofem::FEI2dLineLin::local2global(), oofem::FEI2dLineQuad::local2global(), oofem::FEI2dLineHermite::local2global(), and oofem::FEI2dQuadQuad::local2global().
|
protected |
Definition at line 48 of file feinterpol2d.h.
Referenced by oofem::FEI2dTrConst::edgeComputeLength(), oofem::FEI2dQuadConst::edgeComputeLength(), oofem::FEI2dTrLin::edgeComputeLength(), oofem::FEI2dQuadLin::edgeComputeLength(), oofem::FEI2dLineLin::edgeComputeLength(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dLineQuad::edgeEvaldNds(), oofem::FEI2dLineLin::edgeEvaldNds(), oofem::FEI2dQuadLin::edgeEvalNormal(), oofem::FEI2dTrQuad::edgeEvalNormal(), oofem::FEI2dTrLin::edgeEvalNormal(), oofem::FEI2dLineHermite::edgeEvalNormal(), oofem::FEI2dLineQuad::edgeEvalNormal(), oofem::FEI2dLineLin::edgeEvalNormal(), oofem::FEI2dQuadQuad::edgeEvalNormal(), oofem::FEI2dTrConst::edgeLocal2global(), oofem::FEI2dTrQuad::edgeLocal2global(), oofem::FEI2dQuadLin::edgeLocal2global(), oofem::FEI2dTrLin::edgeLocal2global(), oofem::FEI2dQuadQuad::edgeLocal2global(), oofem::FEI2dTrQuad::evald2Ndx2(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dQuadLin::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dLineHermite::evaldNdx(), oofem::FEI2dTrQuad::evalNXIntegral(), oofem::FEI2dQuadLin::evalNXIntegral(), oofem::FEI2dTrLin::evalNXIntegral(), oofem::FEI2dLineQuad::evalNXIntegral(), oofem::FEI2dLineLin::evalNXIntegral(), oofem::FEI2dQuadQuad::evalNXIntegral(), oofem::FEI2dQuadLin::giveArea(), oofem::FEI2dQuadQuad::giveArea(), giveJacobianMatrixAt(), oofem::FEI2dLineQuad::giveJacobianMatrixAt(), oofem::FEI2dLineHermite::giveLength(), oofem::FEI2dTrConst::giveTransformationJacobian(), oofem::FEI2dTrLin::giveTransformationJacobian(), oofem::FEI2dLineHermite::giveTransformationJacobian(), oofem::FEI2dLineLin::giveTransformationJacobian(), oofem::FEI2dLineQuad::giveTransformationJacobian(), oofem::FEI2dTrConst::global2local(), oofem::FEI2dTrLin::global2local(), oofem::FEI2dLineLin::global2local(), oofem::FEI2dQuadLin::global2local(), oofem::FEI2dLineQuad::global2local(), oofem::FEI2dLineHermite::global2local(), oofem::FEI2dTrConst::local2global(), oofem::FEI2dQuadConst::local2global(), oofem::FEI2dTrLin::local2global(), oofem::FEI2dTrQuad::local2global(), oofem::FEI2dQuadLin::local2global(), oofem::FEI2dLineLin::local2global(), oofem::FEI2dLineQuad::local2global(), oofem::FEI2dLineHermite::local2global(), and oofem::FEI2dQuadQuad::local2global().