dune-functions 2.9.0
|
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. More... | |
using | size_type = std::size_t |
using | Node = BSplineNode< GV > |
using | R = double |
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. More... | |
const GridView & | gridView () const |
Obtain the grid view that the basis is defined on. More... | |
void | update (const GridView &gv) |
Update the stored grid view, to be called if the grid has changed. More... | |
Node | makeNode () const |
Create tree node. More... | |
template<class ST , int i> | |
size_type | size (const Dune::ReservedVector< ST, i > &prefix) const |
Return number of possible values for next position in multi index. More... | |
size_type | dimension () const |
Get the total dimension of the space spanned by this basis. More... | |
size_type | maxNodeSize () const |
Get the maximal number of DOFs associated to node for any element. More... | |
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. More... | |
unsigned int | size () const |
Total number of B-spline basis functions. More... | |
unsigned int | size (size_t d) const |
Number of shape functions in one direction. More... | |
void | evaluateFunction (const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldVector< R, 1 > > &out, const std::array< unsigned, dim > ¤tKnotSpan) const |
Evaluate all B-spline basis functions at a given point. More... | |
void | evaluateJacobian (const FieldVector< typename GV::ctype, dim > &in, std::vector< FieldMatrix< R, 1, dim > > &out, const std::array< unsigned, dim > ¤tKnotSpan) 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 > ¤tKnotSpan) const |
Evaluate Derivatives of all B-spline basis functions. More... | |
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. More... | |
std::array< std::vector< double >, dim > | knotVectors_ |
The knot vectors, one for each space dimension. More... | |
std::array< unsigned, dim > | elements_ |
Number of grid elements in the different coordinate directions. More... | |
GridView | gridView_ |
Static Public Attributes | |
static constexpr size_type | maxMultiIndexSize = 1 |
static constexpr size_type | minMultiIndexSize = 1 |
static constexpr size_type | multiIndexBufferSize = 1 |
Pre-basis for B-spline basis.
GV | The 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.
using Dune::Functions::BSplinePreBasis< GV >::GridView = GV |
The grid view that the FE space is defined on.
using Dune::Functions::BSplinePreBasis< GV >::Node = BSplineNode<GV> |
using Dune::Functions::BSplinePreBasis< GV >::R = double |
using Dune::Functions::BSplinePreBasis< GV >::size_type = std::size_t |
|
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.:
Unfortunately, not all of these conditions can be checked for automatically.
knotVector | A single knot vector, which will be used for all coordinate directions |
order | B-spline order, will be used for all coordinate directions |
makeOpen | If 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. |
|
inline |
Construct a B-spline basis for a given grid view with uniform knot vectors.
The grid must match the knot vectors, i.e.:
Unfortunately, not all of these conditions can be checked for automatically.
gridView | The grid we are defining the basis on |
lowerLeft | Lower left corner of the structured grid |
upperRight | Upper right corner of the structured grid |
elements | Number of elements in each coordinate direction |
order | B-spline order, will be used for all coordinate directions |
makeOpen | If 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. |
|
inline |
Get the total dimension of the space spanned by this basis.
|
inline |
Evaluate Derivatives of all B-spline basis functions.
|
inlinestatic |
Evaluate the second derivatives of all one-dimensional B-spline functions for a given coordinate direction.
in | Scalar(!) coordinate where to evaluate the functions | |
enableEvaluations | switches calculation of desired derivatives on | |
[out] | out | Vector containing the values of all B-spline derivatives at 'in' |
[out] | outJac | Vector containing the first derivatives of all B-spline derivatives at 'in' (only if calculation was switched on by enableEvaluations) |
[out] | outHess | Vector containing the second derivatives of all B-spline derivatives at 'in' (only if calculation was switched on by enableEvaluations) |
|
inline |
Evaluate all B-spline basis functions at a given point.
|
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"
in | Scalar(!) coordinate where to evaluate the functions | |
[out] | out | Vector containing the values of all B-spline functions at 'in' |
|
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"
in | Scalar(!) coordinate where to evaluate the functions | |
[out] | out | Vector containing the values of all B-spline functions at 'in' |
|
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.
|
inlinestatic |
Compute integer element coordinates from the element index.
|
inline |
Obtain the grid view that the basis is defined on.
|
inline |
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
|
inline |
Initialize the global indices.
|
inline |
Create tree node.
|
inline |
Get the maximal number of DOFs associated to node for any element.
|
inline |
Total number of B-spline basis functions.
|
inline |
Return number of possible values for next position in multi index.
|
inline |
Number of shape functions in one direction.
|
inline |
Update the stored grid view, to be called if the grid has changed.
std::array<unsigned,dim> Dune::Functions::BSplinePreBasis< GV >::elements_ |
Number of grid elements in the different coordinate directions.
GridView Dune::Functions::BSplinePreBasis< GV >::gridView_ |
std::array<std::vector<double>, dim> Dune::Functions::BSplinePreBasis< GV >::knotVectors_ |
The knot vectors, one for each space dimension.
|
staticconstexpr |
|
staticconstexpr |
|
staticconstexpr |
std::array<unsigned int, dim> Dune::Functions::BSplinePreBasis< GV >::order_ |
Order of the B-spline for each space dimension.