KASKADE 7 development version
|
A Converter for scalar shape functions that do not change their value under transformation from reference to actual element geometry. More...
#include <converter.hh>
A Converter for scalar shape functions that do not change their value under transformation from reference to actual element geometry.
If grid and world dimension are the same, the restriction of ansatz functions \( \varphi \) on the actual element are defined in terms of the shape functions \( \phi \) on the reference elements as
\[ \varphi(x) = \phi(B^{-1}(x-b)), \quad \xi = B^{-1}(x-b). \]
Therefore the derivatives and gradients are transformed as
\[ \varphi'(x) = \phi'(\xi)B^{-1}, \quad \nabla \varphi(x) = B^{-T}\nabla\phi(\xi). \]
Hessians (derivatives of gradients) are transformed as
\[ H\varphi(x) = B^{-T}H_\phi(\xi)B^{-1}. \]
This is used as a converter for scalar Lagrange and hierarchical spaces.
If the world dimension is larger than the grid dimension, the gradient is defined as the gradient of the function extending as a constant in normal direction. This is just
\[ \varphi'(x) = \phi'(\xi)B^{+}, \quad \nabla \varphi(x) = B^{+T}\nabla\phi(\xi), \]
where \( B \in {\bf R}^{w\times d}\) .
Template parameters:
Cell | the type of codim 0 entity on which the shape functions are defined. |
Definition at line 61 of file converter.hh.
Public Member Functions | |
ScalarConverter () | |
ScalarConverter (Cell const &cell) | |
void | moveTo (Cell const &cell) |
Indicates that following accesses refer to the given cell. More... | |
void | setLocalPosition (Dune::FieldVector< typename Cell::Geometry::ctype, dim > const &xi) |
Dune::FieldMatrix< Scalar, 1, 1 > | global (Dune::FieldMatrix< Scalar, 1, 1 > const &sf) const |
Applies the transformation \( \psi(x) \) to shape function value. More... | |
VariationalArg< Scalar, dimw, 1 > | global (VariationalArg< Scalar, dim, 1 > const &sf, int deriv=1) const |
Applies the transformation \( \psi \) to shape function value, gradient, and Hessian. More... | |
Dune::FieldMatrix< Scalar, 1, 1 > | local (Dune::FieldMatrix< Scalar, 1, 1 > const &glob) const |
|
inline |
Definition at line 67 of file converter.hh.
|
inline |
Definition at line 69 of file converter.hh.
|
inline |
Applies the transformation \( \psi(x) \) to shape function value.
Definition at line 83 of file converter.hh.
Referenced by Kaskade::MorleyMapper< ScalarType, GV >::Combiner::Combiner().
|
inline |
Applies the transformation \( \psi \) to shape function value, gradient, and Hessian.
Definition at line 88 of file converter.hh.
|
inline |
Definition at line 97 of file converter.hh.
|
inline |
Indicates that following accesses refer to the given cell.
Definition at line 74 of file converter.hh.
|
inline |
Definition at line 76 of file converter.hh.
Referenced by Kaskade::MorleyMapper< ScalarType, GV >::Combiner::Combiner().