|
OOFEM 3.0
|
#include <feibspline.h>
Public Member Functions | |
| 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 | initializeFrom (InputRecord &ir, ParameterManager &pm, int elnum, int priority) override |
| Initializes receiver according to object description stored in input record. | |
| 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 |
| 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 Attributes | |
| 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 |
Static Protected Attributes | |
| 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 |
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 |
| 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 |
| 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 |
| 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 |
| 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 |
Interpolation for B-splines.
Definition at line 62 of file feibspline.h.
|
inline |
Definition at line 93 of file feibspline.h.
References oofem::FEInterpolation::FEInterpolation(), and nsd.
Referenced by oofem::NURBSInterpolation::NURBSInterpolation(), and oofem::TSplineInterpolation::TSplineInterpolation().
|
protected |
Evaluates the nonvanishing basis functions of 1d BSpline (algorithm A2.2 from NURBS book)
| span | Knot span index (zero based). |
| u | Value at which to evaluate. |
| p | Degree. |
| U | Knot vector. |
| N | Computed p+1 nonvanishing functions (N_{span-p,p}-N_{span,p}) |
Definition at line 644 of file feibspline.C.
References N.
Referenced by evalN(), oofem::NURBSInterpolation::evalN(), local2global(), and oofem::NURBSInterpolation::local2global().
|
inlineoverridevirtual |
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 146 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 152 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 144 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 148 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 150 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 175 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 177 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 173 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 179 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 181 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 160 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 158 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 162 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 168 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 166 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
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 164 of file feibspline.h.
References OOFEM_ERROR.
|
protected |
Computes nonzero basis functions and their derivatives at u.
For information on the algorithm, see A2.3 on p72 of the NURBS book. The result is stored in the ders matrix, where ders is of size (n+1,p+1) and the derivative
\begin{align*}N(u) &= \mathit{ders}(0,span-p+j) \quad\text{where } j=0...p \ \ N'(u) &= \mathit{ders}(1,span-p+j) \quad\text{where } j=0...p \ \ N''(u) &= \mathit{ders}(2,span-p+j) \quad\text{where } j=0...p \end{align*}
| n | Degree of the derivation. |
| u | Parametric value. |
| span | Knot span index (zero based). |
| p | Degree. |
| U | Knot vector. |
| ders | Matrix containing the derivatives of the basis functions. |
Definition at line 675 of file feibspline.C.
References oofem::FloatMatrix::resize().
Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), giveJacobianMatrixAt(), and oofem::NURBSInterpolation::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. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 221 of file feibspline.C.
References degree, dersBasisFuns(), findSpan(), oofem::FloatMatrix::giveDeterminant(), giveNumberOfKnotSpanBasisFunctions(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, nsd, numberOfControlPoints, OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatArray::zero(), 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. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 176 of file feibspline.C.
References oofem::FloatArray::at(), basisFuns(), degree, findSpan(), giveNumberOfKnotSpanBasisFunctions(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, N, nsd, numberOfControlPoints, OOFEM_ERROR, and oofem::FloatArray::resize().
|
protected |
Determines the knot span index (Algorithm A2.1 from the NURBS book)
Determines the knot span for which there exists non-zero basis functions. The span is the index k for which the parameter u is valid in the (u_k,u_{k+1}] range.
| n | Number of control points - 1. |
| u | Parametric value. |
| p | Degree. |
| U | Knot vector. |
Definition at line 774 of file feibspline.C.
Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), local2global(), oofem::NURBSInterpolation::local2global(), and oofem::TSplineInterpolation::local2global().
|
inlineoverridevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 132 of file feibspline.h.
References oofem::_Line, oofem::_UnknownIntegrationDomain, and nsd.
|
inlineoverridevirtual |
Sets up a suitable integration rule for integrating over the requested boundary. The required polynomial order for the determinant of the jacobian is added automatically.
| order | Polynomial order of the integrand (should NOT including determinant of jacobian). |
| boundary | Boundary number. |
Reimplemented from oofem::FEInterpolation.
Definition at line 209 of file feibspline.h.
References OOFEM_ERROR, and oofem::FEInterpolation::order.
|
inlineoverridevirtual |
Returns boundary geometry type
Implements oofem::FEInterpolation.
Definition at line 109 of file feibspline.h.
|
inlineoverridevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 112 of file feibspline.h.
References oofem::_Line, oofem::_Point, oofem::_Square, oofem::_UnknownIntegrationDomain, and nsd.
|
inlineoverridevirtual |
Sets up a suitable integration rule for integrating over the requested boundary. The required polynomial order for the determinant of the jacobian is added automatically.
| order | Polynomial order of the integrand (should NOT including determinant of jacobian). |
| boundary | Boundary number. |
Reimplemented from oofem::FEInterpolation.
Definition at line 207 of file feibspline.h.
References OOFEM_ERROR, and oofem::FEInterpolation::order.
|
inlineoverridevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 123 of file feibspline.h.
References oofem::_Square, oofem::_UnknownIntegrationDomain, and nsd.
|
inline |
Definition at line 202 of file feibspline.h.
|
inlineoverridevirtual |
Returns the geometry type fo the interpolator.
Implements oofem::FEInterpolation.
Definition at line 108 of file feibspline.h.
|
inlineoverridevirtual |
Returns the integration domain of the interpolator.
Implements oofem::FEInterpolation.
Definition at line 97 of file feibspline.h.
References oofem::_Cube, oofem::_Line, oofem::_Square, oofem::_UnknownIntegrationDomain, and nsd.
|
inlineoverridevirtual |
Sets up a suitable integration rule for numerical integrating over volume. The required polynomial order for the determinant of the jacobian is added automatically.
| order | Polynomial order of integrand (should NOT including determinant of jacobian). |
Reimplemented from oofem::FEInterpolation.
Definition at line 205 of file feibspline.h.
References OOFEM_ERROR, and oofem::FEInterpolation::order.
|
overridevirtual |
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::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 482 of file feibspline.C.
References oofem::FloatArray::add(), degree, dersBasisFuns(), findSpan(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, nsd, numberOfControlPoints, OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().
|
inlineoverridevirtual |
Returns the knot multiplicity of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 190 of file feibspline.h.
|
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::FEInterpolation.
Reimplemented in oofem::TSplineInterpolation.
Definition at line 587 of file feibspline.C.
References oofem::IntArray::at(), degree, giveNumberOfKnotSpanBasisFunctions(), nsd, numberOfControlPoints, OOFEM_ERROR, and oofem::IntArray::resize().
|
inlineoverridevirtual |
Returns the knot values of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 191 of file feibspline.h.
Referenced by oofem::TSplineInterpolation::evaldNdx(), oofem::TSplineInterpolation::evalN(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), and oofem::TSplineInterpolation::local2global().
|
inlineoverridevirtual |
Returns the subdivision of patch parametric space
Reimplemented from oofem::FEInterpolation.
Definition at line 187 of file feibspline.h.
References oofem::FloatArray::data().
|
inlineprotected |
Returns the range of nonzero basis functions for given knot span and given degree.
Definition at line 263 of file feibspline.h.
|
inlineoverridevirtual |
Returns number of spatial dimensions.
Implements oofem::FEInterpolation.
Definition at line 140 of file feibspline.h.
References nsd.
|
inlinevirtual |
Definition at line 186 of file feibspline.h.
References numberOfControlPoints.
Referenced by oofem::BsplinePlaneStressElement::checkConsistency(), oofem::NURBSPlaneStressElement::checkConsistency(), and oofem::NURBSSpace3dElement::checkConsistency().
|
overridevirtual |
Returns the number of nonzero basis functions at individual knot span,
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::TSplineInterpolation.
Definition at line 625 of file feibspline.C.
Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), and giveKnotSpanBasisFuncMask().
|
inlineoverridevirtual |
Returns the number of knot spans of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 185 of file feibspline.h.
References numberOfKnotSpans.
|
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. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 195 of file feibspline.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
Returns true, if receiver is formulated on sub-patch basis.
Reimplemented from oofem::FEInterpolation.
Definition at line 203 of file feibspline.h.
|
overridevirtual |
Initializes receiver according to object description stored in input record.
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 56 of file feibspline.C.
References _IFT_BSplineInterpolation_degree, oofem::IntArray::at(), degree, oofem::IntArray::giveSize(), IPK_BSplineInterpolation_degree, IPK_BSplineInterpolation_knotMultiplicityU, IPK_BSplineInterpolation_knotMultiplicityV, IPK_BSplineInterpolation_knotMultiplicityW, IPK_BSplineInterpolation_knotVectorU, IPK_BSplineInterpolation_knotVectorV, IPK_BSplineInterpolation_knotVectorW, knotMultiplicity, knotValues, nsd, PM_UPDATE_PARAMETER, and PM_UPDATE_PARAMETER_AND_REPORT.
|
overridevirtual |
Evaluates global coordinates from given local ones.
| answer | Contains resulting global coordinates. |
| lcoords | Array containing (local) coordinates. |
| cellgeo | Underlying cell geometry. |
Implements oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.
Definition at line 400 of file feibspline.C.
References oofem::FloatArray::add(), basisFuns(), degree, findSpan(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, N, nsd, numberOfControlPoints, OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::NURBSInterpolation.
Definition at line 86 of file feibspline.C.
References oofem::ComponentInputException::ctElement, degree, IPK_BSplineInterpolation_knotMultiplicityU, IPK_BSplineInterpolation_knotMultiplicityV, IPK_BSplineInterpolation_knotMultiplicityW, IPK_BSplineInterpolation_knotVectorU, IPK_BSplineInterpolation_knotVectorV, IPK_BSplineInterpolation_knotVectorW, knotMultiplicity, knotValues, knotVector, nsd, numberOfControlPoints, numberOfKnotSpans, OOFEM_LOG_RELEVANT, and oofem::sum().
Referenced by oofem::IGATSplineElement::initializeFinish().
|
protected |
Degree in each direction.
Definition at line 68 of file feibspline.h.
Referenced by oofem::TSplineInterpolation::createLocalKnotVector(), evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), initializeFrom(), oofem::TSplineInterpolation::initializeFrom(), local2global(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), and postInitialize().
|
staticprotected |
Definition at line 84 of file feibspline.h.
Referenced by initializeFrom().
|
staticprotected |
Definition at line 88 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
Definition at line 89 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
Definition at line 90 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
Definition at line 85 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
Definition at line 86 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
Definition at line 87 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
protected |
Knot multiplicity [nsd].
Definition at line 72 of file feibspline.h.
Referenced by initializeFrom(), and postInitialize().
|
protected |
Knot values [nsd].
Definition at line 70 of file feibspline.h.
Referenced by oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), initializeFrom(), oofem::TSplineInterpolation::initializeFrom(), and postInitialize().
|
protected |
Knot vectors [nsd][knot_vector_size].
Definition at line 80 of file feibspline.h.
Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), local2global(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), and postInitialize().
|
protected |
Number of spatial directions.
Definition at line 66 of file feibspline.h.
Referenced by BSplineInterpolation(), oofem::TSplineInterpolation::createLocalKnotVector(), evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), giveBoundaryEdgeIntegrationDomain(), giveBoundaryIntegrationDomain(), giveBoundarySurfaceIntegrationDomain(), giveIntegrationDomain(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), giveNsd(), giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), initializeFrom(), oofem::TSplineInterpolation::initializeFrom(), local2global(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), oofem::NURBSInterpolation::NURBSInterpolation(), postInitialize(), and oofem::TSplineInterpolation::TSplineInterpolation().
|
protected |
numberOfControlPoints[nsd]
For TSpline this is filled by values corresponding to case when there are no T-junctions (i.e. TSpline = BSpline)
Definition at line 78 of file feibspline.h.
Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveKnotSpanBasisFuncMask(), giveNumberOfControlPoints(), local2global(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), and postInitialize().
|
protected |
Nonzero spans in each directions [nsd].
Definition at line 82 of file feibspline.h.
Referenced by giveNumberOfKnotSpans(), and postInitialize().