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

A grid manager for asymptotically \( G^1 \)-continuous domains. More...

#include <deforminggridmanager.hh>

Detailed Description

template<class Grid>
class Kaskade::DeformingGridManager< Grid >

A grid manager for asymptotically \( G^1 \)-continuous domains.

On construction, a deformation of the grid domain is computed, such that the deformed domain has a \( G^1 \)-continuous boundary. On mesh refinement, the new vertices are moved such that the smooth boundary is interpolated.

Definition at line 117 of file deforminggridmanager.hh.

Inheritance diagram for Kaskade::DeformingGridManager< Grid >:
Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >

Public Types

using Grid = Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > >
 the type of managed grid More...
 
typedef GridManagerBase< GridSelf
 

Public Member Functions

 DeformingGridManager (std::unique_ptr< Grid > &&grid, ctype featureCurveAngle=1.0, ctype featureEdgeAngle=0.5, bool verbose=false, bool enforceConcurrentReads=false)
 Constructor. More...
 
 DeformingGridManager (std::unique_ptr< Grid > &&grid, GeometricFeatureDetector< Index, ctype, dim > const &featureDetector, bool verbose=false, bool enforceConcurrentReads=false)
 
virtual void update ()
 Must be overloaded by derived classes to adjust internal state immediately after mesh refinement. More...
 
virtual void refineGrid (int refCount)
 
virtual bool adaptGrid ()
 
template<class Scalar >
void correctingDisplacement (FunctionSpaceElement< H1Space< DeformedGrid, Scalar >, dim > &correction)
 correctingDisplacement provides a displacement function which displaces the points on the boundary of the grid on to the actual domain boundary. More...
 
SelfenforceConcurrentReads (bool enforceConcurrentReads)
 Tells the grid manager that concurrent reads of the grid are safe. More...
 
SelfsetVerbosity (const bool verbosity)
 sets the verbosity level of the grid manager More...
 
bool gridIsThreadSafe () const
 Returns true if concurrent read accesses to the grid do not lead to data races. More...
 
CellRanges< typename Grid::LeafGridView > const & cellRanges (typename Grid::LeafGridView const &) const
 Returns a CellRanges object for the given grid view. The cell ranges are created on demand. The reference is valid up to the next mesh modification. This method is thread-safe. More...
 
CellRanges< typename Grid::LevelGridView > const & cellRanges (typename Grid::LevelGridView const &gridView) const
 Returns a CellRanges object for the given grid view. More...
 
EntityNumbering< typename Grid::LeafGridView > const & entityNumbering (typename Grid::LeafGridView const &) const
 Returns an EntityNumbering object for the given grid view. More...
 
Access to the managed grid.
Grid const & grid () const
 Returns a const reference on the owned grid. More...
 
Gridgrid_non_const ()
 Returns a non-const reference on the owned grid. More...
 
std::shared_ptr< Grid const > gridShared () const
 Provides shared ownership management. More...
 
Grid adaptation
bool mark (int refCount, typename Grid::template Codim< 0 >::Entity const &cell)
 Marks the given cell for refinement or coarsening. More...
 
void mark (CellData< Grid, S > const &cellData)
 Marks all Entities of a grid according to the corresponding entries in a CellData. More...
 
size_t countMarked () const
 Reports the number of marked elements. More...
 
bool preAdapt ()
 Calls grid.preAdapt(). More...
 
void postAdapt ()
 Calls grid.postAdapt(). More...
 
bool adaptAtOnce ()
 DEPRECATED Performs grid refinement without prolongating registered FE functions. More...
 
void flushMarks ()
 Remove all cell marks from the grid. More...
 
bool adapt ()
 Refines each marked leaf cell of the grid and prolongates registered FE functions to the new leaf grid view. More...
 
void globalRefine (int refCount)
 Refines each leaf cell of the grid and prolongates registered FE functions to the new leaf grid view. More...
 

Public Attributes

GridSignals signals
 

Protected Member Functions

void beforeRefinement ()
 
void afterRefinement ()
 

Protected Attributes

std::shared_ptr< Gridgridptr
 The grid itself. More...
 
bool verbose
 

Member Typedef Documentation

◆ Grid

using Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::Grid = Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > >
inherited

the type of managed grid

