dune-functions 2.9.0
Public Types | Public Member Functions | Friends | List of all members
Dune::Functions::BSplineLocalBasis< GV, R > Class Template Reference

LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch to a knot span. More...

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

Inheritance diagram for Dune::Functions::BSplineLocalBasis< GV, R >:
Inheritance graph

Public Types

typedef LocalBasisTraits< D, dim, FieldVector< D, dim >, R, 1, FieldVector< R, 1 >, FieldMatrix< R, 1, dim > > Traits
 export type traits for function signature More...
 

Public Member Functions

 BSplineLocalBasis (const BSplinePreBasis< GV > &preBasis, const BSplineLocalFiniteElement< GV, R > &lFE)
 Constructor with a given B-spline patch. More...
 
void evaluateFunction (const FieldVector< D, dim > &in, std::vector< FieldVector< R, 1 > > &out) const
 Evaluate all shape functions. More...
 
void evaluateJacobian (const FieldVector< D, dim > &in, std::vector< FieldMatrix< D, 1, dim > > &out) const
 Evaluate Jacobian of all shape functions. More...
 
template<size_t k>
void evaluate (const typename std::array< int, k > &directions, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
 Evaluate all shape functions and derivatives of any order. More...
 
unsigned int order () const
 Polynomial order of the shape functions. More...
 
std::size_t size () const
 Return the number of basis functions on the current knot span. More...
 

Friends

class BSplineLocalFiniteElement< GV, R >
 

Detailed Description

template<class GV, class R>
class Dune::Functions::BSplineLocalBasis< GV, R >

LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch to a knot span.

Template Parameters
GVGrid view that the basis is defined on
RNumber type used for spline function values

Member Typedef Documentation

◆ Traits

template<class GV , class R >
typedef LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>, FieldMatrix<R,1,dim> > Dune::Functions::BSplineLocalBasis< GV, R >::Traits

export type traits for function signature

Constructor & Destructor Documentation

◆ BSplineLocalBasis()

template<class GV , class R >
Dune::Functions::BSplineLocalBasis< GV, R >::BSplineLocalBasis ( const BSplinePreBasis< GV > &  preBasis,
const BSplineLocalFiniteElement< GV, R > &  lFE 
)
inline

Constructor with a given B-spline patch.

The patch object does all the work.

Member Function Documentation

◆ evaluate()

template<class GV , class R >
template<size_t k>
void Dune::Functions::BSplineLocalBasis< GV, R >::evaluate ( const typename std::array< int, k > &  directions,
const typename Traits::DomainType &  in,
std::vector< typename Traits::RangeType > &  out 
) const
inline

Evaluate all shape functions and derivatives of any order.

◆ evaluateFunction()

template<class GV , class R >
void Dune::Functions::BSplineLocalBasis< GV, R >::evaluateFunction ( const FieldVector< D, dim > &  in,
std::vector< FieldVector< R, 1 > > &  out 
) const
inline

Evaluate all shape functions.

Parameters
inCoordinates where to evaluate the functions, in local coordinates of the current knot span

◆ evaluateJacobian()

template<class GV , class R >
void Dune::Functions::BSplineLocalBasis< GV, R >::evaluateJacobian ( const FieldVector< D, dim > &  in,
std::vector< FieldMatrix< D, 1, dim > > &  out 
) const
inline

Evaluate Jacobian of all shape functions.

Parameters
inCoordinates where to evaluate the Jacobian, in local coordinates of the current knot span

◆ order()

template<class GV , class R >
unsigned int Dune::Functions::BSplineLocalBasis< GV, R >::order ( ) const
inline

Polynomial order of the shape functions.

Unfortunately, the general interface of the LocalBasis class mandates that the 'order' method takes no arguments, and returns a single integer. It therefore cannot reflect that fact that a B-spline basis function can easily have different orders in the different coordinate directions. We therefore take the conservative choice and return the maximum over the orders of all directions.

◆ size()

template<class GV , class R >
std::size_t Dune::Functions::BSplineLocalBasis< GV, R >::size ( ) const
inline

Return the number of basis functions on the current knot span.

Friends And Related Function Documentation

◆ BSplineLocalFiniteElement< GV, R >

template<class GV , class R >
friend class BSplineLocalFiniteElement< GV, R >
friend

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