Dune Core Modules (2.8.0)

Dune::Functions::BSplineLocalBasis< GV, R, MI > 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>

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
 

Public Member Functions

 BSplineLocalBasis (const BSplinePreBasis< GV, MI > &preBasis, const BSplineLocalFiniteElement< GV, R, MI > &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.
 
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.
 

Detailed Description

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

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

Constructor & Destructor Documentation

◆ BSplineLocalBasis()

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

Constructor with a given B-spline patch.

The patch object does all the work.

Member Function Documentation

◆ evaluateFunction()

template<class GV , class R , class MI >
void Dune::Functions::BSplineLocalBasis< GV, R, MI >::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

References Dune::DiagonalMatrix< K, n >::umv().

Referenced by Dune::Functions::BSplineLocalBasis< GV, R, MI >::evaluate().

◆ evaluateJacobian()

template<class GV , class R , class MI >
void Dune::Functions::BSplineLocalBasis< GV, R, MI >::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

References Dune::DiagonalMatrix< K, n >::umv().

◆ order()

template<class GV , class R , class MI >
unsigned int Dune::Functions::BSplineLocalBasis< GV, R, MI >::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.


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.111.3 (Nov 12, 23:30, 2024)