18#include "dune/common/fmatrix.hh"
28 template <
class ScalarA,
class ScalarB,
int r,
int c>
33 res[i][j] =
static_cast<ScalarA
>(m[i][j]);
60 template <
class Cell,
class Scalar>
63 static int const dim = Cell::dimension;
64 static int const dimw = Cell::Geometry::coorddimension;
79 inv = castMatrixEntries<Scalar>(cell_->geometry().jacobianInverseTransposed(xi));
A Converter for scalar shape functions that do not change their value under transformation from refer...
void moveTo(Cell const &cell)
Indicates that following accesses refer to the given cell.
VariationalArg< Scalar, dimw, 1 > global(VariationalArg< Scalar, dim, 1 > const &sf, int deriv=1) const
Applies the transformation to shape function value, gradient, and Hessian.
Dune::FieldMatrix< Scalar, 1, 1 > global(Dune::FieldMatrix< Scalar, 1, 1 > const &sf) const
Applies the transformation to shape function value.
ScalarConverter(Cell const &cell)
Dune::FieldMatrix< Scalar, 1, 1 > local(Dune::FieldMatrix< Scalar, 1, 1 > const &glob) const
void setLocalPosition(Dune::FieldVector< typename Cell::Geometry::ctype, dim > const &xi)
This file contains various utility functions that augment the basic functionality of Dune.
typename GridView::template Codim< 0 >::Entity Cell
The type of cells (entities of codimension 0) in the grid view.
Dune::FieldMatrix< ScalarA, r, c > castMatrixEntries(Dune::FieldMatrix< ScalarB, r, c > const &m)
Copies matrix into another matrix with possibly different entry type (e.g. from double to float).
A class that stores values, gradients, and Hessians of evaluated shape functions.
Dune::FieldMatrix< Scalar, components, dim > derivative
The shape function's spatial derivative, a components x dim matrix.
ValueType value
The shape function's value, a vector of dimension components
Tensor< Scalar, components, dim, dim > hessian
Variables and their descriptions.