dune-mmesh (1.4)
Provides a DUNE grid interface class for the interface of a MMesh interface grid. More...
#include <dune/mmesh/interface/grid.hh>
Public Types | |
typedef MMeshInterfaceGridFamily< MMesh > | GridFamily |
the Traits | |
using | GridImp = typename GridFamily::Traits::Grid |
the grid implementation | |
typedef std::unique_ptr< GridImp > | GridPtrType |
the unique pointer to the grid | |
using | HostGridType = typename MMesh::HostGridType |
the underlying hostgrid | |
using | MMeshType = MMesh |
the underlying mmesh | |
typedef Traits::template Codim< 0 >::LeafIterator | LeafIterator |
the leaf iterator | |
typedef FieldType | ctype |
The type used to store coordinates, inherited from the MMesh. | |
typedef Dune::FieldVector< ctype, dimensionworld > | GlobalCoordinate |
The type used for coordinates. | |
template<int cd> | |
using | MMeshInterfaceEntity = typename HostGridEntityChooser_< HostGridType, dimensionworld, cd+1 >::type |
The type of the underlying entities. | |
using | VertexHandle = MMeshInterfaceEntity< dimension > |
The type of the underlying vertex handle. | |
using | EdgeHandle = MMeshInterfaceEntity< dimension-1 > |
The type of the underlying edge handle. | |
using | ElementOutput = std::list< MMeshInterfaceEntity< 0 > > |
The type of the element output. | |
using | Entity = typename Traits::template Codim< 0 >::Entity |
The type of a codim 0 entity. | |
using | Vertex = typename Traits::template Codim< dimension >::Entity |
The type of a vertex. | |
using | ConnectedComponent = MMeshInterfaceConnectedComponent< const GridImp > |
The type of a connected component. | |
using | RemeshingIndicator = RatioIndicator< GridImp > |
The type of the default remeshing indicator. | |
Public Member Functions | |
MMeshInterfaceGrid (MMesh *mMesh, BoundarySegments boundarySegments={}) | |
Constructor. | |
int | maxLevel () const |
Return maximum level defined in this grid. More... | |
int | size (int level, int codim) const |
Number of grid entities per level and codim. | |
size_t | numBoundarySegments () const |
returns the number of boundary segments within the macro grid | |
int | size (int codim) const |
number of leaf entities per codim in this process | |
int | size (int level, GeometryType type) const |
number of entities per level, codim and geometry type in this process | |
int | size (GeometryType type) const |
number of leaf entities per codim and geometry type in this process | |
const Traits::GlobalIdSet & | globalIdSet () const |
Access to the GlobalIdSet. | |
const Traits::LocalIdSet & | localIdSet () const |
Access to the LocalIdSet. | |
const MMeshInterfaceGridLeafIndexSet< const GridImp > & | levelIndexSet (int level) const |
Access to the LevelIndexSets. | |
const MMeshInterfaceGridLeafIndexSet< const GridImp > & | leafIndexSet () const |
Access to the LeafIndexSet. | |
template<class EntitySeed > | |
Traits::template Codim< EntitySeed::codimension >::Entity | entity (const EntitySeed &seed) const |
Create Entity from EntitySeed. | |
Traits::template Codim< dimension >::Entity | entity (const MMeshInterfaceEntity< dimension > &hostEntity) const |
Create Entity from a host entity. | |
Traits::template Codim< 0 >::Entity | entity (const MMeshInterfaceEntity< 0 > &hostEntity) const |
Create codim 0 entity from a host entity. | |
template<int dimw = dimension+1> | |
std::enable_if_t< dimw==3, typename Traits::template Codim< 1 >::Entity > | entity (const MMeshInterfaceEntity< 1 > &hostEntity) const |
Create codim 0 entity from a host entity. | |
template<int codim> | |
Traits::template Codim< codim >::LevelIterator | lbegin (int level) const |
Iterator to first entity of given codim on level. | |
template<int codim> | |
Traits::template Codim< codim >::LevelIterator | lend (int level) const |
one past the end on this level | |
template<int codim, PartitionIteratorType PiType> | |
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator | lbegin (int level) const |
Iterator to first entity of given codim on level. | |
template<int codim, PartitionIteratorType PiType> | |
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator | lend (int level) const |
one past the end on this level | |
template<int codim> | |
Traits::template Codim< codim >::LeafIterator | leafbegin () const |
Iterator to first leaf entity of given codim. | |
template<int codim> | |
Traits::template Codim< codim >::LeafIterator | leafend () const |
one past the end of the sequence of leaf entities | |
template<int codim, PartitionIteratorType PiType> | |
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator | leafbegin () const |
Iterator to first leaf entity of given codim. | |
template<int codim, PartitionIteratorType PiType> | |
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator | leafend () const |
one past the end of the sequence of leaf entities | |
void | globalRefine (int steps=1) |
Global refine. More... | |
bool | mark (int refCount, const typename Traits::template Codim< 0 >::Entity &e) const |
Mark entity for refinement. More... | |
int | getMark (const typename Traits::template Codim< 0 >::Entity &e) const |
Return refinement mark for entity. More... | |
bool | preAdapt () |
returns false, if at least one entity is marked for adaption | |
bool | adapt () |
Triggers the grid adaptation process. More... | |
template<class GridImp , class DataHandle > | |
bool | adapt (AdaptDataHandleInterface< GridImp, DataHandle > &handle) |
Callback for the grid adaptation process with restrict/prolong. | |
void | postAdapt () |
Clean up refinement markers. | |
unsigned int | overlapSize (int codim) const |
Size of the overlap on the leaf level. | |
unsigned int | ghostSize (int codim) const |
Size of the ghost cell layer on the leaf level. | |
unsigned int | overlapSize (int level, int codim) const |
Size of the overlap on a given level. | |
unsigned int | ghostSize (int level, int codim) const |
Size of the ghost cell layer on a given level. | |
void | loadBalance () |
Distributes this grid over the available nodes in a distributed machine. More... | |
const Communication< Comm > & | comm () const |
Collective communication. | |
bool | isInterface (const MMeshInterfaceEntity< 0 > &segment) const |
Return if interface segment is part of the interface. | |
template<int d = dimension> | |
std::enable_if_t< d==2, bool > | isInterface (const MMeshInterfaceEntity< 1 > &edge) const |
Return if an edge is of the interface. | |
bool | isInterface (const VertexHandle &vertex) const |
Return if vertex is part of the interface. | |
bool | isInterface (const typename Traits::LeafIntersection &intersection) const |
Return if intersection is part of the interface. | |
template<class Entity > | |
bool | isInterface (const Entity &entity) const |
Return if dune entity is part of the interface. | |
template<class Entity > | |
std::size_t | domainMarker (const Entity &entity) const |
Return domain marker of entity. | |
void | markAsRefined (const std::vector< std::vector< std::size_t > > &children, const ConnectedComponent connectedComponent) |
Mark a set of children elements as refinement of a connected component. | |
bool | hasConnectedComponent (const Entity &entity) const |
Return if an entity has a connected component. | |
const ConnectedComponent & | getConnectedComponent (const Entity &entity) const |
Return the connected component for a given entity. | |
ConnectedComponent & | getConnectedComponent (const Entity &entity) |
Return a non-const reference to the connected component for a given entity. | |
template<int d = dimensionworld> | |
std::enable_if_t< d==2, const MMeshInterfaceEntity< 0 > > | mirrorHostEntity (const MMeshInterfaceEntity< 0 > &other) const |
Mirror a host entity (2d) | |
template<int d = dimensionworld> | |
std::enable_if_t< d==3, const MMeshInterfaceEntity< 0 > > | mirrorHostEntity (const MMeshInterfaceEntity< 0 > &other) const |
Mirror a host entity (3d) | |
void | setIds () |
compute the grid ids | |
void | setIndices () |
compute the grid indices | |
const MMesh & | getMMesh () const |
Return reference to MMesh. | |
const HostGridType & | getHostGrid () const |
Return reference to underlying CGAL triangulation. | |
int | sequence () const |
Return sequence number. | |
Protected Attributes | |
MMesh * | mMesh_ |
The host grid which contains the actual grid hierarchy structure. | |
Detailed Description
class Dune::MMeshInterfaceGrid< MMesh >
Provides a DUNE grid interface class for the interface of a MMesh interface grid.
Member Function Documentation
◆ adapt()
|
inline |
Triggers the grid adaptation process.
- Returns
- if triangulation has changed
References Dune::MMesh< HostGrid, dim >::adapt(), and Dune::MMeshInterfaceGrid< MMesh >::mMesh_.
Referenced by Dune::MMeshInterfaceGrid< MMesh >::globalRefine().
◆ getMark()
|
inline |
Return refinement mark for entity.
- Returns
- refinement mark
Referenced by Dune::MMesh< HostGrid, dim >::adapt().
◆ globalRefine()
|
inline |
Global refine.
Marks all elements for refinement and adapts the grid
References Dune::MMeshInterfaceGrid< MMesh >::adapt(), Dune::MMeshInterfaceGrid< MMesh >::mark(), Dune::MMeshInterfaceGrid< MMesh >::postAdapt(), and Dune::MMeshInterfaceGrid< MMesh >::preAdapt().
◆ loadBalance()
|
inline |
Distributes this grid over the available nodes in a distributed machine.
- Parameters
-
minlevel The coarsest grid level that gets distributed maxlevel does currently get ignored
References Dune::MMesh< HostGrid, dim >::loadBalance(), and Dune::MMeshInterfaceGrid< MMesh >::mMesh_.
◆ mark()
|
inline |
Mark entity for refinement.
This only works for entities of codim 0.
Referenced by Dune::MMeshInterfaceGrid< MMesh >::globalRefine().
◆ maxLevel()
|
inline |
Return maximum level defined in this grid.
Levels are numbered 0 ... maxlevel with 0 the coarsest level.
The documentation for this class was generated from the following file:
- dune/mmesh/interface/grid.hh
