KASKADE 7 development version
|
Basic functionality for managing grids and their refinement. More...
#include <gridmanager.hh>
Basic functionality for managing grids and their refinement.
Definition at line 260 of file gridmanager.hh.
Public Types | |
using | Grid = Grd |
the type of managed grid More... | |
typedef GridManagerBase< Grid > | Self |
Public Member Functions | |
GridManagerBase (Dune::GridPtr< Grid > grid_, bool verbose_, bool enforceConcurrentReads=false) | |
GridManagerBase (std::unique_ptr< Grid > &&grid_, bool verbose_, bool enforceConcurrentReads=false) | |
GridManagerBase (Grid *&&grid_, bool verbose_, bool enforceConcurrentReads=false) | |
virtual | ~GridManagerBase () |
Self & | enforceConcurrentReads (bool enforceConcurrentReads) |
Tells the grid manager that concurrent reads of the grid are safe. More... | |
Self & | setVerbosity (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... | |
Grid & | grid_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... | |
template<class S > | |
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... | |
virtual void | update ()=0 |
Must be overloaded by derived classes to adjust internal state immediately after mesh refinement. More... | |
Public Attributes | |
GridSignals | signals |
Protected Member Functions | |
void | beforeRefinement () |
void | afterRefinement () |
Protected Attributes | |
std::shared_ptr< Grid > | gridptr |
The grid itself. More... | |
bool | verbose |
Friends | |
template<class LocalToGlobalMapper > | |
class | FEFunctionSpace |
template<class GridMan , class ErrEst > | |
class | AdaptiveGrid |
using Kaskade::GridManagerBase< Grd >::Grid = Grd |
the type of managed grid
Definition at line 266 of file gridmanager.hh.
typedef GridManagerBase<Grid> Kaskade::GridManagerBase< Grd >::Self |
Definition at line 268 of file gridmanager.hh.
|
inline |
Definition at line 270 of file gridmanager.hh.
|
inline |
Definition at line 274 of file gridmanager.hh.
|
inline |
Definition at line 278 of file gridmanager.hh.
|
inlinevirtual |
Definition at line 284 of file gridmanager.hh.
|
inline |
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.
Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::adapt(), and Kaskade::GridManagerBase< Grd >::adaptAtOnce().
|
inline |
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.
Referenced by Kaskade::coarsening(), Kaskade::markAndRefine(), Kaskade::refineAtBoundary(), Kaskade::Adaptivity::FixedCellFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::FixedErrorFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration2< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration< Grid >::refineGrid_impl(), and Kaskade::Limex< Eq >::step().
|
inlineprotected |
Definition at line 595 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::adapt(), and Kaskade::GridManagerBase< Grd >::globalRefine().
|
inlineprotected |
Definition at line 583 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::adapt(), and Kaskade::GridManagerBase< Grd >::globalRefine().
|
inline |
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.
Referenced by Kaskade::VariationalFunctionalAssembler< F, SparseIndex, BoundaryDetector, QuadRule >::assemble().
|
inline |
Returns a CellRanges object for the given grid view.
gridView | a 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.
|
inline |
Reports the number of marked elements.
Definition at line 359 of file gridmanager.hh.
|
inline |
Tells the grid manager that concurrent reads of the grid are safe.
Definition at line 484 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::enforceConcurrentReads().
|
inline |
Returns an EntityNumbering object for the given grid view.
Definition at line 563 of file gridmanager.hh.
|
inline |
Remove all cell marks from the grid.
Definition at line 402 of file gridmanager.hh.
Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::flushMarks(), and Kaskade::GridManagerBase< Grd >::mark().
|
inline |
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.
|
inline |
Returns a const reference on the owned grid.
Definition at line 292 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::afterRefinement(), Kaskade::AmiraMeshReader::boundaryConditionsToScalarField(), Kaskade::GridManagerBase< Grd >::cellRanges(), Kaskade::coarsening(), LossyStorage< Grid, VariableSet, Space, QuantizationPolicy >::decode(), LossyStorage< Grid, VariableSet, Space, QuantizationPolicy >::encode(), Kaskade::GridManagerBase< Grd >::entityNumbering(), LossyStorage< Grid, VariableSet, Space, QuantizationPolicy >::flush(), Kaskade::Giant< Grid, Equation, VariableSet, Spaces >::giantGbit(), Kaskade::FEFunctionSpace< LocalToGlobalMapper >::grid(), LossyStorage< Grid, VariableSet, Space, QuantizationPolicy >::LossyStorage(), Kaskade::makeAuxiliarySpaceMultigridStack(), Kaskade::markCells(), MultilevelTransfer< Space, Grid >::MultilevelTransfer(), Kaskade::NleqSolver< Grid, Equation, VariableSet, Spaces >::nleqErr(), MultilevelTransfer< Space, Grid >::out(), Kaskade::refineAtBoundary(), Kaskade::Adaptivity::FixedCellFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::FixedErrorFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration2< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration< Grid >::refineGrid_impl(), Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::size(), Kaskade::Limex< Eq >::step(), and Kaskade::CheckPoint::writeGrid().
|
inline |
Returns a non-const reference on the owned grid.
Definition at line 298 of file gridmanager.hh.
|
inline |
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.
Referenced by Kaskade::VariationalFunctionalAssembler< F, SparseIndex, BoundaryDetector, QuadRule >::assemble().
|
inline |
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.
|
inline |
Marks all Entities of a grid according to the corresponding entries in a CellData.
Definition at line 342 of file gridmanager.hh.
|
inline |
Marks the given cell for refinement or coarsening.
refCount | if positive, marks for (multiple) refinement; if negative, for (multiple) coarsening. Not every grid implementation accepts multiple refinement/coarsening. |
cell | the affected cell |
Definition at line 330 of file gridmanager.hh.
Referenced by Kaskade::coarsening(), Kaskade::GridManagerBase< Grd >::flushMarks(), Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::mark(), Kaskade::GridManagerBase< Grd >::mark(), Kaskade::markAndRefine(), Kaskade::markCells(), Kaskade::refineAtBoundary(), Kaskade::Adaptivity::FixedCellFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::FixedErrorFraction< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration2< Grid >::refineGrid_impl(), Kaskade::Adaptivity::ErrorEquilibration< Grid >::refineGrid_impl(), Kaskade::Limex< Eq >::step(), and Kaskade::uniformEmbeddedErrorEstimation().
|
inline |
Calls grid.postAdapt().
Definition at line 374 of file gridmanager.hh.
Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::adapt(), and Kaskade::GridManagerBase< Grd >::adaptAtOnce().
|
inline |
Calls grid.preAdapt().
Definition at line 365 of file gridmanager.hh.
Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::adapt(), and Kaskade::GridManagerBase< Grd >::adaptAtOnce().
|
inline |
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:
Definition at line 499 of file gridmanager.hh.
Referenced by Kaskade::Limex< Eq >::step().
|
pure virtual |
Must be overloaded by derived classes to adjust internal state immediately after mesh refinement.
Implemented in Kaskade::DeformingGridManager< Grid >, Kaskade::GridManager< Grd >, Kaskade::GridManager< Grid >, and Kaskade::GridManager< typename EvolutionEquation::AnsatzVars::Grid >.
Referenced by Kaskade::GridManagerBase< Grd >::adapt(), and Kaskade::GridManagerBase< Grd >::globalRefine().
|
friend |
Definition at line 603 of file gridmanager.hh.
Definition at line 602 of file gridmanager.hh.
|
protected |
The grid itself.
Definition at line 608 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::adapt(), Kaskade::GridManagerBase< Grd >::flushMarks(), Kaskade::GridManagerBase< Grd >::grid(), Kaskade::GridManagerBase< Grd >::grid_non_const(), Kaskade::GridManagerBase< Grd >::gridShared(), Kaskade::GridManagerBase< Grd >::mark(), Kaskade::GridManagerBase< Grd >::postAdapt(), and Kaskade::GridManagerBase< Grd >::preAdapt().
GridSignals Kaskade::GridManagerBase< Grd >::signals |
Definition at line 611 of file gridmanager.hh.
Referenced by Kaskade::Bridge::AdaptiveGrid< GridManager, Estimate >::AdaptiveGrid(), Kaskade::GridManagerBase< Grd >::afterRefinement(), Kaskade::GridManagerBase< Grd >::beforeRefinement(), Kaskade::FEFunctionSpace< LocalToGlobalMapper >::FEFunctionSpace(), HierarchicIndexSet< Grid >::HierarchicIndexSet(), and HierarchicIndexSet< Grid >::operator=().
|
protected |
Definition at line 626 of file gridmanager.hh.
Referenced by Kaskade::GridManagerBase< Grd >::adapt(), Kaskade::GridManagerBase< Grd >::postAdapt(), and Kaskade::GridManagerBase< Grd >::setVerbosity().