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

#include <feitspline.h>

Inheritance diagram for oofem::TSplineInterpolation:
Collaboration diagram for oofem::TSplineInterpolation:

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 FloatArraygiveKnotVector () const override
const IntArraygiveKnotMultiplicity (int dim) const override
const FloatArraygiveKnotValues (int dim) 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
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 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

Detailed Description

Interpolation for T-splines.

Definition at line 62 of file feitspline.h.

Constructor & Destructor Documentation

◆ TSplineInterpolation()

oofem::TSplineInterpolation::TSplineInterpolation ( int nsd)
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.

Member Function Documentation

◆ basisFunction()

double oofem::TSplineInterpolation::basisFunction ( double u,
int p,
const FloatArray & U,
const IntArray & I ) const
protected

Evaluates the middle basis function on local knot vector at u.

Parameters
uValue at which to evaluate.
pDegree.
UGlobal knot values.
ILocal index knot vector.
Returns
N Computed middle basis function.
Warning
Parameter u must be in a valid range.

Definition at line 607 of file feitspline.C.

References createLocalKnotVector(), and N.

Referenced by evalN(), and local2global().

◆ createLocalKnotVector()

void oofem::TSplineInterpolation::createLocalKnotVector ( FloatArray & answer,
int p,
const FloatArray & U,
const IntArray & I,
int & prepend,
int & append ) const
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.

Parameters
pDegree.
UGlobal knot values.
ILocal index knot vector
prependNumber of prepended entries
appendNumber 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().

◆ dersBasisFunction()

void oofem::TSplineInterpolation::dersBasisFunction ( int n,
double u,
int p,
const FloatArray & U,
const IntArray & I,
FloatArray & ders ) const
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*}

Parameters
nDegree of the derivation.
uParametric value.
pDegree.
UGlobal knot values.
ILocal index knot vector.
dersVector containing the derivatives of the middle basis function.
Warning
Parameters n and u must be in a valid range.

Definition at line 624 of file feitspline.C.

References createLocalKnotVector(), and oofem::FloatArray::resize().

Referenced by evaldNdx(), and giveJacobianMatrixAt().

◆ evaldNdx()

double oofem::TSplineInterpolation::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.

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

◆ evalN()

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

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.

◆ giveClassName()

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

Definition at line 90 of file feitspline.h.

◆ giveJacobianMatrixAt()

◆ giveKnotSpanBasisFuncMask() [1/2]

int oofem::TSplineInterpolation::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::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().

◆ giveKnotSpanBasisFuncMask() [2/2]

int oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask ( const IntArray & startKnotSpan,
const IntArray & endKnotSpan,
IntArray & mask ) const
protected

◆ giveNumberOfKnotSpanBasisFunctions() [1/2]

int oofem::TSplineInterpolation::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::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.

◆ giveNumberOfKnotSpanBasisFunctions() [2/2]

int oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions ( const IntArray & startKnotSpan,
const IntArray & endKnotSpan ) const
protected

◆ global2local()

int oofem::TSplineInterpolation::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.

Reimplemented from oofem::BSplineInterpolation.

Definition at line 82 of file feitspline.h.

References OOFEM_ERROR.

◆ initializeFrom()

◆ local2global()

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

◆ setNumberOfControlPoints()

void oofem::TSplineInterpolation::setNumberOfControlPoints ( int num)
inline

Definition at line 78 of file feitspline.h.

Referenced by oofem::IGATSplineElement::initializeFinish().

Member Data Documentation

◆ localIndexKnotVector

std::vector< std::array<IntArray, 3> > oofem::TSplineInterpolation::localIndexKnotVector
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().

◆ totalNumberOfControlPoints

int oofem::TSplineInterpolation::totalNumberOfControlPoints = 0
protected

◆ weights

FloatArray oofem::TSplineInterpolation::weights
protected

Definition at line 68 of file feitspline.h.

Referenced by evaldNdx(), evalN(), giveJacobianMatrixAt(), initializeFrom(), and local2global().


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

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