DUNE-FUNCTIONS (unstable)

Pre-basis for B-spline basis. More...

#include <dune/functions/functionspacebases/bsplinebasis.hh>

Public Types

using GridView = GV
 The grid view that the FE space is defined on.
 

Public Member Functions

 BSplinePreBasis (const GridView &gridView, const std::vector< double > &knotVector, unsigned int order, bool makeOpen=true)
 Construct a B-spline basis for a given grid view and set of knot vectors. More...
 
 BSplinePreBasis (const GridView &gridView, const FieldVector< double, dim > &lowerLeft, const FieldVector< double, dim > &upperRight, const std::array< unsigned int, dim > &elements, unsigned int order, bool makeOpen=true)
 Construct a B-spline basis for a given grid view with uniform knot vectors. More...
 
void initializeIndices ()
 Initialize the global indices.
 
const GridViewgridView () const
 Obtain the grid view that the basis is defined on.
 
void update (const GridView &gv)
 Update the stored grid view, to be called if the grid has changed.
 
Node makeNode () const
 Create tree node.
 
size_type maxNodeSize () const
 Get the maximal number of DOFs associated to node for any element.
 
template<typename It >
It indices (const Node &node, It it) const
 Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
 
unsigned int dimension () const
 Total number of B-spline basis functions.
 
unsigned int size (size_t d) const
 Number of shape functions in one direction.
 
void evaluateFunction (const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > &currentKnotSpan) const
 Evaluate all B-spline basis functions at a given point.
 
void evaluateJacobian (const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldMatrix< R, 1, dim > > &out, const std::array< unsigned, dim > &currentKnotSpan) const
 Evaluate Jacobian of all B-spline basis functions. More...
 
template<size_type k>
void evaluate (const typename std::array< int, k > &directions, const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > &currentKnotSpan) const
 Evaluate Derivatives of all B-spline basis functions.
 
size_type size (const SizePrefix &prefix) const
 Return number of possible values for next position in multi index.
 
size_type size () const
 Get the total dimension of the space spanned by this basis.
 
auto containerDescriptor () const
 Return a flat container-descriptor.
 

Static Public Member Functions

static std::array< unsigned int, dim > getIJK (typename GridView::IndexSet::IndexType idx, std::array< unsigned int, dim > elements)
 Compute integer element coordinates from the element index. More...
 
