KASKADE 7 development version
Public Types | Public Member Functions | Protected Attributes | List of all members
Kaskade::MultiLevelLocalGeometry< Grid > Class Template Reference

Transformation of coordinates between ancestor and child. More...

#include <mllgeometry.hh>

Detailed Description

template<class Grid>
class Kaskade::MultiLevelLocalGeometry< Grid >

Transformation of coordinates between ancestor and child.

Provides a convenient transfer between a cell and its (possibly higher degree) ancestor.

This defines a mapping \( g:L\to G \), where \( L \) is the "local" grid cell and \( G \) the "global" grid cell. The evaluation of the Jacobian \( J = g'(x) \) is also supported.

It is possible to choose, which of child and ancestor are to be considered as "global". Moreover, it is possible to provide a flag that decides if (chained) transformations should be performed via geometry.global/local, or via localGeometry.global/local. The first option is more efficient for high degree ancestors, but it is only sensible, if the child is actually geometrically a subset of the father. This requirement might be violated in some special cases. In general the transformation from child to father is more efficient than the other way round. This has to do with the data structures in the grid.

Todo:
: implement pull-back and pull-forward operation for gradients

Definition at line 50 of file mllgeometry.hh.

Inheritance diagram for Kaskade::MultiLevelLocalGeometry< Grid >:
Kaskade::LocalGeometryInCoarseGridAncestor< Grid >

Public Types

enum  Direction { ChildIsGlobal , FatherIsGlobal }
 
typedef Grid::template Codim< 0 >::Entity Entity
 The entity type between which this geometry object maps. More...
 

Public Member Functions

 MultiLevelLocalGeometry (Grid const &grid, Cell const &child_, Cell const &ancestor_, MultiLevelLocalGeometry< Grid >::Direction direction_, bool geometricallyNested_=(dim==dimw), ctype tol_=-100)
 Constructor. More...
 
Entity const & localEntity () const
 A reference to the local entity. More...
 
void global (std::vector< Dune::FieldVector< ctype, dim > > &x, std::vector< int > &componentsInside) const
 Transformation from local to global. More...
 
void local (std::vector< Dune::FieldVector< ctype, dim > > &x, std::vector< int > &componentsInside) const
 Transformation from global to local. More...
 
Dune::FieldMatrix< ctype, dim, dim > jacobianInverseTransposed (Dune::FieldVector< ctype, dim > const &xLocal) const
 inverse of transposed Jacobian. More...
 

Protected Attributes

Cell child
 
Cell ancestor
 

Member Typedef Documentation

◆ Entity

template<class Grid >
typedef Grid::template Codim<0>::Entity Kaskade::MultiLevelLocalGeometry< Grid >::Entity

The entity type between which this geometry object maps.

Definition at line 74 of file mllgeometry.hh.

Member Enumeration Documentation

◆ Direction

template<class Grid >
enum Kaskade::MultiLevelLocalGeometry::Direction
Enumerator
ChildIsGlobal 
FatherIsGlobal 

Definition at line 68 of file mllgeometry.hh.

Constructor & Destructor Documentation

◆ MultiLevelLocalGeometry()

template<class Grid >
Kaskade::MultiLevelLocalGeometry< Grid >::MultiLevelLocalGeometry ( Grid const &  grid,
Cell const &  child_,
Cell const &  ancestor_,
MultiLevelLocalGeometry< Grid >::Direction  direction_,
bool  geometricallyNested_ = (dim==dimw),
ctype  tol_ = -100 
)
inline

Constructor.

Parameters
gridGrid in which everthing takes place
child_Pointer to the child entity
ancestor_Pointer to the ancestor entity
direction_determine, which of both entities is seen as the global entitiy
geometricallyNested_true: assume that father and child are geometrically nested. This allows sometimes more efficient transformation, but is not alwasy feasible to assume this (for example in boundary adaptation).
tol_geometric tolerance for checking points to be inside cells or not. Use this to make sure that points inside are detected as such, even if by rounding errors they appear to be outside (by some negligible margin). Defaults to \( 100 \epsilon \).
Todo:
: (i) Simplify code by using entity pointer's comparison operator instead of ids (caveat: can ids lead to shorter hierarchies? Remember that entities copied to higher grid levels have the same id. If the comparison operators between entity pointers evaluate to "not equal" in this case, comparing ids may be more efficient. (ii) Check whether spanwidth can be computed as difference of the entity levels instead of using the loop. The same consideration as for (i) applies.

Definition at line 110 of file mllgeometry.hh.

Member Function Documentation

◆ global()

template<class Grid >
void Kaskade::MultiLevelLocalGeometry< Grid >::global ( std::vector< Dune::FieldVector< ctype, dim > > &  x,
std::vector< int > &  componentsInside 
) const
inline

Transformation from local to global.

Definition at line 138 of file mllgeometry.hh.

Referenced by Kaskade::LocalGeometryInCoarseGridAncestor< Grid >::global().

◆ jacobianInverseTransposed()

template<class Grid >
Dune::FieldMatrix< ctype, dim, dim > Kaskade::MultiLevelLocalGeometry< Grid >::jacobianInverseTransposed ( Dune::FieldVector< ctype, dim > const &  xLocal) const
inline

inverse of transposed Jacobian.

The Jacobian is just the derivative of the local to global map \( g \). This methhod returns the transpose of its inverse.

Definition at line 188 of file mllgeometry.hh.

◆ local()

template<class Grid >
void Kaskade::MultiLevelLocalGeometry< Grid >::local ( std::vector< Dune::FieldVector< ctype, dim > > &  x,
std::vector< int > &  componentsInside 
) const
inline

Transformation from global to local.

Only those points are returned that lie inside the local entity. All others are deleted during the process.

Parameters
xa sequence of global positions. On return, only those are retained and transformed to local positions that are inside the local cell.
componentsInsidea sequence of the same length as x. Only the elements corresponding to retained positions are retained.

Usually, componentsInside contains an index list 0,1,...,size(x)-1 on entry. This way, on exit the vector contains the original indices of the retained positions.

Note that "insideness" can be controlled by setting the tol_ parameter on construction.

Definition at line 169 of file mllgeometry.hh.

Referenced by Kaskade::localProlongationMatrix(), and Kaskade::LocalTransfer< Space, CoarseningPolicy >::LocalTransfer().

◆ localEntity()

template<class Grid >
Entity const & Kaskade::MultiLevelLocalGeometry< Grid >::localEntity ( ) const
inline

A reference to the local entity.

Definition at line 134 of file mllgeometry.hh.

Member Data Documentation

◆ ancestor

template<class Grid >
Cell Kaskade::MultiLevelLocalGeometry< Grid >::ancestor
protected

◆ child

template<class Grid >
Cell Kaskade::MultiLevelLocalGeometry< Grid >::child
protected

The documentation for this class was generated from the following file: