Dune Core Modules (2.6.0)
Convert a simple scalar local basis into a global basis. More...
#include <dune/localfunctions/common/localtoglobaladaptors.hh>
Public Member Functions | |
ScalarLocalToGlobalBasisAdaptor (const LocalBasis &localBasis_, const Geometry &geometry_) | |
construct a ScalarLocalToGlobalBasisAdaptor More... | |
std::size_t | order () const |
return maximum polynomial order of the base function More... | |
void | evaluateFunction (const Traits::DomainLocal &in, std::vector< Traits::Range > &out) const |
Evaluate all shape functions at given position. | |
void | evaluateJacobian (const Traits::DomainLocal &in, std::vector< Traits::Jacobian > &out) const |
Evaluate Jacobian of all shape functions at given position. | |
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
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
- 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
-
LocalBasis Type of the local basis to adapt. Geometry Type of the local-to-global transformation.
Constructor & Destructor Documentation
◆ ScalarLocalToGlobalBasisAdaptor()
|
inline |
construct a ScalarLocalToGlobalBasisAdaptor
- Parameters
-
localBasis_ The local basis object to adapt. geometry_ The geometry object to use for adaption.
- Note
- This class stores the references passed here. Any use of this class after these references have become invalid results in undefined behaviour. The exception is that the destructor of this class may still be called.
Member Function Documentation
◆ evaluate()
|
inherited |
Evaluate derivatives of all shape functions at given position.
- Note
- Only required for Traits::diffOrder >= 2
◆ order()
|
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 assumtion for non-affine geometries is that they are still multi-linear.
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::affine().
The documentation for this class was generated from the following file:
- dune/localfunctions/common/localtoglobaladaptors.hh