OOFEM 3.0
Loading...
Searching...
No Matches
oofem::BSplineInterpolation Class Reference

#include <feibspline.h>

Inheritance diagram for oofem::BSplineInterpolation:
Collaboration diagram for oofem::BSplineInterpolation:

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< IntegrationRulegiveBoundarySurfaceIntegrationRule (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 FloatArraygiveKnotVector () const override
const IntArraygiveKnotMultiplicity (int dim) const override
const FloatArraygiveKnotValues (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< IntegrationRulegiveIntegrationRule (int order, const Element_Geometry_Type) const override
std::unique_ptr< IntegrationRulegiveBoundaryIntegrationRule (int order, int boundary, Element_Geometry_Type) const override
std::unique_ptr< IntegrationRulegiveBoundaryEdgeIntegrationRule (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

Detailed Description

Interpolation for B-splines.

Definition at line 62 of file feibspline.h.

Constructor & Destructor Documentation

◆ BSplineInterpolation()

oofem::BSplineInterpolation::BSplineInterpolation ( int nsd)
inline

Member Function Documentation

◆ basisFuns()

void oofem::BSplineInterpolation::basisFuns ( FloatArray & N,
int span,
double u,
int p,
const FloatArray & U ) const
protected

Evaluates the nonvanishing basis functions of 1d BSpline (algorithm A2.2 from NURBS book)

Parameters
spanKnot span index (zero based).
uValue at which to evaluate.
pDegree.
UKnot vector.
NComputed p+1 nonvanishing functions (N_{span-p,p}-N_{span,p})
Warning
Parameter u and span must be in a valid range.

Definition at line 644 of file feibspline.C.

References N.

Referenced by evalN(), oofem::NURBSInterpolation::evalN(), local2global(), and oofem::NURBSInterpolation::local2global().

◆ boundaryEdgeEvalN()

void oofem::BSplineInterpolation::boundaryEdgeEvalN ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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 146 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEdgeEvalNormal()

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

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 152 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEdgeGiveNodes()

IntArray oofem::BSplineInterpolation::boundaryEdgeGiveNodes ( int boundary,
const Element_Geometry_Type ,
bool includeHierarchical = false ) const
inlineoverridevirtual

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 144 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEdgeGiveTransformationJacobian()

double oofem::BSplineInterpolation::boundaryEdgeGiveTransformationJacobian ( int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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 148 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEdgeLocal2Global()

void oofem::BSplineInterpolation::boundaryEdgeLocal2Global ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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

Implements oofem::FEInterpolation.

Definition at line 150 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEvalN()

void oofem::BSplineInterpolation::boundaryEvalN ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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 175 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryEvalNormal()

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

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 177 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryGiveNodes()

IntArray oofem::BSplineInterpolation::boundaryGiveNodes ( int boundary,
const Element_Geometry_Type  ) const
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.

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

Implements oofem::FEInterpolation.

Definition at line 173 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryGiveTransformationJacobian()

double oofem::BSplineInterpolation::boundaryGiveTransformationJacobian ( int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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 179 of file feibspline.h.

References OOFEM_ERROR.

◆ boundaryLocal2Global()

void oofem::BSplineInterpolation::boundaryLocal2Global ( FloatArray & answer,
int boundary,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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.

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

Implements oofem::FEInterpolation.

Definition at line 181 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceEvaldNdx()

void oofem::BSplineInterpolation::boundarySurfaceEvaldNdx ( FloatMatrix & answer,
int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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).

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 160 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceEvalN()

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

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 158 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceEvalNormal()

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

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 162 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceGiveNodes()

IntArray oofem::BSplineInterpolation::boundarySurfaceGiveNodes ( int boundary,
const Element_Geometry_Type ,
bool includeHierarchical = false ) const
inlineoverridevirtual

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 168 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceGiveTransformationJacobian()

double oofem::BSplineInterpolation::boundarySurfaceGiveTransformationJacobian ( int isurf,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
inlineoverridevirtual

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 166 of file feibspline.h.

References OOFEM_ERROR.

◆ boundarySurfaceLocal2global()

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

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 164 of file feibspline.h.

References OOFEM_ERROR.

◆ dersBasisFuns()

void oofem::BSplineInterpolation::dersBasisFuns ( int n,
double u,
int span,
int p,
const FloatArray & U,
FloatMatrix & ders ) const
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*}

Parameters
nDegree of the derivation.
uParametric value.
spanKnot span index (zero based).
pDegree.
UKnot vector.
dersMatrix containing the derivatives of the basis functions.
Warning
Parameters n, u and span must be in a valid range.

Definition at line 675 of file feibspline.C.

References oofem::FloatMatrix::resize().

Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), giveJacobianMatrixAt(), and oofem::NURBSInterpolation::giveJacobianMatrixAt().

◆ evaldNdx()

double oofem::BSplineInterpolation::evaldNdx ( FloatMatrix & answer,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
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)

Parameters
answerContains resulting matrix of derivatives, the member at i,j position contains value of dNi/dxj.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying cell geometry.
Returns
Determinant of the Jacobian.

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

◆ evalN()

void oofem::BSplineInterpolation::evalN ( FloatArray & answer,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

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

Parameters
answerContains resulting array of evaluated interpolation functions.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying 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().

◆ findSpan()

int oofem::BSplineInterpolation::findSpan ( int n,
int p,
double u,
const FloatArray & U ) const
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.

Parameters
nNumber of control points - 1.
uParametric value.
pDegree.
UKnot vector.
Returns
Span index at u (zero based).
Warning
Parameter u must be in a valid range.

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

◆ giveBoundaryEdgeIntegrationDomain()

integrationDomain oofem::BSplineInterpolation::giveBoundaryEdgeIntegrationDomain ( int boundary,
const Element_Geometry_Type  ) const
inlineoverridevirtual

Returns boundary integration domain.

Implements oofem::FEInterpolation.

Definition at line 132 of file feibspline.h.

References oofem::_Line, oofem::_UnknownIntegrationDomain, and nsd.

◆ giveBoundaryEdgeIntegrationRule()

std::unique_ptr< IntegrationRule > oofem::BSplineInterpolation::giveBoundaryEdgeIntegrationRule ( int order,
int boundary,
Element_Geometry_Type egt ) const
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.

Parameters
orderPolynomial order of the integrand (should NOT including determinant of jacobian).
boundaryBoundary number.

Reimplemented from oofem::FEInterpolation.

Definition at line 209 of file feibspline.h.

References OOFEM_ERROR, and oofem::FEInterpolation::order.

◆ giveBoundaryGeometryType()

const Element_Geometry_Type oofem::BSplineInterpolation::giveBoundaryGeometryType ( int boundary) const
inlineoverridevirtual

Returns boundary geometry type

Implements oofem::FEInterpolation.

Definition at line 109 of file feibspline.h.

◆ giveBoundaryIntegrationDomain()

integrationDomain oofem::BSplineInterpolation::giveBoundaryIntegrationDomain ( int boundary,
const Element_Geometry_Type  ) const
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.

◆ giveBoundaryIntegrationRule()

std::unique_ptr< IntegrationRule > oofem::BSplineInterpolation::giveBoundaryIntegrationRule ( int order,
int boundary,
Element_Geometry_Type egt ) const
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.

Parameters
orderPolynomial order of the integrand (should NOT including determinant of jacobian).
boundaryBoundary number.

Reimplemented from oofem::FEInterpolation.

Definition at line 207 of file feibspline.h.

References OOFEM_ERROR, and oofem::FEInterpolation::order.

◆ giveBoundarySurfaceIntegrationDomain()

integrationDomain oofem::BSplineInterpolation::giveBoundarySurfaceIntegrationDomain ( int boundary,
const Element_Geometry_Type  ) const
inlineoverridevirtual

Returns boundary integration domain.

Implements oofem::FEInterpolation.

Definition at line 123 of file feibspline.h.

References oofem::_Square, oofem::_UnknownIntegrationDomain, and nsd.

◆ giveClassName()

const char * oofem::BSplineInterpolation::giveClassName ( ) const
inline

Definition at line 202 of file feibspline.h.

◆ giveGeometryType()

const Element_Geometry_Type oofem::BSplineInterpolation::giveGeometryType ( ) const
inlineoverridevirtual

Returns the geometry type fo the interpolator.

Implements oofem::FEInterpolation.

Definition at line 108 of file feibspline.h.

◆ giveIntegrationDomain()

integrationDomain oofem::BSplineInterpolation::giveIntegrationDomain ( const Element_Geometry_Type ) const
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.

◆ giveIntegrationRule()

std::unique_ptr< IntegrationRule > oofem::BSplineInterpolation::giveIntegrationRule ( int order,
const Element_Geometry_Type egt ) const
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.

Parameters
orderPolynomial 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.

◆ giveJacobianMatrixAt()

void oofem::BSplineInterpolation::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.

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

◆ giveKnotMultiplicity()

const IntArray * oofem::BSplineInterpolation::giveKnotMultiplicity ( int dim) const
inlineoverridevirtual

Returns the knot multiplicity of the receiver.

Reimplemented from oofem::FEInterpolation.

Definition at line 190 of file feibspline.h.

◆ giveKnotSpanBasisFuncMask()

int oofem::BSplineInterpolation::giveKnotSpanBasisFuncMask ( const IntArray & knotSpan,
IntArray & mask ) const
overridevirtual

Returns indices (zero based) of nonzero basis functions for given knot span. The knot span identifies the sub-region of the finite element.

Returns
Nonzero if mask is provided, zero otherwise meaning that all basis functions are generally nonzero.

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

◆ giveKnotValues()

const FloatArray * oofem::BSplineInterpolation::giveKnotValues ( int dim) const
inlineoverridevirtual

◆ giveKnotVector()

const FloatArray * oofem::BSplineInterpolation::giveKnotVector ( ) const
inlineoverridevirtual

Returns the subdivision of patch parametric space

Reimplemented from oofem::FEInterpolation.

Definition at line 187 of file feibspline.h.

References oofem::FloatArray::data().

◆ giveNonzeroBasisFuncInterval()

void oofem::BSplineInterpolation::giveNonzeroBasisFuncInterval ( int span,
int deg,
int & s,
int & e ) const
inlineprotected

Returns the range of nonzero basis functions for given knot span and given degree.

Definition at line 263 of file feibspline.h.

◆ giveNsd()

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

Returns number of spatial dimensions.

Implements oofem::FEInterpolation.

Definition at line 140 of file feibspline.h.

References nsd.

◆ giveNumberOfControlPoints()

virtual int oofem::BSplineInterpolation::giveNumberOfControlPoints ( int dim) const
inlinevirtual

◆ giveNumberOfKnotSpanBasisFunctions()

int oofem::BSplineInterpolation::giveNumberOfKnotSpanBasisFunctions ( const IntArray & knotSpan) const
overridevirtual

Returns the number of nonzero basis functions at individual knot span,

Returns
Zero in case of all basis functions generaedgeEvaldNdslly nonzero, answer otherwise.
Todo
This loop seems meaningless. It just returns degree[nsd-1]+1 in the end ?

Reimplemented from oofem::FEInterpolation.

Reimplemented in oofem::TSplineInterpolation.

Definition at line 625 of file feibspline.C.

References degree, and nsd.

Referenced by evaldNdx(), oofem::NURBSInterpolation::evaldNdx(), evalN(), oofem::NURBSInterpolation::evalN(), and giveKnotSpanBasisFuncMask().

◆ giveNumberOfKnotSpans()

int oofem::BSplineInterpolation::giveNumberOfKnotSpans ( int dim) const
inlineoverridevirtual

Returns the number of knot spans of the receiver.

Reimplemented from oofem::FEInterpolation.

Definition at line 185 of file feibspline.h.

References numberOfKnotSpans.

◆ global2local()

int oofem::BSplineInterpolation::global2local ( FloatArray & answer,
const FloatArray & gcoords,
const FEICellGeometry & cellgeo ) const
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.

Parameters
answerContains evaluated local coordinates.
gcoordsArray containing global coordinates.
cellgeoUnderlying cell geometry.
Returns
Nonzero is returned if point is within the element geometry, zero otherwise.

Implements oofem::FEInterpolation.

Reimplemented in oofem::NURBSInterpolation, and oofem::TSplineInterpolation.

Definition at line 195 of file feibspline.h.

References OOFEM_ERROR.

◆ hasSubPatchFormulation()

bool oofem::BSplineInterpolation::hasSubPatchFormulation ( ) const
inlineoverridevirtual

Returns true, if receiver is formulated on sub-patch basis.

Reimplemented from oofem::FEInterpolation.

Definition at line 203 of file feibspline.h.

◆ initializeFrom()

◆ local2global()

void oofem::BSplineInterpolation::local2global ( FloatArray & answer,
const FloatArray & lcoords,
const FEICellGeometry & cellgeo ) const
overridevirtual

Evaluates global coordinates from given local ones.

Parameters
answerContains resulting global coordinates.
lcoordsArray containing (local) coordinates.
cellgeoUnderlying 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().

◆ postInitialize()

Member Data Documentation

◆ degree

◆ IPK_BSplineInterpolation_degree

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_degree
staticprotected

Definition at line 84 of file feibspline.h.

Referenced by initializeFrom().

◆ IPK_BSplineInterpolation_knotMultiplicityU

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotMultiplicityU
staticprotected

Definition at line 88 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ IPK_BSplineInterpolation_knotMultiplicityV

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotMultiplicityV
staticprotected

Definition at line 89 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ IPK_BSplineInterpolation_knotMultiplicityW

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotMultiplicityW
staticprotected

Definition at line 90 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ IPK_BSplineInterpolation_knotVectorU

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotVectorU
staticprotected

Definition at line 85 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ IPK_BSplineInterpolation_knotVectorV

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotVectorV
staticprotected

Definition at line 86 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ IPK_BSplineInterpolation_knotVectorW

ParamKey oofem::BSplineInterpolation::IPK_BSplineInterpolation_knotVectorW
staticprotected

Definition at line 87 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ knotMultiplicity

std::array<IntArray, 3> oofem::BSplineInterpolation::knotMultiplicity
protected

Knot multiplicity [nsd].

Definition at line 72 of file feibspline.h.

Referenced by initializeFrom(), and postInitialize().

◆ knotValues

◆ knotVector

◆ nsd

◆ numberOfControlPoints

◆ numberOfKnotSpans

std::array<int, 3> oofem::BSplineInterpolation::numberOfKnotSpans
protected

Nonzero spans in each directions [nsd].

Definition at line 82 of file feibspline.h.

Referenced by giveNumberOfKnotSpans(), and postInitialize().


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