KASKADE 7 development version
Classes | Typedefs | Enumerations | Functions | Variables

Tools for managing grids and their refinement, and for iterating over grids. More...

Classes

class  Kaskade::BoundaryEdge< Index, Vector >
 A class storing a boundary edge consisting of two vertex indices and associated boundary curve tangents. More...
 
struct  Kaskade::Corner< Index, ctype, dim >
 A class storing algebraic and geometric information about a vertex-boundaryface incidence in simplicial grids. More...
 
struct  Kaskade::BoundaryStar< Index, ctype, dim >
 A class storing all the corners (vertex-boundaryface incidences) at a boundary vertex. More...
 
class  Kaskade::GeometricFeatureDetector< Index, ctype, dim >
 The interface for specifying feature edges and feature vertices. More...
 
class  Kaskade::AngleFeatureDetector< Index, ctype, dim >
 Specifying geometric features in terms of angles between face normals and edge tangents. More...
 
class  Kaskade::BoundaryTangents< dim, Index, ctype >
 A class providing information about domain boundary tangent vectors. More...
 
class  Kaskade::BoundaryInterpolationDisplacement< GridView >
 The BoundaryInterpolationDisplacement class provides a Kaskade WeakFunctionView on the displacement from grid to actual domain boundary (which is computed by boundary interpolation). More...
 
class  Kaskade::DeformingGridManager< Grid >
 A grid manager for asymptotically \( G^1 \)-continuous domains. More...
 
class  Kaskade::GlobalBarycentricPermutation< dimension >
 A class for computing permutations of local vertex numbers of simplex subentities to a globally unique ordering. More...
 
struct  Kaskade::GridSignals
 A class that provides access to signals that are emitted from the grid manager on various occasions. More...
 
class  Kaskade::CellRanges< GridView >
 A class storing cell ranges of roughly equal size for multithreaded mesh traversal. More...
 
class  Kaskade::GridManagerBase< Grd >
 Basic functionality for managing grids and their refinement. More...
 
class  Kaskade::GridManager< Grd >
 Class that takes care of the transfer of data (of each FunctionSpaceElement) after modification of the grid. More...
 

Typedefs

template<class GridView >
using GridBoundaryTangents = BoundaryTangents< GridView::dimension, typename GridView::IndexSet::IndexType, typename GridView::ctype >
 A convenient type alias for BoundaryTangents if the grid view type is known. More...
 
template<class GridView >
using Kaskade::Cell = typename GridView::template Codim< 0 >::Entity
 The type of cells (entities of codimension 0) in the grid view. More...
 
template<class GridView >
using Kaskade::Face = typename GridView::Intersection
 The type of faces in the grid view. More...
 
template<class GridView >
using Kaskade::GlobalPosition = Dune::FieldVector< typename GridView::ctype, GridView::dimensionworld >
 The type of global positions within the grid view. More...
 
template<class GridView >
using Kaskade::LocalPosition = Dune::FieldVector< typename GridView::ctype, GridView::dimension >
 The type of local positions within cells of the grid view. More...
 
template<class GridView >
using Kaskade::FaceLocalPosition = Dune::FieldVector< typename GridView::ctype, GridView::dimension-1 >
 The type of local positions within faces of the grid view. More...
 

Enumerations

enum  { Kaskade::GridSignals::freeResources =0 , Kaskade::GridSignals::allocResources =1 }
 
enum  Kaskade::GridSignals::Status { Kaskade::GridSignals::BeforeRefinement , Kaskade::GridSignals::BeforeRefinementComplete , Kaskade::GridSignals::AfterRefinement , Kaskade::GridSignals::TransferCompleted }
 The argument type of the signal that is emitted before and after grid adaptation. More...
 

Functions

template<class CoordType , int dim>
Dune::FieldVector< CoordType, dim+1 > Kaskade::barycentric (Dune::FieldVector< CoordType, dim > const &x)
 Computes the barycentric coordinates of a point in the unit simplex. More...
 
template<class CoordType , int dim>
Dune::FieldVector< CoordType, dim+1 > Kaskade::barycentric2 (Dune::FieldVector< CoordType, dim > const &x)
 Computes barycentric coordinates of a point in the unit simplex. More...
 
template<class CoordType , int dim>
Dune::FieldMatrix< CoordType, dim+1, dim > Kaskade::barycentric2Derivative (Dune::FieldVector< CoordType, dim > const &x)
 Computes the derivative of the barycentric coordinates of a point in the unit simplex. More...
 
