dune-localfunctions  2.4.1
Public Member Functions | List of all members
Dune::ScalarLocalToGlobalBasisAdaptor< LocalBasis, Geometry > Class Template Reference

Convert a simple scalar local basis into a global basis. More...

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

Inheritance diagram for Dune::ScalarLocalToGlobalBasisAdaptor< LocalBasis, Geometry >:
Inheritance graph

Public Member Functions

std::size_t size () const
 Number of shape functions. More...
 
std::size_t order () const
 Polynomial order of the shape functions for quadrature. More...
 
void evaluateFunction (const Traits::DomainLocal &in, std::vector< Traits::Range > &out) const
 Evaluate all shape functions at given position. More...
 
void evaluateJacobian (const Traits::DomainLocal &in, std::vector< Traits::Jacobian > &out) const
 Evaluate Jacobian of all shape functions at given position. More...
 
void evaluate (const std::array< std::size_t, Traits::dimDomainGlobal > &directions, const Traits::DomainLocal &in, std::vector< Traits::Range > &out) const
 Evaluate derivatives of all shape functions at given position. More...
 

Detailed Description

template<class LocalBasis, class Geometry>
class Dune::ScalarLocalToGlobalBasisAdaptor< LocalBasis, Geometry >

Convert a simple scalar local basis into a global basis.

The local basis must be scalar, i.e. LocalBasis::Traits::dimRange must be

  1. It's values are not transformed.

For scalar function $f$, the gradient is equivalent to the transposed Jacobian $\nabla f|_x = J_f^T(x)$. The Jacobian is thus transformed using

\[ \nabla f|_{\mu(\hat x)} = \hat J_\mu^{-T}(\hat x) \cdot \hat\nabla\hat f|_{\hat x} \]

Here the hat $\hat{\phantom x}$ denotes local quantities and $\mu$ denotes the local-to-global map of the geometry.

Template Parameters
LocalBasisType of the local basis to adapt.
GeometryType of the local-to-global transformation.

Member Function Documentation

void Dune::BasisInterface::evaluate ( const std::array< std::size_t, Traits::dimDomainGlobal > &  directions,
const Traits::DomainLocal in,
std::vector< Traits::Range > &  out 
) const
inherited

Evaluate derivatives of all shape functions at given position.

Note
Only required for Traits::diffOrder >= 2
void Dune::BasisInterface::evaluateFunction ( const Traits::DomainLocal in,
std::vector< Traits::Range > &  out 
) const
inherited

Evaluate all shape functions at given position.

void Dune::BasisInterface::evaluateJacobian ( const Traits::DomainLocal in,
std::vector< Traits::Jacobian > &  out 
) const
inherited

Evaluate Jacobian of all shape functions at given position.

Note: Only required for Traits::diffOrder >= 1

std::size_t Dune::BasisInterface::order ( ) const
inherited

Polynomial order of the shape functions for quadrature.

std::size_t Dune::BasisInterface::size ( ) const
inherited

Number of shape functions.


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