Definition at line 266 of file gridmanager.hh.

◆ Self

typedef GridManagerBase<Grid> Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::Self
inherited

Definition at line 268 of file gridmanager.hh.

Constructor & Destructor Documentation

◆ DeformingGridManager() [1/2]

template<class Grid >
Kaskade::DeformingGridManager< Grid >::DeformingGridManager ( std::unique_ptr< Grid > &&  grid,
ctype  featureCurveAngle = 1.0,
ctype  featureEdgeAngle = 0.5,
bool  verbose = false,
bool  enforceConcurrentReads = false 
)
inline

Constructor.

See BoundaryTangents for a detailed description of feature edge definition.

Parameters
featureCurveAnglethe angle (in radians) between two feature edges that determines whether they form a differentiable feature curve. Defaults to roughly 57°.
featureEdgeAnglethe angle (in radians) between adjacent facets' normals that determines whether they share a feature edge or a normal ( \( G^1 \) continuous) edge. Defaults to roughly 29°.

Definition at line 141 of file deforminggridmanager.hh.

◆ DeformingGridManager() [2/2]

template<class Grid >
Kaskade::DeformingGridManager< Grid >::DeformingGridManager ( std::unique_ptr< Grid > &&  grid,
GeometricFeatureDetector< Index, ctype, dim > const &  featureDetector,
bool  verbose = false,
bool  enforceConcurrentReads = false 
)
inline

Definition at line 150 of file deforminggridmanager.hh.

Member Function Documentation

◆ adapt()

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::adapt ( )
inlineinherited

Refines each marked leaf cell of the grid and prolongates registered FE functions to the new leaf grid view.

Definition at line 413 of file gridmanager.hh.

◆ adaptAtOnce()

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::adaptAtOnce ( )
inlineinherited

DEPRECATED Performs grid refinement without prolongating registered FE functions.

As this can easily lead to inconsistent states of FE functions, usage is not recommended. Probably this method will be removed soon.

Definition at line 390 of file gridmanager.hh.

◆ adaptGrid()

template<class Grid >
virtual bool Kaskade::DeformingGridManager< Grid >::adaptGrid ( )
inlinevirtual

◆ afterRefinement()

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::afterRefinement ( )
inlineprotectedinherited

Definition at line 595 of file gridmanager.hh.

◆ beforeRefinement()

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::beforeRefinement ( )
inlineprotectedinherited

Definition at line 583 of file gridmanager.hh.

◆ cellRanges() [1/2]

CellRanges< typename Grid::LeafGridView > const & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::cellRanges ( typename Grid::LeafGridView const &  ) const
inlineinherited

Returns a CellRanges object for the given grid view. The cell ranges are created on demand. The reference is valid up to the next mesh modification. This method is thread-safe.

Definition at line 531 of file gridmanager.hh.

◆ cellRanges() [2/2]

CellRanges< typename Grid::LevelGridView > const & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::cellRanges ( typename Grid::LevelGridView const &  gridView) const
inlineinherited

Returns a CellRanges object for the given grid view.

Parameters
gridViewa level grid view of the grid managed by this grid manager

The cell ranges are created on demand. The reference is valid up to the next mesh modification. This method is thread-safe.

Definition at line 547 of file gridmanager.hh.

◆ correctingDisplacement()

template<class Grid >
template<class Scalar >
void Kaskade::DeformingGridManager< Grid >::correctingDisplacement ( FunctionSpaceElement< H1Space< DeformedGrid, Scalar >, dim > &  correction)
inline

correctingDisplacement provides a displacement function which displaces the points on the boundary of the grid on to the actual domain boundary.

This function is zero on all grid vertices, so it is only useful if requested for FE order greater than one. We leave it undefined what values the provided correction has in the interior of the grid.

Parameters
correctionFE function of any order >1, will be the computed boundary correction on exit.

Definition at line 181 of file deforminggridmanager.hh.

◆ countMarked()

size_t Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::countMarked ( ) const
inlineinherited

Reports the number of marked elements.

Definition at line 359 of file gridmanager.hh.

◆ enforceConcurrentReads()

Self & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::enforceConcurrentReads ( bool  enforceConcurrentReads)
inlineinherited

Tells the grid manager that concurrent reads of the grid are safe.

