KASKADE 7 development version
|
Tools for converting finite element functions. More...
Modules | |
Functionviews | |
Predefined function views. | |
Classes | |
class | Kaskade::AdvectedFunctionView< Function, Velocity > |
A weak function view that defines an advected function. More... | |
struct | Kaskade::AdaptationCoarseningPolicy |
policy class for LocalTransfer These two structures are used for determining whether to perform a "normal" grid transfer (default, AdaptionCoarseningPolicy), or to calculate prolongation matrices for a FE function space based on a HierarchicIndexSet (see mgtools.hh). More... | |
struct | Kaskade::MultilevelCoarseningPolicy |
policy class for LocalTransfer These two structures are used for determining whether to perform a "normal" grid transfer (default, AdaptionCoarseningPolicy), or to calculate prolongation matrices for a FE function space based on a HierarchicIndexSet (see mgtools.hh). More... | |
class | Kaskade::LocalTransfer< Space, CoarseningPolicy > |
Class that stores for a cell ("this entity" or "father") a local projection matrix. More... | |
class | Kaskade::TransferData< Space, CoarseningPolicy > |
A class storing data collected before grid adaptation that is necessary to transfer FE functions. This includes global indices of shape functions on leaf cells and restriction matrices to coarser cells in case some leaf cells might be removed during grid adaptation. More... | |
struct | Kaskade::InverseVolume |
A Weighing that associates to each cell its inverse volume. More... | |
struct | Kaskade::Volume |
A Weighing that associates to each cell its volume. More... | |
struct | Kaskade::PlainAverage |
A Weighing that associates to each cell a constant weight of 1. More... | |
class | Kaskade::WeakFunctionViewAdaptor< Functor > |
An adaptor that turns callables (e.g., lambdas) into weak function views. More... | |
class | Kaskade::FunctionViewAdaptor< Space_, Functor > |
An adaptor that allows to turn lambda functions into function views. More... | |
Functions | |
template<class Space , class ShapeFunctionSet > | |
void | Kaskade::evaluateGlobalShapeFunctions (Space const &space, Cell< typename Space::GridView > const &cell, std::vector< LocalPosition< typename Space::GridView > > const &xi, DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, Space::sfComponents, 1 > > &sfValues, ShapeFunctionSet const &shapeFunctionSet) |
Computes the values of global shape functions at the given points (which are supposed to be local coordinates inside the cell to which the evaluator is currently moved). More... | |
template<class Space , class ShapeFunctionSet > | |
void | Kaskade::approximateGlobalValues (Space const &space, typename Space::Grid::template Codim< 0 >::Entity const &cell, DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, Space::sfComponents, 1 > > &globalValues, DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, 1, 1 > > &coeff, ShapeFunctionSet const &sfs) |
Computes global shape function coefficients such that the given set of global function values at the interpolation nodes of the shape function set are approximated by the resulting linear combinations of global shape functions. More... | |
template<class ChildSpace , class FatherSpace > | |
void | Kaskade::localProlongationMatrix (ChildSpace const &childSpace, Cell< typename ChildSpace::GridView > const &child, FatherSpace const &fatherSpace, Cell< typename FatherSpace::GridView > const &father, DynamicMatrix< Dune::FieldMatrix< typename ChildSpace::Scalar, 1, 1 > > &prolongation) |
Computes a local prolongation matrix from father cell to child cell. More... | |
template<typename Weighing = PlainAverage, typename FSElement , typename Function > | |
void | Kaskade::interpolateGlobally (FSElement &fse, Function const &fu) |
Interpolates FunctionSpaceElement to FunctionSpaceElement. More... | |
template<typename Weighing = PlainAverage, typename Target , typename Source > | |
void | Kaskade::interpolateGloballyWeak (Target &target, Source const &fu) |
Interpolates WeakFunctionViews to FunctionSpaceElement. More... | |
template<class FromSpace , class ToSpace > | |
Prolongation | Kaskade::twoGridProlongation (FromSpace const &from, ToSpace const &to) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different grids. More... | |
template<class FromSpace , class ToSpace > | |
Prolongation | Kaskade::nonNestedProlongation (FromSpace const &from, ToSpace const &to, double tolTrunc=0.1) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different (possibly non-nested) grids. More... | |
std::vector< Prolongation > | Kaskade::nonNestedProlognationStack (std::vector< Prolongation > ps) |
Computes a full column rank prolongation stack from a colum rank deficient prolongation stack. More... | |
template<class Space > | |
std::vector< Prolongation > | Kaskade::makeNonNestedPstack (std::initializer_list< Space > spacelist, double tolTrunc=0.1) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different (possibly non-nested) grids. More... | |
template<class Space1 , class Space2 > | |
std::vector< Prolongation > | Kaskade::makeJointPstack (Space1 space1, Space2 space2, int coarseLevel=0) |
Computes a (joint) prolongation stack for a geometry defined by two grids. More... | |
template<class Function , class Velocity > | |
auto | makeAdvectedFunctionView (Function const &f, Velocity const &v, double tau, int n=1) |
Creates a weak function view for advecting FE functions. More... | |
template<class Functor > | |
auto | makeWeakFunctionView (Functor const &g) |
A convenience function turning callables (e.g., lambda functions) into weak function views. More... | |
template<class Space , class Functor > | |
auto | makeFunctionView (Space const &space, Functor const &g) |
A convenience functor that supports the easy creation of function view adaptors. More... | |
Tools for converting finite element functions.
void Kaskade::approximateGlobalValues | ( | Space const & | space, |
typename Space::Grid::template Codim< 0 >::Entity const & | cell, | ||
DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, Space::sfComponents, 1 > > & | globalValues, | ||
DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, 1, 1 > > & | coeff, | ||
ShapeFunctionSet const & | sfs | ||
) |
Computes global shape function coefficients such that the given set of global function values at the interpolation nodes of the shape function set are approximated by the resulting linear combinations of global shape functions.
This realizes the matrix product \( \Phi^{+} \Psi^{-1} \) ( \( \Psi \) is the converter, and \( \Phi \) the values of the shape functions at the interpolation nodes, see LocalToGlobalMapperConcept).
On exit, the value of globalValues is undefined.
The actual meaning of "approximate" (formally denoted by a pseudoinverse symbol here) is defined by the specific FE space. However, for any space it holds that after executing
the matrix coeff is the identity matrix.
[in] | globalValues | global function values \( f(\xi_i) \) at interpolation nodes \( \xi_i \) of the shape function set (will be overwritten and invalidated) |
[out] | coeff | coefficients \( C \) of shape function linear combination such that the global values at the interpolation points are approximated: \( f(\xi_i) = \sum_j C_{ij} \phi_j(\xi_i) \) |
sfs | the space's shape function set on the given cell (note that the cell needs not be contained in the space's grid view at all, e.g. when creating multilevel prolongations) |
Definition at line 133 of file fetransfer.hh.
Referenced by Kaskade::approximateGlobalValues(), Kaskade::localProlongationMatrix(), Kaskade::LocalTransfer< Space, CoarseningPolicy >::LocalTransfer(), and Kaskade::prolongation().
void Kaskade::evaluateGlobalShapeFunctions | ( | Space const & | space, |
Cell< typename Space::GridView > const & | cell, | ||
std::vector< LocalPosition< typename Space::GridView > > const & | xi, | ||
DynamicMatrix< Dune::FieldMatrix< typename Space::Scalar, Space::sfComponents, 1 > > & | sfValues, | ||
ShapeFunctionSet const & | shapeFunctionSet | ||
) |
Computes the values of global shape functions at the given points (which are supposed to be local coordinates inside the cell to which the evaluator is currently moved).
This realizes the matrix product \(\Psi\Phi \) ( \( \Phi \) is the matrix of shape function values at the given points \( \xi \) in the reference element and \(\Psi\) is the converter, see LocalToGlobalMapperConcept). On return,
See also approximateGlobalValues.
Definition at line 69 of file fetransfer.hh.
Referenced by Kaskade::evaluateGlobalShapeFunctions(), Kaskade::localProlongationMatrix(), Kaskade::LocalTransfer< Space, CoarseningPolicy >::LocalTransfer(), and Kaskade::prolongation().
void Kaskade::interpolateGlobally | ( | FSElement & | fse, |
Function const & | fu | ||
) |
Interpolates FunctionSpaceElement to FunctionSpaceElement.
This is particulary useful together with FunctionViews and makes it possible to transform a FunctionView to a FunctionSpaceElement
fse | function space element which is set to (approximately) the value of fu |
fu | function space element of a (probably different) FE space |
Weighing | has to be default constructible |
FSElement | type of the target finite element function |
Function | type of the source function |
Definition at line 855 of file fetransfer.hh.
Referenced by Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator=(), and Kaskade::spaceTransfer().
void Kaskade::interpolateGloballyWeak | ( | Target & | target, |
Source const & | fu | ||
) |
Interpolates WeakFunctionViews to FunctionSpaceElement.
This is particulary useful together with WeakFunctionViews and makes it possible to transform a WeakFunctionView to a FunctionSpaceElement.
Weighing | the type of averaging |
Target | the target function type |
Source | the source weak function view type |
Remember that WeakFunctionView s support the evaluation by cell and local coordinate. If the source function is a finite element function, prefer interpolateGlobally, since this is a little bit faster.
Definition at line 946 of file fetransfer.hh.
Referenced by Kaskade::DeformingGridManager< Grid >::correctingDisplacement(), and Kaskade::FunctionSpaceElement< FunctionSpace, m >::operator=().
void Kaskade::localProlongationMatrix | ( | ChildSpace const & | childSpace, |
Cell< typename ChildSpace::GridView > const & | child, | ||
FatherSpace const & | fatherSpace, | ||
Cell< typename FatherSpace::GridView > const & | father, | ||
DynamicMatrix< Dune::FieldMatrix< typename ChildSpace::Scalar, 1, 1 > > & | prolongation | ||
) |
Computes a local prolongation matrix from father cell to child cell.
child | a pointer to a cell contained in the child space's grid view |
father | a pointer to a cell contained in the father space's grid view |
Definition at line 232 of file fetransfer.hh.
|
related |
A convenience functor that supports the easy creation of function view adaptors.
Space | a FEFunctionSpace type |
Functor | a callable object with operator()(Space::Evaluator eval) |
Definition at line 1107 of file fetransfer.hh.
std::vector< Prolongation > Kaskade::makeJointPstack | ( | Space1 | space1, |
Space2 | space2, | ||
int | coarseLevel = 0 |
||
) |
Computes a (joint) prolongation stack for a geometry defined by two grids.
Space1 | FEFunctionSpace built on gridManager of first geometry |
Space2 | FEFunctionSpace built on gridManager of second geometry |
space1 | space of first geometry |
space2 | space of second geometry |
This function computes a prolongation stack for a hierarchy based on two (hierarchical) spaces. At least one of the spaces should be hierarchical, i.e. have more than one level. If the number of levels of the two spaces differ, the hierarchy of the space with fewer levels is "augmeted" by identity prolongation matrices. For each of the spaces this function calls makeDeepProlongationStack(...) and concatenates the two prolongation stacks.
std::vector< Prolongation > Kaskade::makeNonNestedPstack | ( | std::initializer_list< Space > | spacelist, |
double | tolTrunc = 0.1 |
||
) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different (possibly non-nested) grids.
FromSpace | a FEFunctionSpace of continuous piecewise linear elements |
ToSpace | a FEFunctionSpace of continuous piecewise linear elements |
from | the coarse space |
to | the fine space |
The computed interpolation matrix transfers coefficient vectors of the coarse space to coefficient vectors of the fine space such that the function is (approximately) the same. Thus, it is a representation of (approximate) identity w.r.t. the Lagrangian bases of the involved spaces. Since the spaces might be non-nested, there are two undesider case that can come up:
|
related |
A convenience function turning callables (e.g., lambda functions) into weak function views.
This is a model of the WeakFunctionView concept.
Functor | a callable object with operator()(Cell, LocalPosition) |
Definition at line 1059 of file fetransfer.hh.
std::vector< Prolongation > Kaskade::nonNestedProlognationStack | ( | std::vector< Prolongation > | ps | ) |
Computes a full column rank prolongation stack from a colum rank deficient prolongation stack.
ps | vector containing the prolognation stack. |
The prolognation stack obtained from non-nested meshes may contain zero columns (in particular if a coarse grid cell is not contained in the fine mesh). To avoid this, we "erase" the corresponding DOFs on an algebraic level, i.e. by deleting the corresponding column in the prolognation matrix. In the prolognation matrix on the next coarser level in the stack we need to delete the corresponding rows, too. Doing this on each of the levels we obtain a full column rank prolognation stack.
Definition at line 225 of file amg.hh.
Referenced by Kaskade::makeNonNestedPstack().
Prolongation Kaskade::nonNestedProlongation | ( | FromSpace const & | from, |
ToSpace const & | to, | ||
double | tolTrunc = 0.1 |
||
) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different (possibly non-nested) grids.
FromSpace | a FEFunctionSpace of continuous piecewise linear elements |
ToSpace | a FEFunctionSpace of continuous piecewise linear elements |
from | the coarse space |
to | the fine space |
tolTrunc | tolerance for truncation of P-matrix entries (entries below tol*rowtotal will be deleted) |
The computed interpolation matrix transfers coefficient vectors of the coarse space to coefficient vectors of the fine space such that the function is (approximately) the same. Thus, it is a representation of (approximate) identity w.r.t. the Lagrangian bases of the involved spaces. Since the spaces might be non-nested, there are two undesired cases that can come up:
To retain good convergence of MG methods based on this prolongation matix, a truncation of "small" entries is required. In each row, we delete all entries, which are smaller than tolTrunc*rowtotal, with rowtotal being the sum of the absolute values of the row entries. Afterwards we scale each row to retain its original row total. The value of 0.2 seems to be a good guess, which balances the two competing goals of operator complexity and a good convergence rate.
Cf. Dickopf, Krause - "Numerical Study of the Almost Nested Case in a Multlevel Method Based on Non-nested Meshes" (2013)
Definition at line 125 of file amg.hh.
Referenced by Kaskade::makeNonNestedPstack().
Prolongation Kaskade::twoGridProlongation | ( | FromSpace const & | from, |
ToSpace const & | to | ||
) |
Computes an interpolation matrix for transfer between P1 finite element spaces on two different grids.
FromSpace | a FEFunctionSpace of continuous piecewise linear elements |
ToSpace | a FEFunctionSpace of continuous piecewise linear elements |
from | the coarse space |
to | the fine space |
The computed interpolation matrix transfers coefficient vectors of the coarse space to coefficient vectors of the fine space such that the function is (approximately) the same. Thus, it is a representation of (approximate) identity w.r.t. the Lagrangian bases of the involved spaces.
An example use is in auxiliary space methods, where
Definition at line 51 of file amg.hh.
Referenced by Kaskade::makeAuxiliarySpaceMultigridStack().