Dune::C1SimpleBasis< LB, Geo > Class Template Reference

Basis implementation where local and global values are equal. More...

#include <dune/localfunctions/common/simplebasis.hh>

Inheritance diagram for Dune::C1SimpleBasis< LB, Geo >:
Inheritance graph

List of all members.

Public Types

typedef BasisTraits
< C1SimpleBasis< LB, Geo > > 
Traits
 Export type traits.

Public Member Functions

 C1SimpleBasis (const LB &lb_, const Geo &geo_)
 construct a C1SimpleBasis from a localbasis
unsigned int size () const
void evaluateFunction (const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
unsigned int order () const
void evaluateJacobian (const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
void evaluateJacobianCoeffs (const typename Traits::DomainType &in, const std::vector< C > &coeffs) const
 Evaluate jacobian with given coefficients at given position.
template<typename C >
Traits::RangeType evaluateCoeffs (const typename Traits::DomainType &in, const std::vector< C > &coeffs) const
 Evaluate basis with coefficients at given position.

Detailed Description

template<typename LB, typename Geo>
class Dune::C1SimpleBasis< LB, Geo >

Basis implementation where local and global values are equal.

It simply takes a localbasis and returns whatever that returns. This is most appropriate for single component functions. For multicomponent functions there is typically some transformation needed. Which transformation to use depends very much on the particular base function.

Template Parameters:
LB Type of the local basis.
Geo Type of the geometry. Its interface should conform to that of Dune::Geometry.

Member Typedef Documentation

template<typename LB , typename Geo >
typedef BasisTraits<C1SimpleBasis<LB, Geo> > Dune::C1SimpleBasis< LB, Geo >::Traits

Export type traits.

Reimplemented from Dune::C1BasisInterface< C1SimpleBasis< LB, Geo > >.


Constructor & Destructor Documentation

template<typename LB , typename Geo >
Dune::C1SimpleBasis< LB, Geo >::C1SimpleBasis ( const LB &  lb_,
const Geo &  geo_ 
) [inline]

construct a C1SimpleBasis from a localbasis

Parameters:
[in] lb_ The localbasis to wrap. A reference to this object is stored, that means, the object must live at least as long as this simple basis is used for evaluation.

Member Function Documentation

template<typename Imp, typename T = BasisTraits<Imp>>
template<typename C >
Traits::RangeType Dune::C0BasisInterface< Imp, T >::evaluateCoeffs ( const typename Traits::DomainType &  in,
const std::vector< C > &  coeffs 
) const [inline, inherited]

Evaluate basis with coefficients at given position.

Evaluates all shape functions at the given (local) position and return the (global) value of the sum weighted with some coefficients.

Note:
This is a default implementation which works by perfoming the weighted sum of the values returned by evaluateFunction().
Template Parameters:
C The type of the coefficients
Parameters:
[in] in Where to evaluate in local (reference element) coordinates.
[in] coeffs The coefficients.
Returns:
The resulting global value.
template<typename LB , typename Geo >
void Dune::C1SimpleBasis< LB, Geo >::evaluateFunction ( const typename Traits::DomainType in,
std::vector< typename Traits::RangeType > &  out 
) const [inline]

Evaluate all basis function at given position. Evaluates all shape functions at the given (local) position and returns the (global) values in a vector.

Parameters:
[in] in Where to evaluate in local (reference element) coordinates.
[out] out The resulting global values, one per shape function.

Reimplemented from Dune::C0BasisInterface< Imp, T >.

template<typename LB , typename Geo >
void Dune::C1SimpleBasis< LB, Geo >::evaluateJacobian ( const typename Traits::DomainType in,
std::vector< typename Traits::JacobianType > &  out 
) const [inline]

Evaluate jacobian of all shape functions at given position. out[i][j][k] is $\partial_k \hat\phi_j^i $, when $\hat\phi^i $ is the i'th shape function. Note that this are the derivatives of the global values by the global coordinates, evaluated at local coordinates.

Parameters:
[in] in Where to evaluate in local coordinates.
[out] out The result, one jacobian per base function.
This implementation is for single component base functions. For multicomponent base functions, it will transform the gradient of each component independently, as if each component was an independent scalar valued base function.

The Jacobian is taken as a vector of the gradients of each component:

\[ J_f=\begin{pmatrix} (\nabla f_0)^T\\ (\nabla f_1)^T\\ \vdots \end{pmatrix} \]

Each gradient is transformed according to

\[ \nabla\phi^i_\alpha|_{\hat x} = J^{-T}_g|_{\hat x} \hat\nabla\hat\phi^i_\alpha|_{\hat x}, \]

where $\phi^i_\alpha$ is component $\alpha$ of the $i$'th shape function, $J_g$ is the Jacobian of the local-to-global map of the geometry and the hat $\hat{}$ denotes local coordinates ($\hat x$) or functions of the local basis (as opposed to this "global" basis).

Reimplemented from Dune::C1BasisInterface< C1SimpleBasis< LB, Geo > >.

References Dune::C0LocalBasisTraits< DF, n, D, RF, m, R >::dimRange.

void Dune::C1BasisInterface< C1SimpleBasis< LB, Geo > , BasisTraits<C1SimpleBasis< LB, Geo > > >::evaluateJacobianCoeffs ( const typename Traits::DomainType &  in,
const std::vector< C > &  coeffs 
) const [inline, inherited]

Evaluate jacobian with given coefficients at given position.

out[j][k] is $\partial_k \hat f_j $, where $f$ is the function given by this basis and the coefficients. Note that this are the derivatives of the global values by the global coordinates, evaluated at local coordinates.

Note:
This is a default implementation which works by perfoming the weighted sum of the values returned by evaluateJacobian().
Template Parameters:
C The type of the coefficients
Parameters:
[in] in Where to evaluate in local (reference element) coordinates.
[in] coeffs The coefficients.
Returns:
The resulting global jacobian.
template<typename LB , typename Geo >
unsigned int Dune::C1SimpleBasis< LB, Geo >::order (  )  const [inline]

Polynomial order of the shape functions.

Todo:
Gurke!

Reimplemented from Dune::C0BasisInterface< Imp, T >.

template<typename LB , typename Geo >
unsigned int Dune::C1SimpleBasis< LB, Geo >::size (  )  const [inline]

Number of shape functions.

Reimplemented from Dune::C0BasisInterface< Imp, T >.


The documentation for this class was generated from the following file:
Generated on Sat Apr 24 11:15:36 2010 for dune-localfunctions by  doxygen 1.6.3