template<class GridView >
std::vector< Dune::FieldVector< typename GridView::ctype, GridView::dimension > > Kaskade::interpolateVertexPositions (GridView const &gv, GridBoundaryTangents< GridView > const &tangents, bool computeDisplacements=false)
 Computes a new position for each grid vertex based on boundary interpolation. More...
 
template<class GridView , class Functor >
void Kaskade::forEachBoundaryFace (GridView const &gridView, Functor functor)
 iterates over each boundary face and applies functor to face. Each boundary face is visited exactly once. More...
 
template<class Cell , class ctype >
std::pair< Cell, Dune::FieldVector< ctype, Cell::dimension > > Kaskade::findCellOnLevel (Cell cell, Dune::FieldVector< ctype, Cell::dimension > xi, int level=std::numeric_limits< int >::max())
 Returns a descendant or ancestor cell of the given cell containing the specified local coordinate. More...
 
template<class GridView >
Cell< GridView > Kaskade::findCell (GridView const &gv, GlobalPosition< GridView > global)
 Returns a cell of the given grid containing the specified global coordinate. More...
 
template<class GridView >
std::pair< Dune::FieldVector< typename GridView::ctype, GridView::dimensionworld >, Dune::FieldVector< typename GridView::ctype, GridView::dimensionworld > > Kaskade::boundingBox (GridView const &gv)
 Computes a the bounding box of the given grid view. More...
 
template<class GridView >
std::unique_ptr< typename GridView::Grid > Kaskade::createCoverGrid (GridView const &gv, double volumeRatio=10)
 Creates a cartesian structured grid occupying the bounding box of the given grid view. More...
 
template<class Grid >
void Kaskade::refineAtBoundary (GridManagerBase< Grid > &gridManager)
 A convenience routine that refines all cells incident to a boundary face. More...
 
template<class VolumeDisplacement >
auto makeCompositeBoundaryInterpolation (VolumeDisplacement const &u, GeometricFeatureDetector< typename VolumeDisplacement::Space::Grid::LevelGridView::IndexSet::IndexType, typename VolumeDisplacement::Space::Grid::ctype, VolumeDisplacement::Space::Grid::dimension > const &features)
 Creates a composite boundary interpolation displacement from a volume displacement. More...
 

Variables

boost::signals2::signal< void(Status)> Kaskade::GridSignals::informAboutRefinement
 A signal that is emitted thrice on grid modifications, once before adaptation takes place and twice after it is completed. More...
 

Detailed Description

Tools for managing grids and their refinement, and for iterating over grids.

Typedef Documentation

◆ Cell

template<class GridView >
using Kaskade::Cell = typedef typename GridView::template Codim<0>::Entity

The type of cells (entities of codimension 0) in the grid view.

Template Parameters
GridViewthe type of grid view used (a grid type works as well)

Definition at line 35 of file gridBasics.hh.

◆ Face

template<class GridView >
using Kaskade::Face = typedef typename GridView::Intersection

The type of faces in the grid view.

Definition at line 42 of file gridBasics.hh.

◆ FaceLocalPosition

template<class GridView >
using Kaskade::FaceLocalPosition = typedef Dune::FieldVector<typename GridView::ctype,GridView::dimension-1>

The type of local positions within faces of the grid view.

Template Parameters
GridViewa grid view type (a grid type works as well)

Definition at line 67 of file gridBasics.hh.

◆ GlobalPosition

template<class GridView >
using Kaskade::GlobalPosition = typedef Dune::FieldVector<typename GridView::ctype,GridView::dimensionworld>

The type of global positions within the grid view.

Template Parameters
GridViewa grid view type (a grid type works as well)

Definition at line 50 of file gridBasics.hh.

◆ GridBoundaryTangents

template<int dim, class Index , class ctype = double>
template<class GridView >
using GridBoundaryTangents = BoundaryTangents<GridView::dimension,typename GridView::IndexSet::IndexType,typename GridView::ctype>
related

A convenient type alias for BoundaryTangents if the grid view type is known.

Definition at line 589 of file boundaryInterpolation.hh.

◆ LocalPosition

template<class GridView >
using Kaskade::LocalPosition = typedef Dune::FieldVector<typename GridView::ctype,GridView::dimension>

The type of local positions within cells of the grid view.

Template Parameters
GridViewa grid view type (a grid type works as well)

Definition at line 59 of file gridBasics.hh.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum

Signal group names for use with the informAboutRefinement signal. Freeing resources should occur before allocating new resources.

Enumerator
freeResources 
allocResources 

Definition at line 54 of file gridmanager.hh.

◆ Status

