OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Interpolation for B-splines. More...
#include <feibspline.h>
Public Member Functions | |
BSplineInterpolation (int nsd) | |
virtual | ~BSplineInterpolation () |
virtual integrationDomain | giveIntegrationDomain () const |
Returns the integration domain of the interpolator. More... | |
virtual Element_Geometry_Type | giveGeometryType () const |
Returns the geometry type fo the interpolator. More... | |
virtual integrationDomain | giveBoundaryIntegrationDomain (int ib) const |
Returns boundary integration domain. More... | |
virtual integrationDomain | giveBoundarySurfaceIntegrationDomain (int isurf) const |
Returns boundary integration domain. More... | |
virtual integrationDomain | giveBoundaryEdgeIntegrationDomain (int iedge) const |
Returns boundary integration domain. More... | |
virtual int | giveNsd () |
Returns number of spatial dimensions. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual void | boundaryEdgeGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
virtual void | boundaryEdgeEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the basis functions on the requested boundary. More... | |
virtual double | boundaryEdgeGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation Jacobian on the requested boundary. More... | |
virtual void | boundaryEdgeLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Maps the local boundary coordinates to global. More... | |
virtual void | boundaryGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
virtual void | boundaryEvalN (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the basis functions on the requested boundary. More... | |
virtual double | boundaryEvalNormal (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the normal on the requested boundary. More... | |
virtual double | boundaryGiveTransformationJacobian (int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation Jacobian on the requested boundary. More... | |
virtual void | boundaryLocal2Global (FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Maps the local boundary coordinates to global. More... | |
virtual int | giveNumberOfKnotSpans (int dim) |
Returns the number of knot spans of the receiver. More... | |
virtual int | giveNumberOfControlPoints (int dim) |
virtual const double *const * | giveKnotVector () |
Returns the subdivision of patch parametric space. More... | |
virtual const IntArray * | giveKnotMultiplicity (int dim) |
Returns the knot multiplicity of the receiver. More... | |
virtual const FloatArray * | giveKnotValues (int dim) |
Returns the knot values of the receiver. More... | |
virtual void | evalN (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the array of interpolation functions (shape functions) at given point. More... | |
virtual double | evaldNdx (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | local2global (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates global coordinates from given local ones. More... | |
virtual int | global2local (FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates local coordinates from given global ones. More... | |
virtual void | giveJacobianMatrixAt (FloatMatrix &jacobianMatrix, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Gives the jacobian matrix at the local coordinates. More... | |
virtual int | giveKnotSpanBasisFuncMask (const IntArray &knotSpan, IntArray &mask) |
Returns indices (zero based) of nonzero basis functions for given knot span. More... | |
virtual int | giveNumberOfKnotSpanBasisFunctions (const IntArray &knotSpan) |
Returns the number of nonzero basis functions at individual knot span,. More... | |
virtual const char * | giveClassName () const |
virtual bool | hasSubPatchFormulation () |
Returns true, if receiver is formulated on sub-patch basis. More... | |
virtual IntegrationRule * | giveIntegrationRule (int order) |
Sets up a suitable integration rule for numerical integrating over volume. More... | |
virtual IntegrationRule * | giveBoundaryIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
virtual IntegrationRule * | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
Surface interpolation services | |
virtual void | boundarySurfaceEvalN (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the array of edge interpolation functions (shape functions) at given point. More... | |
virtual void | boundarySurfaceEvaldNdx (FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point. More... | |
virtual double | boundarySurfaceEvalNormal (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the normal out of the surface at given point. More... | |
virtual void | boundarySurfaceLocal2global (FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates edge global coordinates from given local ones. More... | |
virtual double | boundarySurfaceGiveTransformationJacobian (int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the edge jacobian of transformation between local and global coordinates. More... | |
virtual void | boundarySurfaceGiveNodes (IntArray &answer, int boundary) |
Gives the boundary nodes for requested boundary number. More... | |
Public Member Functions inherited from oofem::FEInterpolation | |
FEInterpolation (int o) | |
virtual | ~FEInterpolation () |
virtual int | giveNumberOfNodes () const |
Returns the number of geometric nodes of the receiver. More... | |
std::string | errorInfo (const char *func) const |
int | giveInterpolationOrder () |
Returns the interpolation order. More... | |
virtual void | evald2Ndx2 (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of second derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | evaldNdxi (FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point. More... | |
virtual void | giveLocalNodeCoords (FloatMatrix &answer) |
Returns a matrix containing the local coordinates for each node corresponding to the interpolation. More... | |
virtual double | giveTransformationJacobian (const FloatArray &lcoords, const FEICellGeometry &cellgeo) |
Evaluates the determinant of the transformation. More... | |
virtual IntegrationRule * | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
Sets up a suitable integration rule for integrating over the requested boundary. More... | |
virtual double | evalNXIntegral (int boundary, const FEICellGeometry &cellgeo) |
Computes the integral . More... | |
virtual int | giveNumberOfEdges () const |
Returns number of edges. More... | |
Protected Member Functions | |
void | basisFuns (FloatArray &N, int span, double u, int p, const double *U) |
Evaluates the nonvanishing basis functions of 1d BSpline (algorithm A2.2 from NURBS book) More... | |
void | dersBasisFuns (int n, double u, int span, int p, double *const U, FloatMatrix &ders) |
Computes nonzero basis functions and their derivatives at u. More... | |
int | findSpan (int n, int p, double u, const double *U) const |
Determines the knot span index (Algorithm A2.1 from the NURBS book) More... | |
void | giveNonzeroBasisFuncInterval (int span, int deg, int &s, int &e) |
Returns the range of nonzero basis functions for given knot span and given degree. More... | |
Protected Attributes | |
int | nsd |
Number of spatial directions. More... | |
int * | degree |
Degree in each direction. More... | |
FloatArray * | knotValues |
Knot values [nsd]. More... | |
IntArray * | knotMultiplicity |
Knot multiplicity [nsd]. More... | |
int * | numberOfControlPoints |
numberOfControlPoints[nsd] More... | |
double ** | knotVector |
Knot vectors [nsd][knot_vector_size]. More... | |
int * | numberOfKnotSpans |
Nonzero spans in each directions [nsd]. More... | |
Protected Attributes inherited from oofem::FEInterpolation | |
int | order |
Interpolation for B-splines.
Definition at line 60 of file feibspline.h.
|
inline |
Definition at line 82 of file feibspline.h.
|
virtual |
Definition at line 42 of file feibspline.C.
References degree, knotMultiplicity, knotValues, knotVector, nsd, numberOfControlPoints, and numberOfKnotSpans.
|
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 660 of file feibspline.C.
References N, and oofem::FloatArray::resize().
Referenced by oofem::TSplineInterpolation::basisFunction(), oofem::NURBSInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::local2global(), and local2global().
|
inlinevirtual |
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 135 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
Gives the boundary nodes for requested boundary number.
answer | Array to be filled with the boundary nodes. |
boundary | Boundary number. |
Implements oofem::FEInterpolation.
Definition at line 133 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 137 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 140 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 169 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 171 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 167 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 174 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 177 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 147 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 146 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 150 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
Gives the boundary nodes for requested boundary number.
answer | Array to be filled with the boundary nodes. |
boundary | Boundary number. |
Implements oofem::FEInterpolation.
Definition at line 159 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 156 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
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 153 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
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 692 of file feibspline.C.
References oofem::FloatMatrix::resize().
Referenced by oofem::TSplineInterpolation::dersBasisFunction(), oofem::NURBSInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), and giveJacobianMatrixAt().
|
virtual |
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::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 228 of file feibspline.C.
References oofem::FloatArray::at(), 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().
|
virtual |
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::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 181 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 792 of file feibspline.C.
Referenced by oofem::TSplineInterpolation::basisFunction(), oofem::TSplineInterpolation::dersBasisFunction(), oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), and local2global().
|
inlinevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 120 of file feibspline.h.
References oofem::_Line, and oofem::_UnknownIntegrationDomain.
|
inlinevirtual |
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 208 of file feibspline.h.
References N, and OOFEM_ERROR.
|
inlinevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 100 of file feibspline.h.
References oofem::_Line, oofem::_Point, oofem::_Square, and oofem::_UnknownIntegrationDomain.
|
inlinevirtual |
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 205 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
Returns boundary integration domain.
Implements oofem::FEInterpolation.
Definition at line 111 of file feibspline.h.
References oofem::_Square, and oofem::_UnknownIntegrationDomain.
|
inlinevirtual |
Reimplemented in oofem::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 199 of file feibspline.h.
|
inlinevirtual |
Returns the geometry type fo the interpolator.
Implements oofem::FEInterpolation.
Definition at line 98 of file feibspline.h.
|
inlinevirtual |
Returns the integration domain of the interpolator.
Implements oofem::FEInterpolation.
Definition at line 87 of file feibspline.h.
References oofem::_Cube, oofem::_Line, oofem::_Square, and oofem::_UnknownIntegrationDomain.
|
inlinevirtual |
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 202 of file feibspline.h.
References OOFEM_ERROR.
|
virtual |
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::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 495 of file feibspline.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), degree, dersBasisFuns(), findSpan(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, nsd, numberOfControlPoints, OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().
|
inlinevirtual |
Returns the knot multiplicity of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 186 of file feibspline.h.
|
virtual |
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 602 of file feibspline.C.
References oofem::IntArray::at(), degree, giveNumberOfKnotSpanBasisFunctions(), nsd, numberOfControlPoints, OOFEM_ERROR, and oofem::IntArray::resize().
|
inlinevirtual |
Returns the knot values of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 187 of file feibspline.h.
Referenced by oofem::TSplineInterpolation::evaldNdx(), oofem::TSplineInterpolation::evalN(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), and oofem::TSplineInterpolation::local2global().
|
inlinevirtual |
Returns the subdivision of patch parametric space.
Reimplemented from oofem::FEInterpolation.
Definition at line 183 of file feibspline.h.
|
inlineprotected |
Returns the range of nonzero basis functions for given knot span and given degree.
Definition at line 263 of file feibspline.h.
|
inlinevirtual |
Returns number of spatial dimensions.
Implements oofem::FEInterpolation.
Definition at line 130 of file feibspline.h.
|
inlinevirtual |
Definition at line 182 of file feibspline.h.
Referenced by oofem::BsplinePlaneStressElement::checkConsistency(), oofem::NURBSPlaneStressElement::checkConsistency(), and oofem::NURBSSpace3dElement::checkConsistency().
|
virtual |
Returns the number of nonzero basis functions at individual knot span,.
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::TSplineInterpolation.
Definition at line 641 of file feibspline.C.
Referenced by oofem::NURBSInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), evalN(), and giveKnotSpanBasisFuncMask().
|
inlinevirtual |
Returns the number of knot spans of the receiver.
Reimplemented from oofem::FEInterpolation.
Definition at line 181 of file feibspline.h.
|
inlinevirtual |
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::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 191 of file feibspline.h.
References OOFEM_ERROR.
|
inlinevirtual |
Returns true, if receiver is formulated on sub-patch basis.
Reimplemented from oofem::FEInterpolation.
Definition at line 200 of file feibspline.h.
|
virtual |
Initializes receiver according to object description stored in input record.
Reimplemented from oofem::FEInterpolation.
Reimplemented in oofem::TSplineInterpolation.
Definition at line 59 of file feibspline.C.
References _IFT_BSplineInterpolation_degree, _IFT_BSplineInterpolation_knotMultiplicityU, _IFT_BSplineInterpolation_knotMultiplicityV, _IFT_BSplineInterpolation_knotMultiplicityW, _IFT_BSplineInterpolation_knotVectorU, _IFT_BSplineInterpolation_knotVectorV, _IFT_BSplineInterpolation_knotVectorW, oofem::IntArray::at(), oofem::FloatArray::at(), degree, oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_BAD_FORMAT, oofem::IRRT_OK, knotMultiplicity, knotValues, knotVector, nsd, numberOfControlPoints, numberOfKnotSpans, OOFEM_LOG_RELEVANT, OOFEM_WARNING, and oofem::IntArray::resize().
Referenced by oofem::TSplineInterpolation::initializeFrom().
|
virtual |
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::TSplineInterpolation, and oofem::NURBSInterpolation.
Definition at line 411 of file feibspline.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), basisFuns(), degree, findSpan(), oofem::FEICellGeometry::giveVertexCoordinates(), oofem::FEIIGAElementGeometryWrapper::knotSpan, knotVector, N, nsd, numberOfControlPoints, OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
protected |
Degree in each direction.
Definition at line 66 of file feibspline.h.
Referenced by oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::initializeFrom(), initializeFrom(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), local2global(), and ~BSplineInterpolation().
|
protected |
Knot multiplicity [nsd].
Definition at line 70 of file feibspline.h.
Referenced by initializeFrom(), and ~BSplineInterpolation().
|
protected |
Knot values [nsd].
Definition at line 68 of file feibspline.h.
Referenced by oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::initializeFrom(), initializeFrom(), and ~BSplineInterpolation().
|
protected |
Knot vectors [nsd][knot_vector_size].
Definition at line 78 of file feibspline.h.
Referenced by oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), initializeFrom(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), local2global(), and ~BSplineInterpolation().
|
protected |
Number of spatial directions.
Definition at line 64 of file feibspline.h.
Referenced by oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveKnotSpanBasisFuncMask(), giveKnotSpanBasisFuncMask(), oofem::TSplineInterpolation::giveNumberOfKnotSpanBasisFunctions(), giveNumberOfKnotSpanBasisFunctions(), oofem::TSplineInterpolation::initializeFrom(), initializeFrom(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), local2global(), ~BSplineInterpolation(), 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 76 of file feibspline.h.
Referenced by oofem::NURBSInterpolation::evaldNdx(), oofem::TSplineInterpolation::evaldNdx(), evaldNdx(), oofem::NURBSInterpolation::evalN(), oofem::TSplineInterpolation::evalN(), evalN(), oofem::NURBSInterpolation::giveJacobianMatrixAt(), oofem::TSplineInterpolation::giveJacobianMatrixAt(), giveJacobianMatrixAt(), giveKnotSpanBasisFuncMask(), initializeFrom(), oofem::NURBSInterpolation::local2global(), oofem::TSplineInterpolation::local2global(), local2global(), ~BSplineInterpolation(), and oofem::TSplineInterpolation::~TSplineInterpolation().
|
protected |
Nonzero spans in each directions [nsd].
Definition at line 80 of file feibspline.h.
Referenced by initializeFrom(), and ~BSplineInterpolation().