static void evaluateFunction (const typename GV::ctype &in, std::vector< R > &out, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
 Evaluate all one-dimensional B-spline functions for a given coordinate direction. More...
 
static void evaluateFunctionFull (const typename GV::ctype &in, DynamicMatrix< R > &out, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
 Evaluate all one-dimensional B-spline functions for a given coordinate direction. More...
 
static void evaluateAll (const typename GV::ctype &in, std::vector< R > &out, bool evaluateJacobian, std::vector< R > &outJac, bool evaluateHessian, std::vector< R > &outHess, const std::vector< R > &knotVector, unsigned int order, unsigned int currentKnotSpan)
 Evaluate the second derivatives of all one-dimensional B-spline functions for a given coordinate direction. More...
 

Public Attributes

std::array< unsigned int, dim > order_
 Order of the B-spline for each space dimension.
 
std::array< std::vector< double >, dim > knotVectors_
 The knot vectors, one for each space dimension.
 
std::array< unsigned, dim > elements_
 Number of grid elements in the different coordinate directions.
 

Static Public Attributes

static constexpr size_type maxMultiIndexSize
 Maximal length of global multi-indices.
 
static constexpr size_type minMultiIndexSize
 Minimal length of global multi-indices.
 
static constexpr size_type multiIndexBufferSize
 Size required temporarily when constructing global multi-indices.
 

Detailed Description

template<typename GV>
class Dune::Functions::BSplinePreBasis< GV >

Pre-basis for B-spline basis.

Template Parameters
GVThe GridView that the space is defined on

The BSplinePreBasis can be used to embed a BSplineBasis in a larger basis for the construction of product spaces.

Constructor & Destructor Documentation

◆ BSplinePreBasis() [1/2]

template<typename GV >
Dune::Functions::BSplinePreBasis< GV >::BSplinePreBasis ( const GridView gridView,
const std::vector< double > &  knotVector,
unsigned int  order,
bool  makeOpen = true 
)
inline

Construct a B-spline basis for a given grid view and set of knot vectors.

The grid must match the knot vectors, i.e.:

  • The grid must be structured and Cartesian, and have cube elements only
  • The number of elements in each direction must match the number of knot spans in that direction
  • In fact, the element spacing in any direction must match the knot spacing in that direction (disregarding knot multiplicities)
  • When ordering the grid elements according to their indices, the resulting order must be lexicographical, with the x-index increasing fastest.

Unfortunately, not all of these conditions can be checked for automatically.

Parameters
knotVectorA single knot vector, which will be used for all coordinate directions
orderB-spline order, will be used for all coordinate directions
makeOpenIf this is true, then knots are prepended and appended to the knot vector to make the knot vector 'open'. i.e., start and end with 'order+1' identical knots. Basis functions from such knot vectors are interpolatory at the end of the parameter interval.
Todo:
maybe test whether the knot vector is already open?

References Dune::Functions::BSplinePreBasis< GV >::elements_, Dune::Functions::BSplinePreBasis< GV >::knotVectors_, and Dune::Functions::BSplinePreBasis< GV >::order_.

◆ BSplinePreBasis() [2/2]

template<typename GV >
Dune::Functions::BSplinePreBasis< GV >::BSplinePreBasis ( const GridView gridView,
const FieldVector< double, dim > &  lowerLeft,
const FieldVector< double, dim > &  upperRight,
const std::array< unsigned int, dim > &  elements,
unsigned int  order,
bool  makeOpen = true 
)
inline

Construct a B-spline basis for a given grid view with uniform knot vectors.

The grid must match the knot vectors, i.e.:

  • The grid must be structured and Cartesian, and have cube elements only
  • Bounding box and number of elements of the grid must match the corresponding arguments given to this constructor.
  • The element spacing must be uniform
  • When ordering the grid elements according to their indices, the resulting order must be lexicographical, with the x-index increasing fastest.

Unfortunately, not all of these conditions can be checked for automatically.

Parameters
gridViewThe grid we are defining the basis on
lowerLeftLower left corner of the structured grid
upperRightUpper right corner of the structured grid
elementsNumber of elements in each coordinate direction
orderB-spline order, will be used for all coordinate directions
makeOpenIf this is true, then knots are prepended and appended to the knot vector to make the knot vector 'open'. i.e., start and end with 'order+1' identical knots. Basis functions from such knot vectors are interpolatory at the end of the parameter interval.
Todo:
maybe test whether the knot vector is already open?

References Dune::Functions::BSplinePreBasis< GV >::elements_, Dune::Functions::BSplinePreBasis< GV >::knotVectors_, and Dune::Functions::BSplinePreBasis< GV >::order_.

Member Function Documentation

◆ evaluateAll()

template<typename GV >
static void Dune::Functions::BSplinePreBasis< GV >::evaluateAll ( const typename GV::ctype &  in,
std::vector< R > &  out,
bool  evaluateJacobian,
std::vector< R > &  outJac,
bool  evaluateHessian,
std::vector< R > &  outHess,
const std::vector< R > &  knotVector,
unsigned int  order,
unsigned int  currentKnotSpan 
)
inlinestatic

Evaluate the second derivatives of all one-dimensional B-spline functions for a given coordinate direction.

Parameters
inScalar(!) coordinate where to evaluate the functions
enableEvaluationsswitches calculation of desired derivatives on
[out]outVector containing the values of all B-spline derivatives at 'in'
[out]outJacVector containing the first derivatives of all B-spline derivatives at 'in' (only if calculation was switched on by enableEvaluations)
[out]outHessVector containing the second derivatives of all B-spline derivatives at 'in' (only if calculation was switched on by enableEvaluations)

References Dune::Functions::BSplinePreBasis< GV >::evaluateFunctionFull(), and Dune::Functions::BSplinePreBasis< GV >::evaluateJacobian().

Referenced by Dune::Functions::BSplinePreBasis< GV >::evaluate().

◆ evaluateFunction()

template<typename GV >
static void Dune::Functions::BSplinePreBasis< GV >::evaluateFunction ( const typename GV::ctype &  in,
std::vector< R > &  out,
const std::vector< R > &  knotVector,
unsigned int  order,
unsigned int  currentKnotSpan 
)
inlinestatic

Evaluate all one-dimensional B-spline functions for a given coordinate direction.

This implementations was based on the explanations in the book of Cottrell, Hughes, Bazilevs, "Isogeometric Analysis"

Parameters
inScalar(!) coordinate where to evaluate the functions
[out]outVector containing the values of all B-spline functions at 'in'
Todo:
We only hand out function values for those basis functions whose support overlaps the current knot span. However, in the preceding loop we still computed all values_. This won't scale.

◆ evaluateFunctionFull()

template<typename GV >
static void Dune::Functions::BSplinePreBasis< GV >::evaluateFunctionFull ( const typename GV::ctype &  in,
DynamicMatrix< R > &  out,
const std::vector< R > &  knotVector,
unsigned int  order,
unsigned int  currentKnotSpan 
)
inlinestatic

Evaluate all one-dimensional B-spline functions for a given coordinate direction.

This implementations was based on the explanations in the book of Cottrell, Hughes, Bazilevs, "Isogeometric Analysis"

Todo:
This method is a hack! I computes the derivatives of ALL B-splines, even the ones that are zero on the current knot span. I need it as an intermediate step to get the derivatives working. It will/must be removed as soon as possible.
Parameters
inScalar(!) coordinate where to evaluate the functions
[out]outVector containing the values of all B-spline functions at 'in'

Referenced by Dune::Functions::BSplinePreBasis< GV >::evaluateAll(), and Dune::Functions::BSplinePreBasis< GV >::evaluateJacobian().

◆ evaluateJacobian()

template<typename GV >
void Dune::Functions::BSplinePreBasis< GV >::evaluateJacobian ( const FieldVector< typename GV::ctype, dim > &  in,
std::vector< FieldMatrix< R, 1, dim > > &  out,
const std::array< unsigned, dim > &  currentKnotSpan 
) const
inline

Evaluate Jacobian of all B-spline basis functions.

In theory this is easy: just look up the formula in a B-spline text of your choice. The challenge is compute only the values needed for the current knot span.

References Dune::Functions::BSplinePreBasis< GV >::evaluateFunctionFull(), Dune::Functions::BSplinePreBasis< GV >::knotVectors_, Dune::Functions::BSplinePreBasis< GV >::order_, and Dune::Functions::LeafPreBasisMixin< BSplinePreBasis< GV > >::size().

Referenced by Dune::Functions::BSplinePreBasis< GV >::evaluateAll().

◆ getIJK()

template<typename GV >
static std::array<unsigned int,dim> Dune::Functions::BSplinePreBasis< GV >::getIJK ( typename GridView::IndexSet::IndexType  idx,
std::array< unsigned int, dim >  elements 
)
inlinestatic

Compute integer element coordinates from the element index.

Warning
This method makes strong assumptions about the grid, namely that it is structured, and that indices are given in a x-fastest fashion.

Referenced by Dune::Functions::BSplinePreBasis< GV >::indices().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 3, 22:32, 2024)