The argument type of the signal that is emitted before and after grid adaptation.

  • BeforeRefinement: start preparing everything needed for adaptation (can be asynchronously)
  • BeforeRefinementComplete: finish preparations (blocking)
  • AfterRefinement: start transfering any data to the new grid (can be asynchronously)
  • TransferCompleted: finish transfer (blocking)
Enumerator
BeforeRefinement 
BeforeRefinementComplete 
AfterRefinement 
TransferCompleted 

Definition at line 64 of file gridmanager.hh.

Function Documentation

◆ barycentric()

template<class CoordType , int dim>
Dune::FieldVector< CoordType, dim+1 > Kaskade::barycentric ( Dune::FieldVector< CoordType, dim > const &  x)

◆ barycentric2()

template<class CoordType , int dim>
Dune::FieldVector< CoordType, dim+1 > Kaskade::barycentric2 ( Dune::FieldVector< CoordType, dim > const &  x)

Computes barycentric coordinates of a point in the unit simplex.

Here, the barycentric coordinates of \( x \) are \( (1-\sum_{i=0}^d xi,x_0,\dots,x_{d-1}) \). This corresponds to vertex zero being the cartesian origin.

Points extremely close to the unit simplex boundary (less than \( 100 \epsilon \) distance) are snapped to the boundary. This is intended to prevent Dirichlet penalty factors with huge penalty factors from affecting the computation due to roundoff.

See also
barycentric2Derivative

Definition at line 90 of file barycentric.hh.

Referenced by Kaskade::bezier(), Kaskade::BoundaryInterpolationDisplacement< GridView >::derivative(), Kaskade::BoundaryInterpolationDetail::getChildVertexDisplacements(), and Kaskade::BoundaryInterpolationDisplacement< GridView >::value().

◆ barycentric2Derivative()

template<class CoordType , int dim>
Dune::FieldMatrix< CoordType, dim+1, dim > Kaskade::barycentric2Derivative ( Dune::FieldVector< CoordType, dim > const &  x)

Computes the derivative of the barycentric coordinates of a point in the unit simplex.

Here, the barycentric coordinates of \( x \) are \( (1-\sum_{i=0}^d xi,x_0,\dots,x_{d-1}) \). This corresponds to vertex zero being the cartesian origin.

See also
barycentric2

Definition at line 122 of file barycentric.hh.

Referenced by Kaskade::BoundaryInterpolationDisplacement< GridView >::derivative().

◆ boundingBox()

template<class GridView >
std::pair< Dune::FieldVector< typename GridView::ctype, GridView::dimensionworld >, Dune::FieldVector< typename GridView::ctype, GridView::dimensionworld > > Kaskade::boundingBox ( GridView const &  gv)

Computes a the bounding box of the given grid view.

Returns
a pair (low,high) of vectors denoting the smallest and largest corner of the bounding box Keep in mind that level grid views need not cover the whole domain.

Definition at line 35 of file gridCover.hh.

Referenced by Kaskade::createCoverGrid(), and Kaskade::AnsysMeshReader::readGrid().

◆ createCoverGrid()

template<class GridView >
std::unique_ptr< typename GridView::Grid > Kaskade::createCoverGrid ( GridView const &  gv,
double  volumeRatio = 10 
)

Creates a cartesian structured grid occupying the bounding box of the given grid view.

Starting from a uniform grid, the cover grid is locally refined until its cells are not larger than the provided grid view cells by a certain factor.

Parameters
gvthe provided grid view
volumeRatiothe maximum local ratio of cell volumes of the cover grid and the provided grid

Definition at line 70 of file gridCover.hh.

Referenced by Kaskade::makeAuxiliarySpaceMultigridStack().

◆ findCell()

template<class GridView >
Cell< GridView > Kaskade::findCell ( GridView const &  gv,
GlobalPosition< GridView >  global 
)

Returns a cell of the given grid containing the specified global coordinate.

This method makes use of the grid hierarchy, and is (on refined meshes) significantly faster than a linear search. If no cell containing the given point is found, a LookupException is thrown.

Uses now checkInside from fixdune.hh and a (fixed) tolerance.

Parameters
gvthe grid view on which to operate (TODO: does this make sense? We only use the grid itself!)
globalthe global position for which a cell shall be found

Definition at line 213 of file functionspace.hh.

Referenced by Kaskade::createCoverGrid(), Kaskade::FunctionSpaceElement< FunctionSpace, m >::derivative(), Kaskade::FunctionSpaceElement< FunctionSpace, m >::hessian(), Kaskade::FEFunctionSpace< LocalToGlobalMapper >::linearCombination(), Kaskade::twoGridProlongation(), and Kaskade::FunctionSpaceElement< FunctionSpace, m >::value().

