KASKADE 7 development version
converter.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2019 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef CONVERTER_HH
14#define CONVERTER_HH
15
16#include <cassert>
17
18#include "dune/common/fmatrix.hh"
19
20#include "fem/fixdune.hh"
21#include "fem/variables.hh"
22
23namespace Kaskade
24{
28 template <class ScalarA, class ScalarB, int r, int c>
31 for(int i=0;i<r;++i)
32 for(int j=0;j<r;++j)
33 res[i][j] = static_cast<ScalarA>(m[i][j]);
34 return res;
35 }
36
60 template <class Cell, class Scalar>
62 {
63 static int const dim = Cell::dimension;
64 static int const dimw = Cell::Geometry::coorddimension;
65
66 public:
67 ScalarConverter(): cell_(0) {}
68
69 ScalarConverter(Cell const& cell): cell_(&cell) {}
70
74 void moveTo(Cell const& cell) { cell_ = &cell; }
75
77 {
78 assert(cell_);
79 inv = castMatrixEntries<Scalar>(cell_->geometry().jacobianInverseTransposed(xi));
80 }
81
84
89 {
91 phi.value = sf.value;
92 if (deriv >= 1) phi.derivative[0] = inv * sf.derivative[0];
93 if (deriv >= 2) phi.hessian[0] = inv * static_cast<Dune::FieldMatrix<Scalar,dim,dim>>(sf.hessian[0]) * transpose(inv);
94 return phi;
95 }
96
98
99
100 private:
101 Cell const* cell_;
103 };
104} /* end of namespace Kaskade */
105
106
107#endif
A Converter for scalar shape functions that do not change their value under transformation from refer...
Definition: converter.hh:62
void moveTo(Cell const &cell)
Indicates that following accesses refer to the given cell.
Definition: converter.hh:74
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.
Definition: converter.hh:88
Dune::FieldMatrix< Scalar, 1, 1 > global(Dune::FieldMatrix< Scalar, 1, 1 > const &sf) const
Applies the transformation to shape function value.
Definition: converter.hh:83
ScalarConverter(Cell const &cell)
Definition: converter.hh:69
Dune::FieldMatrix< Scalar, 1, 1 > local(Dune::FieldMatrix< Scalar, 1, 1 > const &glob) const
Definition: converter.hh:97
void setLocalPosition(Dune::FieldVector< typename Cell::Geometry::ctype, dim > const &xi)
Definition: converter.hh:76
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.
Definition: gridBasics.hh:35
T transpose(T x)
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).
Definition: converter.hh:29
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.