Definition at line 484 of file gridmanager.hh.

◆ entityNumbering()

EntityNumbering< typename Grid::LeafGridView > const & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::entityNumbering ( typename Grid::LeafGridView const &  ) const
inlineinherited

Returns an EntityNumbering object for the given grid view.

Todo:
allow to specify the sorting

Definition at line 563 of file gridmanager.hh.

◆ flushMarks()

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::flushMarks ( )
inlineinherited

Remove all cell marks from the grid.

Definition at line 402 of file gridmanager.hh.

◆ globalRefine()

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::globalRefine ( int  refCount)
inlineinherited

Refines each leaf cell of the grid and prolongates registered FE functions to the new leaf grid view.

Definition at line 449 of file gridmanager.hh.

◆ grid()

Grid const & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::grid ( ) const
inlineinherited

Returns a const reference on the owned grid.

Definition at line 292 of file gridmanager.hh.

◆ grid_non_const()

Grid & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::grid_non_const ( )
inlineinherited

Returns a non-const reference on the owned grid.

Definition at line 298 of file gridmanager.hh.

◆ gridIsThreadSafe()

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::gridIsThreadSafe ( ) const
inlineinherited

Returns true if concurrent read accesses to the grid do not lead to data races.

This is either the case if the grid implementation itself reports that it's thread-safe, or if the user has explicitly stated that concurrent reads are safe (using enforceConcurrentReads) - even if they actually are not.

The thread-safety reported by this method is used e.g. by the assembler to decide how many threads to use during assembly.

If the grid access may not be thread-safe, consider using GridLocking.

Definition at line 516 of file gridmanager.hh.

◆ gridShared()

std::shared_ptr< Grid const > Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::gridShared ( ) const
inlineinherited

Provides shared ownership management.

Note that only an immutable grid is provided, such that the responsibility for grid modifications remains with the GridManager object. This is necessary for the GridManager to trigger prolongation/interpolation on grid modification.

Definition at line 310 of file gridmanager.hh.

◆ mark() [1/2]

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::mark ( CellData< Grid, S > const &  cellData)
inlineinherited

Marks all Entities of a grid according to the corresponding entries in a CellData.

Definition at line 342 of file gridmanager.hh.

◆ mark() [2/2]

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::mark ( int  refCount,
typename Grid::template Codim< 0 >::Entity const &  cell 
)
inlineinherited

Marks the given cell for refinement or coarsening.

Parameters
refCountif positive, marks for (multiple) refinement; if negative, for (multiple) coarsening. Not every grid implementation accepts multiple refinement/coarsening.
cellthe affected cell

Definition at line 330 of file gridmanager.hh.

◆ postAdapt()

void Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::postAdapt ( )
inlineinherited

Calls grid.postAdapt().

Definition at line 374 of file gridmanager.hh.

◆ preAdapt()

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::preAdapt ( )
inlineinherited

Calls grid.preAdapt().

Definition at line 365 of file gridmanager.hh.

◆ refineGrid()

template<class Grid >
virtual void Kaskade::DeformingGridManager< Grid >::refineGrid ( int  refCount)
inlinevirtual

◆ setVerbosity()

Self & Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::setVerbosity ( const bool  verbosity)
inlineinherited

sets the verbosity level of the grid manager

The grid manager prints status and progress information to std::cout. The amount of information depends on the verbosity:

  • 0: no messages are printed
  • 1: mesh size change on refinement is reported

Definition at line 499 of file gridmanager.hh.

◆ update()

template<class Grid >
virtual void Kaskade::DeformingGridManager< Grid >::update ( )
inlinevirtual

Must be overloaded by derived classes to adjust internal state immediately after mesh refinement.

Implements Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >.

Definition at line 158 of file deforminggridmanager.hh.

Member Data Documentation

◆ gridptr

std::shared_ptr<Grid> Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::gridptr
protectedinherited

The grid itself.

Definition at line 608 of file gridmanager.hh.

◆ signals

GridSignals Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::signals
inherited

Definition at line 611 of file gridmanager.hh.

◆ verbose

bool Kaskade::GridManagerBase< Dune::GeometryGrid< Grid, DeformingGridManagerDetail::Deformation< Grid > > >::verbose
protectedinherited

Definition at line 626 of file gridmanager.hh.


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