◆ findCellOnLevel()

template<class Cell , class ctype >
std::pair< Cell, Dune::FieldVector< ctype, Cell::dimension > > Kaskade::findCellOnLevel ( Cell  cell,
Dune::FieldVector< ctype, Cell::dimension >  xi,
int  level = std::numeric_limits<int>::max() 
)

Returns a descendant or ancestor cell of the given cell containing the specified local coordinate.

This method makes use of the grid hierarchy, and is significantly faster than a linear search.

Parameters
cellthe start cell
xithe local position in the start cell for which a descendant cell shall be found
levelgrid level on which the target cell is searched (nonnegative)

If level <= cell.level(), we go downwards to coarser levels, until the requested level has been reached. If level > cell.level(), we go upwards until the requested level or the leaf level is reached.

Returns
a pair (cell,xi)

Postcondition: return.level()<=level && checkInside(return.type(),return.geometry().local(cell.geometry().global(xi)))

Definition at line 159 of file functionspace.hh.

Referenced by Kaskade::BoundaryInterpolationDisplacement< GridView >::applyVolumeDisplacement().

◆ forEachBoundaryFace()

template<class GridView , class Functor >
void Kaskade::forEachBoundaryFace ( GridView const &  gridView,
Functor  functor 
)

iterates over each boundary face and applies functor to face. Each boundary face is visited exactly once.

Parameters
functora callable object taking a single intersection as argument

Definition at line 88 of file forEach.hh.

Referenced by Kaskade::computeBoundaryStars(), Kaskade::integrateOverBoundary(), and Kaskade::BoundaryMapper< Implementation_, Scalar_, GridView_ >::numFaces().

◆ interpolateVertexPositions()

template<class GridView >
std::vector< Dune::FieldVector< typename GridView::ctype, GridView::dimension > > Kaskade::interpolateVertexPositions ( GridView const &  gv,
GridBoundaryTangents< GridView > const &  tangents,
bool  computeDisplacements = false 
)

Computes a new position for each grid vertex based on boundary interpolation.

The grid has to contain only simplicial cells. If a cell of the grid has more than one boundary face, these need to be separated by a feature edge - otherwise the grid deformation will necessarily create degenerated flat cells on refinement. Cells without boundary face shall have at most one boundary edge (only possible in 3D).

Parameters
gvthe grid view for which the vertex positions are computes, usually a leaf grid view
tangentsprovides the surface tangent vectors corresponding to boundary edges
computeDisplacementsif true, only the vertex displacements are computed, otherwise their new positions will be returned
Returns
a vector container of the displaced vertex positions (or the displacements, depending on the computeDisplacements parameter)

Definition at line 749 of file boundaryInterpolation.hh.

◆ makeCompositeBoundaryInterpolation()

template<class VolumeDisplacement >
auto makeCompositeBoundaryInterpolation ( VolumeDisplacement const &  u,
GeometricFeatureDetector< typename VolumeDisplacement::Space::Grid::LevelGridView::IndexSet::IndexType, typename VolumeDisplacement::Space::Grid::ctype, VolumeDisplacement::Space::Grid::dimension > const &  features 
)
related

Creates a composite boundary interpolation displacement from a volume displacement.

Definition at line 984 of file boundaryInterpolation.hh.

◆ refineAtBoundary()

template<class Grid >
void Kaskade::refineAtBoundary ( GridManagerBase< Grid > &  gridManager)

A convenience routine that refines all cells incident to a boundary face.

This is useful to create a priori boundary-refined meshes.

Definition at line 736 of file gridmanager.hh.

Variable Documentation

◆ informAboutRefinement

boost::signals2::signal<void (Status)> Kaskade::GridSignals::informAboutRefinement

A signal that is emitted thrice on grid modifications, once before adaptation takes place and twice after it is completed.

The current state of grid refinement and FE coefficient transfer is signaled by the status argument.

Definition at line 74 of file gridmanager.hh.

Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::AdaptiveGrid(), Kaskade::GridManagerBase< Grd >::afterRefinement(), Kaskade::GridManagerBase< Grd >::beforeRefinement(), Kaskade::CachingBoundaryDetector< GridView >::CachingBoundaryDetector(), Kaskade::FEFunctionSpace< LocalToGlobalMapper >::FEFunctionSpace(), HierarchicIndexSet< Grid >::HierarchicIndexSet(), and HierarchicIndexSet< Grid >::operator=().