KASKADE 7 development version
Classes | Public Member Functions | List of all members
Kaskade::ShapeFunctionCache< G, T, ComponentsEnd > Class Template Reference

This class caches values and derivatives of shape functions at quadrature points. More...

#include <shapefunctioncache.hh>

Detailed Description

template<class G, class T, int ComponentsEnd = 4>
class Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >

This class caches values and derivatives of shape functions at quadrature points.

Limitations are:

Template Parameters
Gthe grid type
Tthe scalar field type
ComponentsEndupper bound on the number of components of the shape functions
Todo:
make shape function cache independent of grid type - only coordinate type and dimension are required

Definition at line 115 of file shapefunctioncache.hh.

Classes

struct  DataType
 Defines the type returned by evaluate(sfs,qr,subId). More...
 
struct  LocalDataType
 Defines the type returned by evaluate(sfs,qr,ip,subId). More...
 

Public Member Functions

template<int nComp, int subDim>
Dune::FieldVector< T, nComp > evaluateFunction (ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &sfs, int isf, Dune::QuadratureRule< typename G::ctype, subDim > const &qr, int ip, int subId)
 Evaluate comp-th component of isf-th shape function of the shape function set sfs ath the ip-th integration point of quadrature rule qr. More...
 
template<int nComp, int subDim>
Dune::FieldMatrix< T, nComp, G::dimension > const & evaluateDerivative (ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &sfs, int isf, Dune::QuadratureRule< typename G::ctype, subDim > const &qr, int ip, int subId)
 Evaluate derivative of comp-th component of isf-th shape function of the shape function set sfs ath the ip-th integration point of quadrature rule qr. More...
 
template<int nComp, int subDim>
LocalDataType< nComp >::type evaluate (ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &sfs, Dune::QuadratureRule< typename G::ctype, subDim > const &qr, int ip, int subId)
 Returns the values of all shape functions at given integration point. More...
 
template<int nComp, int subDim>
DataType< nComp >::type const & evaluate (ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &sfs, Dune::QuadratureRule< typename G::ctype, subDim > const &qr, int subId)
 Returns the values of all shape functions. More...
 
int maximumSubentityIndex () const
 Reports the maximum subentity number that is allowed. More...
 

Member Function Documentation

◆ evaluate() [1/2]

template<class G , class T , int ComponentsEnd = 4>
template<int nComp, int subDim>
LocalDataType< nComp >::type Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >::evaluate ( ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &  sfs,
Dune::QuadratureRule< typename G::ctype, subDim > const &  qr,
int  ip,
int  subId 
)
inline

Returns the values of all shape functions at given integration point.

The return value can be indexed by shape function number, yielding a pair of value and derivative of the shape function's component.

E.g., the value can be obtained by evaluate(sfs,qr,ip,subId)[isf].value

This method is explicitly NOT thread safe. Use one cache object per thread.

Parameters
subIdthe index of the subentity in the cell (has to be between 0 and the maximum supported subentity number)

Definition at line 220 of file shapefunctioncache.hh.

Referenced by Kaskade::FEFunctionSpace< LocalToGlobalMapper >::Evaluator::evaluateAt(), and Kaskade::FEFunctionSpace< LocalToGlobalMapper >::Evaluator::useQuadratureRule().

◆ evaluate() [2/2]

template<class G , class T , int ComponentsEnd = 4>
template<int nComp, int subDim>
DataType< nComp >::type const & Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >::evaluate ( ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &  sfs,
Dune::QuadratureRule< typename G::ctype, subDim > const &  qr,
int  subId 
)
inline

Returns the values of all shape functions.

The return value can be indexed by shape function number, yielding a pair of value and derivative of the shape function's component.

E.g., the value can be obtained by evaluate(sfs,qr,subId)[ip][isf].value

This method is explicitly NOT thread safe. Use one cache object per thread.

Parameters
subIdthe index of the subentity in the cell (has to be between 0 and the maximum supported subentity number)

Definition at line 245 of file shapefunctioncache.hh.

◆ evaluateDerivative()

template<class G , class T , int ComponentsEnd = 4>
template<int nComp, int subDim>
Dune::FieldMatrix< T, nComp, G::dimension > const & Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >::evaluateDerivative ( ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &  sfs,
int  isf,
Dune::QuadratureRule< typename G::ctype, subDim > const &  qr,
int  ip,
int  subId 
)
inline

Evaluate derivative of comp-th component of isf-th shape function of the shape function set sfs ath the ip-th integration point of quadrature rule qr.

This method is explicitly NOT thread safe. Use one cache object per thread.

Definition at line 196 of file shapefunctioncache.hh.

◆ evaluateFunction()

template<class G , class T , int ComponentsEnd = 4>
template<int nComp, int subDim>
Dune::FieldVector< T, nComp > Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >::evaluateFunction ( ShapeFunctionSet< typename G::ctype, G::dimension, T, nComp > const &  sfs,
int  isf,
Dune::QuadratureRule< typename G::ctype, subDim > const &  qr,
int  ip,
int  subId 
)
inline

Evaluate comp-th component of isf-th shape function of the shape function set sfs ath the ip-th integration point of quadrature rule qr.

This method is explicitly NOT thread safe. Use one cache object per thread.

Definition at line 180 of file shapefunctioncache.hh.

◆ maximumSubentityIndex()

template<class G , class T , int ComponentsEnd = 4>
int Kaskade::ShapeFunctionCache< G, T, ComponentsEnd >::maximumSubentityIndex ( ) const
inline

Reports the maximum subentity number that is allowed.

Definition at line 256 of file shapefunctioncache.hh.


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