|
OOFEM 3.0
|
#include <feinterpol2d.h>
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< IntegrationRule > | giveIntegrationRule (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< IntegrationRule > | giveBoundaryEdgeIntegrationRule (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< IntegrationRule > | giveBoundarySurfaceIntegrationRule (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< IntegrationRule > | giveBoundaryIntegrationRule (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 FloatArray * | giveKnotVector () const |
| virtual int | giveNumberOfKnotSpans (int dim) const |
| virtual const FloatArray * | giveKnotValues (int dim) const |
| virtual const IntArray * | giveKnotMultiplicity (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 |
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.
References oofem::FEInterpolation::FEInterpolation(), xind, and yind.
Referenced by oofem::FEI2dLineConst::FEI2dLineConst(), oofem::FEI2dLineHermite::FEI2dLineHermite(), oofem::FEI2dLineLin::FEI2dLineLin(), oofem::FEI2dLineQuad::FEI2dLineQuad(), oofem::FEI2dQuadConst::FEI2dQuadConst(), oofem::FEI2dQuadLin::FEI2dQuadLin(), oofem::FEI2dQuadQuad::FEI2dQuadQuad(), oofem::FEI2dTrConst::FEI2dTrConst(), oofem::FEI2dTrLin::FEI2dTrLin(), and oofem::FEI2dTrQuad::FEI2dTrQuad().
|
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.
| 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 46 of file feinterpol2d.C.
References edgeEvalN().
|
overridevirtual |
Evaluates the normal out of the edge 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 61 of file feinterpol2d.C.
References edgeEvalNormal().
|
overridevirtual |
Gives the boundary nodes for requested boundary number.
| answer | Array to be filled with the boundary nodes. |
| boundary | Boundary number. |
| includeHierarchical | If true, include hierarchical nodes, introduced by interpolations on universal cells (mpm) |
Implements oofem::FEInterpolation.
Definition at line 41 of file feinterpol2d.C.
References computeLocalEdgeMapping().
|
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.
| boundary | Boundary number. |
| lcoords | The local coordinates (on the boundary local coordinate system). |
| cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Definition at line 51 of file feinterpol2d.C.
References edgeGiveTransformationJacobian().
|
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.
| 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 56 of file feinterpol2d.C.
References edgeLocal2global().
|
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.
| 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 149 of file feinterpol2d.C.
References edgeEvalN().
|
overridevirtual |
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 154 of file feinterpol2d.C.
References edgeEvalNormal().
|
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.
| answer | Array to be filled with the boundary nodes. |
| boundary | Boundary number. |
Implements oofem::FEInterpolation.
Definition at line 144 of file feinterpol2d.C.
References computeLocalEdgeMapping().
Referenced by oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement().
|
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.
| 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 edgeGiveTransformationJacobian().
|
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.
| 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 164 of file feinterpol2d.C.
References edgeLocal2global().
|
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).
| 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 191 of file feinterpol2d.C.
References oofem::FEInterpolation::evaldNdx().
|
overridevirtual |
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 186 of file feinterpol2d.C.
References oofem::FEInterpolation::evalN().
|
overridevirtual |
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 197 of file feinterpol2d.C.
References oofem::FEInterpolation::giveTransformationJacobian(), and oofem::Vec3().
|
overridevirtual |
Gives the boundary nodes for requested boundary number.
| answer | Array to be filled with the boundary nodes. |
| boundary | Boundary number. |
| includeHierarchical | If 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().
|
overridevirtual |
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 241 of file feinterpol2d.C.
References oofem::FEInterpolation::giveTransformationJacobian().
|
overridevirtual |
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 235 of file feinterpol2d.C.
References oofem::FEInterpolation::local2global().
| 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().
|
pure virtual |
Implemented in oofem::FEI2dLineConst, oofem::FEI2dLineHermite, oofem::FEI2dLineLin, oofem::FEI2dLineQuad, oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrConst, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.
References computeEdgeMapping().
Referenced by boundaryEdgeGiveNodes(), boundaryGiveNodes(), computeEdgeMapping(), oofem::XfemElementInterface::partitionEdgeSegment(), and surfaceEvalBaseVectorsAt().
|
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 \( \frac{\mathrm{d}N_i}{\mathrm{d}s} \). |
| iedge | Determines the edge number. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying 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().
|
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::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().
|
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::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().
|
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::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().
|
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::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().
|
virtual |
Computes the exact area.
| cellgeo | Cell geometry for the element. |
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().
|
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.
References giveArea().
Referenced by global2local().
|
overridevirtual |
Gives the jacobian matrix at the local coordinates.
| jacobianMatrix | The requested matrix. |
| lcoords | Local coordinates. |
| cellgeo | Element 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().
|
inlineoverridevirtual |
Returns number of spatial dimensions.
Implements oofem::FEInterpolation.
Definition at line 53 of file feinterpol2d.h.
|
overridevirtual |
Default implementation using Newton's method to find the local coordinates. Can be overloaded if desired.
Implements oofem::FEInterpolation.
Definition at line 74 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(), oofem::Vec2(), and oofem::FloatArray::zero().
|
virtual |
Reimplemented in oofem::FEI2dQuadConst, oofem::FEI2dQuadLin, oofem::FEI2dQuadQuad, oofem::FEI2dTrLin, and oofem::FEI2dTrQuad.
Definition at line 139 of file feinterpol2d.C.
References OOFEM_ERROR.
Referenced by oofem::FEI2dTrConst::global2local(), and global2local().
|
virtual |
Evaluates the tangent vectors of the surface at given point.
| isurf | Determines the surface number. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying cell geometry. |
Definition at line 205 of file feinterpol2d.C.
References computeLocalEdgeMapping(), oofem::FEICellGeometry::giveVertexCoordinates(), and surfaceEvaldNdxi().
Referenced by oofem::FEContactSurface::computeContactPointLocalCoordinates_2d().
|
overridevirtual |
Evaluates the matrix of second derivatives of surface interpolation functions (shape functions) wrt parametric coordinates at given point.
| answer | Contains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi. |
| lcoords | Array 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().
|
overridevirtual |
Evaluates the matrix of derivatives of surface interpolation functions (shape functions) wrt parametric coordinates at given point.
| answer | Contains resulting matrix of derivatives, the member at i,j position contains value of dNj/dxi. |
| lcoords | Array 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().
|
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().
|
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().