DUNE PDELab (git)

Interpolation for Brezzi-Douglas-Fortin-Marini shape functions on cubes. More...

#include <dune/localfunctions/brezzidouglasfortinmarini/cube/localinterpolation.hh>

Public Member Functions

 BDFMCubeLocalInterpolation ()
 Standard constructor.
 
 BDFMCubeLocalInterpolation (std::bitset< numFaces > s)
 Make set number s, where 0 <= s < 2^(2*dim) More...
 
template<class F , class C >
void interpolate (const F &f, C &out) const
 Interpolate a given function with shape functions. More...
 
template<class F , class C >
void trace (unsigned int face, const F &f, C &out) const
 Interpolate a given function with shape functions on a face of the reference element. More...
 
template<class F , class C >
void interior (const F &f, C &out) const
 Interpolate a given function with shape functions in the interior of the reference element. More...
 

Detailed Description

template<class D, class R, unsigned int dim, unsigned int order>
class Dune::BDFMCubeLocalInterpolation< D, R, dim, order >

Interpolation for Brezzi-Douglas-Fortin-Marini shape functions on cubes.

Template Parameters
DType of represent the field in the domain.
RType of represent the field in the domain.
dimdimension of the reference element, must be >= 2.
orderorder of the element, must be >= 1.

Constructor & Destructor Documentation

◆ BDFMCubeLocalInterpolation()

template<class D , class R , unsigned int dim, unsigned int order>
Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::BDFMCubeLocalInterpolation ( std::bitset< numFaces >  s)
inline

Make set number s, where 0 <= s < 2^(2*dim)

Parameters
sEdge orientation indicator

Member Function Documentation

◆ interior()

template<class D , class R , unsigned int dim, unsigned int order>
template<class F , class C >
void Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interior ( const F &  f,
C &  out 
) const
inline

Interpolate a given function with shape functions in the interior of the reference element.

Template Parameters
FFunction type for function which should be interpolated
CCoefficient vector type
Parameters
ffunction which should be interpolated
outreturn value, vector of coefficients

References Dune::GeometryTypes::cube(), Dune::Hybrid::max, and Dune::QuadratureRules< ctype, dim >::rule().

Referenced by Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interpolate().

◆ interpolate()

template<class D , class R , unsigned int dim, unsigned int order>
template<class F , class C >
void Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interpolate ( const F &  f,
C &  out 
) const
inline

Interpolate a given function with shape functions.

Template Parameters
FFunction type for function which should be interpolated
CCoefficient vector type
Parameters
ffunction which should be interpolated
outreturn value, vector of coefficients

References Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interior(), and Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::trace().

◆ trace()

template<class D , class R , unsigned int dim, unsigned int order>
template<class F , class C >
void Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::trace ( unsigned int  face,
const F &  f,
C &  out 
) const
inline

Interpolate a given function with shape functions on a face of the reference element.

Template Parameters
FFunction type for function which should be interpolated
CCoefficient vector type
Parameters
faceindex of the face on which to interpolate
ffunction which should be interpolated
outreturn value, vector of coefficients

References Dune::GeometryTypes::cube(), and Dune::QuadratureRules< ctype, dim >::rule().

Referenced by Dune::BDFMCubeLocalInterpolation< D, R, dim, order >::interpolate().


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 23, 23:29, 2024)