dune-localfunctions 2.10
|
Convert a simple scalar local basis into a global basis. More...
#include <dune/localfunctions/common/localtoglobaladaptors.hh>
Public Types | |
typedef LocalToGlobalBasisAdaptorTraits< typename LocalBasis::Traits, Geometry::coorddimension > | Traits |
Public Member Functions | |
ScalarLocalToGlobalBasisAdaptor (const LocalBasis &localBasis_, const Geometry &geometry_) | |
construct a ScalarLocalToGlobalBasisAdaptor | |
std::size_t | size () const |
std::size_t | order () const |
return maximum polynomial order of the base function | |
void | evaluateFunction (const typename Traits::DomainLocal &in, std::vector< typename Traits::Range > &out) const |
void | evaluateJacobian (const typename Traits::DomainLocal &in, std::vector< typename Traits::Jacobian > &out) const |
void | evaluateFunction (const Traits::DomainType &in, std::vector< Traits::RangeType > &out) const |
Evaluate all shape functions at given position. | |
void | evaluateJacobian (const Traits::DomainType &in, std::vector< Traits::Jacobian > &out) const |
Evaluate Jacobian of all shape functions at given position. | |
void | partial (const std::array< unsigned int, Traits::dimDomain > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const |
Evaluate partial derivatives of any order of all shape functions. | |
Convert a simple scalar local basis into a global basis.
The local basis must be scalar, i.e. LocalBasis::Traits::dimRange must be
For scalar function
Here the hat
LocalBasis | Type of the local basis to adapt. |
Geometry | Type of the local-to-global transformation. |
typedef LocalToGlobalBasisAdaptorTraits<typename LocalBasis::Traits, Geometry::coorddimension> Dune::ScalarLocalToGlobalBasisAdaptor< LocalBasis, Geometry >::Traits |
|
inline |
construct a ScalarLocalToGlobalBasisAdaptor
localBasis_ | The local basis object to adapt. |
geometry_ | The geometry object to use for adaption. |
|
inherited |
Evaluate all shape functions at given position.
|
inline |
|
inherited |
Evaluate Jacobian of all shape functions at given position.
|
inline |
|
inline |
return maximum polynomial order of the base function
This is to determine the required quadrature order. For an affine geometry this is the same order as for the local basis. For other geometries this returns the order of the local basis plus the global dimension minus 1. The assumption for non-affine geometries is that they are still multi-linear.
|
inherited |
Evaluate partial derivatives of any order of all shape functions.
order | Order of the partial derivatives, in the classic multi-index notation | |
in | Position where to evaluate the derivatives | |
[out] | out | Return value: the desired partial derivatives |
|
inline |