|
OOFEM 3.0
|
#include <feitspline.h>
Public Member Functions | |
| TSplineInterpolation (int nsd) | |
| void | initializeFrom (InputRecord &ir, ParameterManager &pm, int elnum, int priority) override |
| Initializes receiver according to object description stored in input record. | |
| void | setNumberOfControlPoints (int num) |
| void | evalN (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| double | evaldNdx (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| void | local2global (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| int | global2local (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| void | giveJacobianMatrixAt (FloatMatrix &jacobianMatrix, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| int | giveKnotSpanBasisFuncMask (const IntArray &knotSpan, IntArray &mask) const override |
| int | giveNumberOfKnotSpanBasisFunctions (const IntArray &knotSpan) const override |
| const char * | giveClassName () const |
| Public Member Functions inherited from oofem::BSplineInterpolation | |
| BSplineInterpolation (int nsd) | |
| integrationDomain | giveIntegrationDomain (const Element_Geometry_Type egt) const override |
| const Element_Geometry_Type | giveGeometryType () const override |
| const Element_Geometry_Type | giveBoundaryGeometryType (int boundary) const override |
| integrationDomain | giveBoundaryIntegrationDomain (int ib, const Element_Geometry_Type) const override |
| Returns boundary integration domain. | |
| integrationDomain | giveBoundarySurfaceIntegrationDomain (int isurf, const Element_Geometry_Type) const override |
| Returns boundary integration domain. | |
| integrationDomain | giveBoundaryEdgeIntegrationDomain (int iedge, const Element_Geometry_Type) const override |
| Returns boundary integration domain. | |
| int | giveNsd (const Element_Geometry_Type) const override |
| void | postInitialize (ParameterManager &pm, int elnum) override |
| 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 | 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 |
| double | boundaryEdgeEvalNormal (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const override |
| 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 |
| 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 |
| int | giveNumberOfKnotSpans (int dim) const override |
| virtual int | giveNumberOfControlPoints (int dim) const |
| const FloatArray * | giveKnotVector () const override |
| const IntArray * | giveKnotMultiplicity (int dim) const override |
| const FloatArray * | giveKnotValues (int dim) const override |
| const char * | giveClassName () const |
| bool | hasSubPatchFormulation () const override |
| std::unique_ptr< IntegrationRule > | giveIntegrationRule (int order, const Element_Geometry_Type) const override |
| std::unique_ptr< IntegrationRule > | giveBoundaryIntegrationRule (int order, int boundary, Element_Geometry_Type) const override |
| std::unique_ptr< IntegrationRule > | giveBoundaryEdgeIntegrationRule (int order, int boundary, Element_Geometry_Type) const override |
| Public Member Functions inherited from oofem::FEInterpolation | |
| FEInterpolation (int o) | |
| virtual | ~FEInterpolation ()=default |
| 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 | 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 double | giveTransformationJacobian (const FloatArray &lcoords, const FEICellGeometry &cellgeo) const |
| 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 void | surfaceEvaldNdxi (FloatMatrix &answer, const FloatArray &lcoords) const |
| virtual void | surfaceEvald2Ndxi2 (FloatMatrix &answer, const FloatArray &lcoords) 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 Member Functions | |
| double | basisFunction (double u, int p, const FloatArray &U, const IntArray &I) const |
| void | dersBasisFunction (int n, double u, int p, const FloatArray &U, const IntArray &I, FloatArray &ders) const |
| void | createLocalKnotVector (FloatArray &answer, int p, const FloatArray &U, const IntArray &I, int &prepend, int &append) const |
| int | giveKnotSpanBasisFuncMask (const IntArray &startKnotSpan, const IntArray &endKnotSpan, IntArray &mask) const |
| int | giveNumberOfKnotSpanBasisFunctions (const IntArray &startKnotSpan, const IntArray &endKnotSpan) const |
| void | basisFuns (FloatArray &N, int span, double u, int p, const FloatArray &U) const |
| void | dersBasisFuns (int n, double u, int span, int p, const FloatArray &U, FloatMatrix &ders) const |
| int | findSpan (int n, int p, double u, const FloatArray &U) const |
| void | giveNonzeroBasisFuncInterval (int span, int deg, int &s, int &e) const |
Protected Attributes | |
| std::vector< std::array< IntArray, 3 > > | localIndexKnotVector |
| Local index knot vector of the dimensions [totalNumberOfControlPoints][nsd][degree+2]. | |
| int | totalNumberOfControlPoints = 0 |
| FloatArray | weights |
| Protected Attributes inherited from oofem::BSplineInterpolation | |
| int | nsd |
| Number of spatial directions. | |
| std::array< int, 3 > | degree |
| Degree in each direction. | |
| std::array< FloatArray, 3 > | knotValues |
| Knot values [nsd]. | |
| std::array< IntArray, 3 > | knotMultiplicity |
| Knot multiplicity [nsd]. | |
| std::array< int, 3 > | numberOfControlPoints |
| numberOfControlPoints[nsd] | |
| std::array< FloatArray, 3 > | knotVector |
| Knot vectors [nsd][knot_vector_size]. | |
| std::array< int, 3 > | numberOfKnotSpans |
| Nonzero spans in each directions [nsd]. | |
| Protected Attributes inherited from oofem::FEInterpolation | |
| int | order = 0 |
Additional Inherited Members | |
| Static Protected Attributes inherited from oofem::BSplineInterpolation | |
| static ParamKey | IPK_BSplineInterpolation_degree |
| static ParamKey | IPK_BSplineInterpolation_knotVectorU |
| static ParamKey | IPK_BSplineInterpolation_knotVectorV |
| static ParamKey | IPK_BSplineInterpolation_knotVectorW |
| static ParamKey | IPK_BSplineInterpolation_knotMultiplicityU |
| static ParamKey | IPK_BSplineInterpolation_knotMultiplicityV |
| static ParamKey | IPK_BSplineInterpolation_knotMultiplicityW |
Interpolation for T-splines.
Definition at line 62 of file feitspline.h.
|
inline |
Temporary open local knot vector to enable use of BSpline algorithms (common for all directions) [3*max_degree+2].
Definition at line 74 of file feitspline.h.
References oofem::BSplineInterpolation::BSplineInterpolation(), and oofem::BSplineInterpolation::nsd.
|
protected |
Evaluates the middle basis function on local knot vector at u.
| u | Value at which to evaluate. |
| p | Degree. |
| U | Global knot values. |
| I | Local index knot vector. |
Definition at line 607 of file feitspline.C.
References createLocalKnotVector(), and N.
Referenced by evalN(), and local2global().
|
protected |
Creates local open knot vector. This is generally done extracting knot values from global knot vector using the local index knot vector and by prepending p times the first knot and appending p times the last knot. However, existing knot multiplicity at the start and end must be accounted for.
| p | Degree. |
| U | Global knot values. |
| I | Local index knot vector |
| prepend | Number of prepended entries |
| append | Number of appended entries |
Definition at line 643 of file feitspline.C.
References oofem::FloatArray::at(), oofem::BSplineInterpolation::degree, oofem::BSplineInterpolation::nsd, and oofem::FloatArray::resize().
Referenced by basisFunction(), and dersBasisFunction().
|
protected |
Computes the middle basis function and it derivatives on local knot vector at u. The result is stored in the ders vector
\begin{align*}N(u) &= ders(0) N'(u) &= ders(1) N''(u)&= ders(2) \end{align*}
| n | Degree of the derivation. |
| u | Parametric value. |
| p | Degree. |
| U | Global knot values. |
| I | Local index knot vector. |
| ders | Vector containing the derivatives of the middle basis function. |
Definition at line 624 of file feitspline.C.
References createLocalKnotVector(), and oofem::FloatArray::resize().
Referenced by evaldNdx(), and giveJacobianMatrixAt().
|
overridevirtual |
Evaluates the matrix of derivatives of 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 dNi/dxj. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying cell geometry. |
Reimplemented from oofem::BSplineInterpolation.
Definition at line 157 of file feitspline.C.
References oofem::BSplineInterpolation::degree, dersBasisFunction(), oofem::BSplineInterpolation::findSpan(), oofem::FloatMatrix::giveDeterminant(), giveKnotSpanBasisFuncMask(), oofem::BSplineInterpolation::giveKnotValues(), oofem::IntArray::giveSize(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, oofem::BSplineInterpolation::numberOfControlPoints, OOFEM_ERROR, oofem::product(), oofem::FloatMatrix::resize(), weights, and oofem::FloatMatrix::zero().
|
overridevirtual |
Evaluates the array of interpolation functions (shape functions) at given point.
| answer | Contains resulting array of evaluated interpolation functions. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying cell geometry. |
Reimplemented from oofem::BSplineInterpolation.
Definition at line 116 of file feitspline.C.
References basisFunction(), oofem::BSplineInterpolation::degree, oofem::BSplineInterpolation::findSpan(), giveKnotSpanBasisFuncMask(), oofem::BSplineInterpolation::giveKnotValues(), oofem::IntArray::giveSize(), oofem::FEIIGAElementGeometryWrapper::knotSpan, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, N, oofem::BSplineInterpolation::nsd, oofem::BSplineInterpolation::numberOfControlPoints, OOFEM_ERROR, oofem::FloatArray::resize(), oofem::sum(), oofem::FloatArray::times(), and weights.
|
inline |
Definition at line 90 of file feitspline.h.
|
overridevirtual |
Gives the jacobian matrix at the local coordinates.
| jacobianMatrix | The requested matrix. |
| lcoords | Local coordinates. |
| cellgeo | Element geometry. |
Reimplemented from oofem::BSplineInterpolation.
Definition at line 340 of file feitspline.C.
References oofem::BSplineInterpolation::degree, dersBasisFunction(), oofem::BSplineInterpolation::findSpan(), giveKnotSpanBasisFuncMask(), oofem::BSplineInterpolation::giveKnotValues(), oofem::IntArray::giveSize(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, oofem::BSplineInterpolation::numberOfControlPoints, OOFEM_ERROR, oofem::product(), oofem::FloatMatrix::resize(), weights, and oofem::FloatMatrix::zero().
|
overridevirtual |
Returns indices (zero based) of nonzero basis functions for given knot span. The knot span identifies the sub-region of the finite element.
Reimplemented from oofem::BSplineInterpolation.
Definition at line 452 of file feitspline.C.
References oofem::BSplineInterpolation::degree, oofem::IntArray::followedBy(), oofem::BSplineInterpolation::knotValues, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, OOFEM_ERROR, oofem::IntArray::preallocate(), and totalNumberOfControlPoints.
Referenced by evaldNdx(), evalN(), giveJacobianMatrixAt(), and local2global().
|
protected |
Returns indices (zero based) of nonzero basis functions for given knot span interval.
Definition at line 531 of file feitspline.C.
References oofem::BSplineInterpolation::degree, oofem::IntArray::followedBy(), oofem::BSplineInterpolation::knotValues, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, OOFEM_ERROR, oofem::IntArray::preallocate(), and totalNumberOfControlPoints.
|
overridevirtual |
Returns the number of nonzero basis functions at individual knot span,
Reimplemented from oofem::BSplineInterpolation.
Definition at line 497 of file feitspline.C.
References oofem::BSplineInterpolation::degree, oofem::BSplineInterpolation::knotValues, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, and totalNumberOfControlPoints.
|
protected |
Returns the number of nonzero basis functions at given knot span interval.
Definition at line 576 of file feitspline.C.
References oofem::BSplineInterpolation::degree, oofem::BSplineInterpolation::knotValues, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, oofem::BSplineInterpolation::nsd, and totalNumberOfControlPoints.
|
inlineoverridevirtual |
Evaluates local coordinates from given global ones. If local coordinates cannot be found (generate elements, or point far outside geometry, then the center coordinate will be used as a last resort, and the return value will be zero.
| answer | Contains evaluated local coordinates. |
| gcoords | Array containing global coordinates. |
| cellgeo | Underlying cell geometry. |
Reimplemented from oofem::BSplineInterpolation.
Definition at line 82 of file feitspline.h.
References OOFEM_ERROR.
|
overridevirtual |
Initializes receiver according to object description stored in input record.
Reimplemented from oofem::BSplineInterpolation.
Definition at line 46 of file feitspline.C.
References _IFT_TSplineInterpolation_localIndexKnotVectorU, _IFT_TSplineInterpolation_localIndexKnotVectorV, _IFT_TSplineInterpolation_localIndexKnotVectorW, _IFT_TSplineInterpolation_weights, oofem::IntArray::clear(), oofem::BSplineInterpolation::degree, oofem::IntArray::giveSize(), IR_GIVE_FIELD, oofem::BSplineInterpolation::knotValues, localIndexKnotVector, oofem::BSplineInterpolation::nsd, totalNumberOfControlPoints, and weights.
Referenced by oofem::IGATSplineElement::initializeFrom().
|
overridevirtual |
Evaluates global coordinates from given local ones.
| answer | Contains resulting global coordinates. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying cell geometry. |
Reimplemented from oofem::BSplineInterpolation.
Definition at line 289 of file feitspline.C.
References basisFunction(), oofem::BSplineInterpolation::degree, oofem::BSplineInterpolation::findSpan(), giveKnotSpanBasisFuncMask(), oofem::BSplineInterpolation::giveKnotValues(), oofem::IntArray::giveSize(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, oofem::BSplineInterpolation::knotVector, localIndexKnotVector, N, oofem::BSplineInterpolation::nsd, oofem::BSplineInterpolation::numberOfControlPoints, OOFEM_ERROR, oofem::product(), oofem::FloatArray::resize(), oofem::FloatArray::times(), weights, and oofem::FloatArray::zero().
|
inline |
Definition at line 78 of file feitspline.h.
Referenced by oofem::IGATSplineElement::initializeFinish().
|
protected |
Local index knot vector of the dimensions [totalNumberOfControlPoints][nsd][degree+2].
Definition at line 66 of file feitspline.h.
Referenced by evaldNdx(), evalN(), giveJacobianMatrixAt(), giveKnotSpanBasisFuncMask(), giveKnotSpanBasisFuncMask(), giveNumberOfKnotSpanBasisFunctions(), giveNumberOfKnotSpanBasisFunctions(), initializeFrom(), and local2global().
|
protected |
Definition at line 67 of file feitspline.h.
Referenced by giveKnotSpanBasisFuncMask(), giveKnotSpanBasisFuncMask(), giveNumberOfKnotSpanBasisFunctions(), giveNumberOfKnotSpanBasisFunctions(), and initializeFrom().
|
protected |
Definition at line 68 of file feitspline.h.
Referenced by evaldNdx(), evalN(), giveJacobianMatrixAt(), initializeFrom(), and local2global().