OOFEM 3.0
Loading...
Searching...
No Matches
oofem::FEInterpolation2d Class Referenceabstract

#include <feinterpol2d.h>

Inheritance diagram for oofem::FEInterpolation2d:
Collaboration diagram for oofem::FEInterpolation2d:

Public Member Functions

 FEInterpolation2d (int o, int ind1, int ind2)
int giveNsd (const Element_Geometry_Type) const override
virtual double giveArea (const FEICellGeometry &cellgeo) const
virtual double giveCharacteristicLength (const FEICellGeometry &cellgeo) const
int global2local (FloatArray &answer, const FloatArray &gcoords, const FEICellGeometry &cellgeo) const override
void giveJacobianMatrixAt (FloatMatrix &jacobianMatrix, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
virtual bool inside (const FloatArray &lcoords) const
Boundary interpolation services.

Boundary is defined as entity of one dimension lower than the interpolation represents

IntArray boundaryEdgeGiveNodes (int boundary, const Element_Geometry_Type, bool includeHierarchical=false) const override
void boundaryEdgeEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundaryEdgeEvalNormal (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundaryEdgeGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
void boundaryEdgeLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
IntArray boundaryGiveNodes (int boundary, const Element_Geometry_Type) const override
void boundaryEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundaryEvalNormal (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundaryGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
void boundaryLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
Surface interpolation services
void boundarySurfaceEvalN (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
void boundarySurfaceEvaldNdx (FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundarySurfaceEvalNormal (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
void boundarySurfaceLocal2global (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
double boundarySurfaceGiveTransformationJacobian (int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override
IntArray boundarySurfaceGiveNodes (int boundary, const Element_Geometry_Type, bool includeHierarchical=false) const override
virtual FloatArrayF< 2 > surfaceEvalBaseVectorsAt (int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
virtual void surfaceEvaldNdxi (FloatMatrix &answer, const FloatArray &lcoords) const override
virtual void surfaceEvald2Ndxi2 (FloatMatrix &answer, const FloatArray &lcoords) const override
Edge interpolation services.
virtual IntArray computeLocalEdgeMapping (int iedge) const =0
IntArray computeEdgeMapping (const IntArray &elemNodes, int iedge) const
virtual void edgeEvalN (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double edgeEvalNormal (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void edgeEvaldNds (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void edgeLocal2global (FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double edgeGiveTransformationJacobian (int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
Public Member Functions inherited from oofem::FEInterpolation
 FEInterpolation (int o)
virtual ~FEInterpolation ()=default
virtual void initializeFrom (InputRecord &ir, ParameterManager &pm, int elnum, int priority)
 Initializes receiver according to object description stored in input record.
virtual void postInitialize (ParameterManager &pm, int elnum)
virtual integrationDomain giveIntegrationDomain (const Element_Geometry_Type) const =0
virtual const Element_Geometry_Type giveGeometryType () const =0
int giveInterpolationOrder () const
virtual void giveCellDofMans (IntArray &nodes, IntArray &internalDofMans, Element *elem) const
 Returns list of element nodes (and list of internal dof managers) (including on edges and surfaces) defining the approximation.
virtual void evalN (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double evaldNdx (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual void evald2Ndx2 (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
virtual void evaldNdxi (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
virtual void giveLocalNodeCoords (FloatMatrix &answer, const Element_Geometry_Type) const
virtual void local2global (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
virtual double giveTransformationJacobian (const FloatArray &lcoords, const FEICellGeometry &cellgeo) const
virtual std::unique_ptr< IntegrationRulegiveIntegrationRule (int order, const Element_Geometry_Type) const
virtual integrationDomain giveBoundaryEdgeIntegrationDomain (int boundary, const Element_Geometry_Type) const =0
 Returns boundary integration domain.
virtual std::unique_ptr< IntegrationRulegiveBoundaryEdgeIntegrationRule (int order, int boundary, const Element_Geometry_Type) const
virtual integrationDomain giveBoundarySurfaceIntegrationDomain (int boundary, const Element_Geometry_Type) const =0
 Returns boundary integration domain.
virtual std::unique_ptr< IntegrationRulegiveBoundarySurfaceIntegrationRule (int order, int boundary, const Element_Geometry_Type) const
virtual double evalNXIntegral (int boundary, const FEICellGeometry &cellgeo) const
virtual integrationDomain giveBoundaryIntegrationDomain (int boundary, const Element_Geometry_Type) const =0
 Returns boundary integration domain.
virtual std::unique_ptr< IntegrationRulegiveBoundaryIntegrationRule (int order, int boundary, const Element_Geometry_Type) const
virtual const Element_Geometry_Type giveBoundaryGeometryType (int boundary) const =0
virtual int giveKnotSpanBasisFuncMask (const IntArray &knotSpan, IntArray &mask) const
virtual int giveNumberOfKnotSpanBasisFunctions (const IntArray &knotSpan) const
virtual bool hasSubPatchFormulation () const
virtual const FloatArraygiveKnotVector () const
virtual int giveNumberOfKnotSpans (int dim) const
virtual const FloatArraygiveKnotValues (int dim) const
virtual const IntArraygiveKnotMultiplicity (int dim) const
virtual int giveNumberOfEdges (const Element_Geometry_Type) const
virtual int giveNumberOfNodes (const Element_Geometry_Type) const
virtual void initializeCell (Element *e) const
std::string errorInfo (const char *func) const

Protected Attributes

int xind
int yind
Protected Attributes inherited from oofem::FEInterpolation
int order = 0

Detailed Description

Class representing a general abstraction for surface finite element interpolation class.

Definition at line 45 of file feinterpol2d.h.

Constructor & Destructor Documentation

◆ FEInterpolation2d()

Member Function Documentation

◆ boundaryEdgeEvalN()

void oofem::FEInterpolation2d::boundaryEdgeEvalN ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
answerBasis functions Array to be filled with the boundary nodes.
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.
Todo

Implements oofem::FEInterpolation.

Definition at line 46 of file feinterpol2d.C.

References edgeEvalN().

◆ boundaryEdgeEvalNormal()

double oofem::FEInterpolation2d::boundaryEdgeEvalNormal ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates the normal out of the edge at given point.

Parameters
answerContains resulting normal vector.
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Surface mapping jacobian.

Implements oofem::FEInterpolation.

Definition at line 61 of file feinterpol2d.C.

References edgeEvalNormal().

◆ boundaryEdgeGiveNodes()

IntArray oofem::FEInterpolation2d::boundaryEdgeGiveNodes ( int boundary,
const Element_Geometry_Type ,
bool includeHierarchical = false ) const
overridevirtual

Gives the boundary nodes for requested boundary number.

Parameters
answerArray to be filled with the boundary nodes.
boundaryBoundary number.
includeHierarchicalIf true, include hierarchical nodes, introduced by interpolations on universal cells (mpm)

Implements oofem::FEInterpolation.

Definition at line 41 of file feinterpol2d.C.

References computeLocalEdgeMapping().

◆ boundaryEdgeGiveTransformationJacobian()

double oofem::FEInterpolation2d::boundaryEdgeGiveTransformationJacobian ( int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.
Returns
The determinant of the boundary transformation Jacobian.

Implements oofem::FEInterpolation.

Definition at line 51 of file feinterpol2d.C.

References edgeGiveTransformationJacobian().

◆ boundaryEdgeLocal2Global()

void oofem::FEInterpolation2d::boundaryEdgeLocal2Global ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
answerGlobal coordinates.
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 56 of file feinterpol2d.C.

References edgeLocal2global().

◆ boundaryEvalN()

void oofem::FEInterpolation2d::boundaryEvalN ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
answerBasis functions Array to be filled with the boundary nodes.
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 149 of file feinterpol2d.C.

References edgeEvalN().

◆ boundaryEvalNormal()

double oofem::FEInterpolation2d::boundaryEvalNormal ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates the normal on the requested boundary.

Parameters
answerThe evaluated normal.
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.
Returns
The boundary transformation Jacobian.

Implements oofem::FEInterpolation.

Definition at line 154 of file feinterpol2d.C.

References edgeEvalNormal().

◆ boundaryGiveNodes()

IntArray oofem::FEInterpolation2d::boundaryGiveNodes ( int boundary,
const Element_Geometry_Type  ) const
overridevirtual

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.

Parameters
answerArray to be filled with the boundary nodes.
boundaryBoundary number.

Implements oofem::FEInterpolation.

Definition at line 144 of file feinterpol2d.C.

References computeLocalEdgeMapping().

Referenced by oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement().

◆ boundaryGiveTransformationJacobian()

double oofem::FEInterpolation2d::boundaryGiveTransformationJacobian ( int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.
Returns
The determinant of the boundary transformation Jacobian.

Implements oofem::FEInterpolation.

Definition at line 159 of file feinterpol2d.C.

References edgeGiveTransformationJacobian().

◆ boundaryLocal2Global()

void oofem::FEInterpolation2d::boundaryLocal2Global ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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.

Parameters
answerGlobal coordinates.
boundaryBoundary number.
lcoordsThe local coordinates (on the boundary local coordinate system).
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 164 of file feinterpol2d.C.

References edgeLocal2global().

◆ boundarySurfaceEvaldNdx()

void oofem::FEInterpolation2d::boundarySurfaceEvaldNdx ( FloatMatrix & answer,
int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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).

Parameters
answerContains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi.
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 191 of file feinterpol2d.C.

References oofem::FEInterpolation::evaldNdx().

◆ boundarySurfaceEvalN()

void oofem::FEInterpolation2d::boundarySurfaceEvalN ( FloatArray & answer,
int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates the array of edge interpolation functions (shape functions) at given point.

Parameters
answerContains resulting array of evaluated interpolation functions.
isurfSurface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 186 of file feinterpol2d.C.

References oofem::FEInterpolation::evalN().

◆ boundarySurfaceEvalNormal()

double oofem::FEInterpolation2d::boundarySurfaceEvalNormal ( FloatArray & answer,
int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates the normal out of the surface at given point.

Parameters
answerContains resulting normal vector.
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Surface mapping jacobian.

Implements oofem::FEInterpolation.

Definition at line 197 of file feinterpol2d.C.

References oofem::FEInterpolation::giveTransformationJacobian(), and oofem::Vec3().

◆ boundarySurfaceGiveNodes()

IntArray oofem::FEInterpolation2d::boundarySurfaceGiveNodes ( int boundary,
const Element_Geometry_Type ,
bool includeHierarchical = false ) const
overridevirtual

Gives the boundary nodes for requested boundary number.

Parameters
answerArray to be filled with the boundary nodes.
boundaryBoundary number.
includeHierarchicalIf true, include hierarchical nodes, introduced by interpolations on universal cells (mpm)

Implements oofem::FEInterpolation.

Definition at line 247 of file feinterpol2d.C.

References oofem::IntArray::enumerate(), and oofem::FEInterpolation::giveNumberOfNodes().

◆ boundarySurfaceGiveTransformationJacobian()

double oofem::FEInterpolation2d::boundarySurfaceGiveTransformationJacobian ( int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates the edge jacobian of transformation between local and global coordinates.

Parameters
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Determinant of the transformation.

Implements oofem::FEInterpolation.

Definition at line 241 of file feinterpol2d.C.

References oofem::FEInterpolation::giveTransformationJacobian().

◆ boundarySurfaceLocal2global()

void oofem::FEInterpolation2d::boundarySurfaceLocal2global ( FloatArray & answer,
int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates edge global coordinates from given local ones. These derivatives are in global coordinate system (where the nodal coordinates are defined).

Parameters
answerContains resulting global coordinates.
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implements oofem::FEInterpolation.

Definition at line 235 of file feinterpol2d.C.

References oofem::FEInterpolation::local2global().

◆ computeEdgeMapping()

IntArray oofem::FEInterpolation2d::computeEdgeMapping ( const IntArray & elemNodes,
int iedge ) const

Definition at line 169 of file feinterpol2d.C.

References oofem::IntArray::at(), and computeLocalEdgeMapping().

Referenced by computeLocalEdgeMapping().

◆ computeLocalEdgeMapping()

◆ edgeEvaldNds()

virtual void oofem::FEInterpolation2d::edgeEvaldNds ( FloatArray & answer,
int iedge,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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).

Parameters
answerContains resulting array of derivatives, the member at i position contains value of \( \frac{\mathrm{d}N_i}{\mathrm{d}s} \).
iedgeDetermines the edge number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrConst, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.

Referenced by oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement(), and oofem::SurfaceTensionBoundaryCondition::computeTangentFromElement().

◆ edgeEvalN()

virtual void oofem::FEInterpolation2d::edgeEvalN ( FloatArray & answer,
int iedge,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
pure virtual

Evaluates the array of edge interpolation functions (shape functions) at given point.

Parameters
answerContains resulting array of evaluated interpolation functions.
iedgeEdge number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrConst, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.

Referenced by boundaryEdgeEvalN(), and boundaryEvalN().

◆ edgeEvalNormal()

virtual double oofem::FEInterpolation2d::edgeEvalNormal ( FloatArray & answer,
int iedge,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
pure virtual

Evaluates the normal on the given edge.

Parameters
answerContains the evaluated normal.
iedgeDetermines the edge number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrConst, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.

Referenced by boundaryEdgeEvalNormal(), boundaryEvalNormal(), and edgeGiveTransformationJacobian().

◆ edgeGiveTransformationJacobian()

double oofem::FEInterpolation2d::edgeGiveTransformationJacobian ( int iedge,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
virtual

Evaluates the edge Jacobian of transformation between local and global coordinates.

Parameters
iedgeDetermines edge number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Determinant of the mapping on the given edge.

Reimplemented in oofem::FEI2dQuadLinAxi, oofem::FEI2dQuadQuadAxi, and oofem::FEI2dTrLinAxi.

Definition at line 180 of file feinterpol2d.C.

References edgeEvalNormal().

Referenced by boundaryEdgeGiveTransformationJacobian(), boundaryGiveTransformationJacobian(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dQuadLinAxi::edgeGiveTransformationJacobian(), oofem::FEI2dQuadQuadAxi::edgeGiveTransformationJacobian(), oofem::FEI2dTrLinAxi::edgeGiveTransformationJacobian(), and edgeLocal2global().

◆ edgeLocal2global()

virtual void oofem::FEInterpolation2d::edgeLocal2global ( FloatArray & answer,
int iedge,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
pure virtual

Evaluates edge global coordinates from given local ones. These derivatives are in global coordinate system (where the nodal coordinates are defined).

Parameters
answerContains resulting global coordinates.
iedgeDetermines edge number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.

Implemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrConst, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.

References edgeGiveTransformationJacobian().

Referenced by boundaryEdgeLocal2Global(), and boundaryLocal2Global().

◆ giveArea()

double oofem::FEInterpolation2d::giveArea ( const FEICellGeometry & cellgeo) const
virtual

Computes the exact area.

Parameters
cellgeoCell geometry for the element.
Returns
Area of geometry.

Reimplemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.

Definition at line 67 of file feinterpol2d.C.

References OOFEM_ERROR.

Referenced by oofem::Element::computeArea(), and giveCharacteristicLength().

◆ giveCharacteristicLength()

virtual double oofem::FEInterpolation2d::giveCharacteristicLength ( const FEICellGeometry & cellgeo) const
inlinevirtual

Returns a characteristic length of the geometry, typically a diagonal or edge length.

Parameters
cellgeoUnderlying cell geometry.
Returns
Square root of area.

Reimplemented in oofem::FEI2dQuadQuad.

Definition at line 67 of file feinterpol2d.h.

References giveArea().

Referenced by global2local().

◆ giveJacobianMatrixAt()

void oofem::FEInterpolation2d::giveJacobianMatrixAt ( FloatMatrix & jacobianMatrix,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Gives the jacobian matrix at the local coordinates.

Parameters
jacobianMatrixThe requested matrix.
lcoordsLocal coordinates.
cellgeoElement geometry.

Reimplemented from oofem::FEInterpolation.

Definition at line 118 of file feinterpol2d.C.

References oofem::FloatMatrix::at(), oofem::FEInterpolation::evaldNdxi(), oofem::FloatMatrix::giveNumberOfRows(), oofem::FloatMatrix::resize(), xind, and yind.

Referenced by global2local().

◆ giveNsd()

int oofem::FEInterpolation2d::giveNsd ( const Element_Geometry_Type ) const
inlineoverridevirtual

Returns number of spatial dimensions.

Implements oofem::FEInterpolation.

Definition at line 53 of file feinterpol2d.h.

◆ global2local()

int oofem::FEInterpolation2d::global2local ( FloatArray & answer,
const FloatArray & gcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

◆ inside()

bool oofem::FEInterpolation2d::inside ( const FloatArray & lcoords) const
virtual

◆ surfaceEvalBaseVectorsAt()

FloatArrayF< 2 > oofem::FEInterpolation2d::surfaceEvalBaseVectorsAt ( int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
virtual

Evaluates the tangent vectors of the surface at given point.

Parameters
isurfDetermines the surface number.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Surface mapping jacobian .
tangent vector

Definition at line 205 of file feinterpol2d.C.

References computeLocalEdgeMapping(), oofem::FEICellGeometry::giveVertexCoordinates(), and surfaceEvaldNdxi().

Referenced by oofem::FEContactSurface::computeContactPointLocalCoordinates_2d().

◆ surfaceEvald2Ndxi2()

void oofem::FEInterpolation2d::surfaceEvald2Ndxi2 ( FloatMatrix & answer,
const FloatArray & lcoords ) const
overridevirtual

Evaluates the matrix of second derivatives of surface interpolation functions (shape functions) wrt parametric coordinates at given point.

Parameters
answerContains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi.
lcoordsArray containing (local) coordinates.

Reimplemented from oofem::FEInterpolation.

Reimplemented in oofem::FEI2dLineLin.

Definition at line 226 of file feinterpol2d.C.

References OOFEM_ERROR.

Referenced by oofem::FEContactSurface::computeContactPointLocalCoordinates_2d().

◆ surfaceEvaldNdxi()

void oofem::FEInterpolation2d::surfaceEvaldNdxi ( FloatMatrix & answer,
const FloatArray & lcoords ) const
overridevirtual

Evaluates the matrix of derivatives of surface interpolation functions (shape functions) wrt parametric coordinates at given point.

Parameters
answerContains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi.
lcoordsArray containing (local) coordinates.

Reimplemented from oofem::FEInterpolation.

Reimplemented in oofem::FEI2dLineLin.

Definition at line 221 of file feinterpol2d.C.

References OOFEM_ERROR.

Referenced by surfaceEvalBaseVectorsAt().

Member Data Documentation

◆ xind

int oofem::FEInterpolation2d::xind
protected

Definition at line 48 of file feinterpol2d.h.

Referenced by oofem::FEI2dQuadBiQuad::_evaldNdx(), oofem::FEI2dLineConst::edgeComputeLength(), oofem::FEI2dLineLin::edgeComputeLength(), oofem::FEI2dQuadConst::edgeComputeLength(), oofem::FEI2dQuadLin::edgeComputeLength(), oofem::FEI2dTrConst::edgeComputeLength(), oofem::FEI2dTrLin::edgeComputeLength(), oofem::FEI2dLineLin::edgeEvaldNds(), oofem::FEI2dLineQuad::edgeEvaldNds(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dLineConst::edgeEvalNormal(), oofem::FEI2dLineHermite::edgeEvalNormal(), oofem::FEI2dLineLin::edgeEvalNormal(), oofem::FEI2dLineQuad::edgeEvalNormal(), oofem::FEI2dQuadLin::edgeEvalNormal(), oofem::FEI2dQuadQuad::edgeEvalNormal(), oofem::FEI2dTrLin::edgeEvalNormal(), oofem::FEI2dTrQuad::edgeEvalNormal(), oofem::FEI2dQuadLin::edgeLocal2global(), oofem::FEI2dQuadQuad::edgeLocal2global(), oofem::FEI2dTrConst::edgeLocal2global(), oofem::FEI2dTrLin::edgeLocal2global(), oofem::FEI2dTrQuad::edgeLocal2global(), oofem::FEI2dTrQuad::evald2Ndx2(), oofem::FEI2dLineHermite::evaldNdx(), oofem::FEI2dQuadLin::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dLineConst::evalNXIntegral(), oofem::FEI2dLineLin::evalNXIntegral(), oofem::FEI2dLineQuad::evalNXIntegral(), oofem::FEI2dQuadLin::evalNXIntegral(), oofem::FEI2dQuadQuad::evalNXIntegral(), oofem::FEI2dTrLin::evalNXIntegral(), oofem::FEI2dTrQuad::evalNXIntegral(), FEInterpolation2d(), oofem::FEI2dQuadLin::giveArea(), oofem::FEI2dQuadQuad::giveArea(), oofem::FEI2dTrLin::giveArea(), oofem::FEI2dTrQuad::giveArea(), oofem::FEI2dLineQuad::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::FEI2dLineHermite::giveLength(), oofem::FEI2dLineConst::giveTransformationJacobian(), oofem::FEI2dLineHermite::giveTransformationJacobian(), oofem::FEI2dLineLin::giveTransformationJacobian(), oofem::FEI2dLineQuad::giveTransformationJacobian(), oofem::FEI2dTrConst::giveTransformationJacobian(), oofem::FEI2dTrLin::giveTransformationJacobian(), oofem::FEI2dLineConst::global2local(), oofem::FEI2dLineHermite::global2local(), oofem::FEI2dLineLin::global2local(), oofem::FEI2dLineQuad::global2local(), oofem::FEI2dQuadLin::global2local(), oofem::FEI2dTrConst::global2local(), oofem::FEI2dTrLin::global2local(), oofem::FEI2dLineConst::local2global(), oofem::FEI2dLineHermite::local2global(), oofem::FEI2dLineLin::local2global(), oofem::FEI2dLineQuad::local2global(), oofem::FEI2dQuadConst::local2global(), oofem::FEI2dQuadLin::local2global(), oofem::FEI2dQuadQuad::local2global(), oofem::FEI2dTrConst::local2global(), oofem::FEI2dTrLin::local2global(), and oofem::FEI2dTrQuad::local2global().

◆ yind

int oofem::FEInterpolation2d::yind
protected

Definition at line 48 of file feinterpol2d.h.

Referenced by oofem::FEI2dQuadBiQuad::_evaldNdx(), oofem::FEI2dLineConst::edgeComputeLength(), oofem::FEI2dLineLin::edgeComputeLength(), oofem::FEI2dQuadConst::edgeComputeLength(), oofem::FEI2dQuadLin::edgeComputeLength(), oofem::FEI2dTrConst::edgeComputeLength(), oofem::FEI2dTrLin::edgeComputeLength(), oofem::FEI2dLineLin::edgeEvaldNds(), oofem::FEI2dLineQuad::edgeEvaldNds(), oofem::FEI2dTrQuad::edgeEvaldNds(), oofem::FEI2dLineConst::edgeEvalNormal(), oofem::FEI2dLineHermite::edgeEvalNormal(), oofem::FEI2dLineLin::edgeEvalNormal(), oofem::FEI2dLineQuad::edgeEvalNormal(), oofem::FEI2dQuadLin::edgeEvalNormal(), oofem::FEI2dQuadQuad::edgeEvalNormal(), oofem::FEI2dTrLin::edgeEvalNormal(), oofem::FEI2dTrQuad::edgeEvalNormal(), oofem::FEI2dQuadLin::edgeLocal2global(), oofem::FEI2dQuadQuad::edgeLocal2global(), oofem::FEI2dTrConst::edgeLocal2global(), oofem::FEI2dTrLin::edgeLocal2global(), oofem::FEI2dTrQuad::edgeLocal2global(), oofem::FEI2dTrQuad::evald2Ndx2(), oofem::FEI2dLineHermite::evaldNdx(), oofem::FEI2dQuadLin::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dQuadQuad::evaldNdx(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrLin::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dTrQuad::evaldNdx(), oofem::FEI2dLineConst::evalNXIntegral(), oofem::FEI2dLineLin::evalNXIntegral(), oofem::FEI2dLineQuad::evalNXIntegral(), oofem::FEI2dQuadLin::evalNXIntegral(), oofem::FEI2dQuadQuad::evalNXIntegral(), oofem::FEI2dTrLin::evalNXIntegral(), oofem::FEI2dTrQuad::evalNXIntegral(), FEInterpolation2d(), oofem::FEI2dQuadLin::giveArea(), oofem::FEI2dQuadQuad::giveArea(), oofem::FEI2dTrLin::giveArea(), oofem::FEI2dTrQuad::giveArea(), oofem::FEI2dLineQuad::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::FEI2dLineHermite::giveLength(), oofem::FEI2dLineConst::giveTransformationJacobian(), oofem::FEI2dLineHermite::giveTransformationJacobian(), oofem::FEI2dLineLin::giveTransformationJacobian(), oofem::FEI2dLineQuad::giveTransformationJacobian(), oofem::FEI2dTrConst::giveTransformationJacobian(), oofem::FEI2dTrLin::giveTransformationJacobian(), oofem::FEI2dLineConst::global2local(), oofem::FEI2dLineHermite::global2local(), oofem::FEI2dLineLin::global2local(), oofem::FEI2dLineQuad::global2local(), oofem::FEI2dQuadLin::global2local(), oofem::FEI2dTrConst::global2local(), oofem::FEI2dTrLin::global2local(), oofem::FEI2dLineConst::local2global(), oofem::FEI2dLineHermite::local2global(), oofem::FEI2dLineLin::local2global(), oofem::FEI2dLineQuad::local2global(), oofem::FEI2dQuadConst::local2global(), oofem::FEI2dQuadLin::local2global(), oofem::FEI2dQuadQuad::local2global(), oofem::FEI2dTrConst::local2global(), oofem::FEI2dTrLin::local2global(), and oofem::FEI2dTrQuad::local2global().


The documentation for this class was generated from the following files:

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011