Dune Core Modules (2.10.0)
Front-end for the grid manager of the finite element toolbox UG3. More...
#include <dune/grid/uggrid.hh>
Public Types | |
typedef UGGridFamily< dim > | GridFamily |
type of the used GridFamily for this grid | |
typedef UG::DOUBLE | ctype |
The type used to store coordinates. | |
typedef unsigned int | Rank |
The type used for process ranks. | |
Exported types | |
typedef GridFamily::Traits::LeafGridView | LeafGridView |
type of view for leaf grid | |
typedef GridFamily::Traits::LevelGridView | LevelGridView |
type of view for level grid | |
typedef GridFamily::Traits::LeafIntersection | LeafIntersection |
A type that is a model of Dune::Intersection, an intersections of two codimension 1 of two codimension 0 entities in the leaf view. | |
typedef GridFamily::Traits::LevelIntersection | LevelIntersection |
A type that is a model of Dune::Intersection, an intersections of two codimension 1 of two codimension 0 entities in a level view. | |
typedef GridFamily::Traits::LeafIntersectionIterator | LeafIntersectionIterator |
A type that is a model of Dune::IntersectionIterator which is an iterator that allows to examine, but not to modify, the intersections of codimension 1 of an leaf element (entity of codimension 0) with other leaf elements. | |
typedef GridFamily::Traits::LevelIntersectionIterator | LevelIntersectionIterator |
A type that is a model of Dune::IntersectionIterator which is an iterator that allows to examine, but not to modify, the intersections of codimension 1 of an element (entity of codimension 0) with other elements on the same level. | |
typedef GridFamily::Traits::HierarchicIterator | HierarchicIterator |
A type that is a model of Dune::HierarchicIterator A type of iterator that allows to examine, but not to modify, entities of codimension 0 that result from refinement of an entity of codimension 0. | |
typedef GridFamily::Traits::LevelIndexSet | LevelIndexSet |
A type that is a model of Dune::IndexSet which provides a consecutive, but non persistent, numbering for entities on a grid level. | |
typedef GridFamily::Traits::LeafIndexSet | LeafIndexSet |
A type that is a model of Dune::IndexSet which provides a consecutive, but non persistent, numbering for entities in the leaf grid. | |
typedef GridFamily::Traits::GlobalIdSet | GlobalIdSet |
A type that is a model of Dune::IdSet which provides a unique and persistent numbering for all entities in the grid. The numbering is unique over all processes over which the grid is partitioned. The numbering is not necessarily consecutive. | |
typedef GridFamily::Traits::LocalIdSet | LocalIdSet |
A type that is a model of Dune::IdSet which provides a unique and persistent numbering for all entities in the grid. The numbering is only unique in a single process and it is not necessarily consecutive. | |
using | Communication = typename GridFamily::Traits::Communication |
A type that is a model of Dune::Communication. It provides a portable way for communication on the set of processes used by the grid. | |
Public Member Functions | |
UGGrid (UGCommunication comm={}) | |
Default constructor. | |
~UGGrid () noexcept(false) | |
Destructor. | |
int | maxLevel () const |
template<typename Seed > | |
Traits::template Codim< Seed::codimension >::Entity | entity (const Seed &seed) const |
Create an Entity from an EntitySeed. | |
int | size (int level, int codim) const |
Number of grid entities per level and codim. | |
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 and geometry type in this process | |
int | size (GeometryType type) const |
number of leaf entities per geometry type in this process | |
size_t | numBoundarySegments () const |
Return the number of boundary segments. | |
const Traits::GlobalIdSet & | globalIdSet () const |
Access to the GlobalIdSet. | |
const Traits::LocalIdSet & | localIdSet () const |
Access to the LocalIdSet. | |
const Traits::LevelIndexSet & | levelIndexSet (int level) const |
Access to the LevelIndexSets. | |
const Traits::LeafIndexSet & | leafIndexSet () const |
Access to the LeafIndexSet. | |
Traits::LevelGridView | levelGridView (int level) const |
View for a grid level for All_Partition. | |
Traits::LeafGridView | leafGridView () const |
View for the leaf grid for All_Partition. | |
Traits::LeafIntersectionIterator | ileafbegin (const typename Traits::template Codim< 0 >::Entity &entity) const |
obtain begin intersection iterator on the leaf level | |
Traits::LeafIntersectionIterator | ileafend (const typename Traits::template Codim< 0 >::Entity &entity) const |
obtain end intersection iterator on the leaf level | |
Traits::LevelIntersectionIterator | ilevelbegin (const typename Traits::template Codim< 0 >::Entity &entity) const |
obtain begin intersection iterator on the entity level | |
Traits::LevelIntersectionIterator | ilevelend (const typename Traits::template Codim< 0 >::Entity &entity) const |
obtain end intersection iterator on the entity level | |
bool | mark (int refCount, const typename Traits ::template Codim< 0 >::Entity &e) |
Marks an entity to be refined/coarsened in a subsequent adapt. More... | |
bool | loadBalance () |
default implementation of load balance does nothing and returns false | |
Adaptivity and grid refinement | |
bool | mark (int refCount, const typename Codim< 0 >::Entity &e) |
Marks an entity to be refined/coarsened in a subsequent adapt. More... | |
int | getMark (const typename Codim< 0 >::Entity &e) const |
returns adaptation mark for given entity More... | |
Static Public Attributes | |
Exported constants | |
static constexpr int | dimension |
The dimension of the grid. | |
static constexpr int | dimensionworld |
The dimension of the world the grid lives in. | |
Friends | |
class | UGGrid< 2 > |
UGGrid is only implemented for 2 and 3 dimension. | |
Grid Refinement Methods | |
enum | RefinementType { LOCAL , COPY } |
The different forms of grid refinement that UG supports. More... | |
enum | ClosureType { GREEN , NONE } |
Decide whether to add a green closure to locally refined grid sections or not. More... | |
bool | mark (int refCount, const typename Traits::template Codim< 0 >::Entity &e) |
Mark element for refinement. More... | |
bool | mark (const typename Traits::template Codim< 0 >::Entity &e, typename UG_NS< dim >::RefinementRule rule, int side=0) |
Mark method accepting a UG refinement rule. More... | |
int | getMark (const typename Traits::template Codim< 0 >::Entity &e) const |
Query whether element is marked for refinement. | |
bool | preAdapt () |
returns true, if some elements might be coarsend during grid adaption, here always returns true | |
bool | adapt () |
Triggers the grid refinement process. | |
void | postAdapt () |
Clean up refinement markers. | |
template<class DataHandle > | |
bool | loadBalance (DataHandle &dataHandle) |
Distributes the grid and some data over the available nodes in a distributed machine. More... | |
bool | loadBalance (int minlevel=0) |
Distributes this grid over the available nodes in a distributed machine. More... | |
bool | loadBalance (const std::vector< Rank > &targetProcessors, unsigned int fromLevel) |
Distribute this grid over a distributed machine. More... | |
template<class DataHandle > | |
bool | loadBalance (const std::vector< Rank > &targetProcessors, unsigned int fromLevel, DataHandle &dataHandle) |
Distributes the grid over the processes of a parallel machine, and sends data along with it. More... | |
const UGCommunication & | comm () const |
void | getChildrenOfSubface (const typename Traits::template Codim< 0 >::Entity &e, int elementSide, int maxl, std::vector< typename Traits::template Codim< 0 >::Entity > &childElements, std::vector< unsigned char > &childElementSides) const |
Rudimentary substitute for a hierarchic iterator on faces. More... | |
void | setRefinementType (RefinementType type) |
Sets the type of grid refinement. | |
void | setClosureType (ClosureType type) |
Sets the type of grid refinement closure. | |
void | setPosition (const typename Traits::template Codim< dim >::Entity &e, const FieldVector< double, dim > &pos) |
Sets a vertex to a new position. More... | |
void | globalRefine (int n) |
Does uniform refinement. More... | |
void | saveState (const std::string &filename) const |
Save entire grid hierarchy to disk. More... | |
void | loadState (const std::string &filename) |
Read entire grid hierarchy from disk. More... | |
Parallel data distribution and communication | |
Codim< EntitySeed::codimension >::Entity | entity (const EntitySeed &seed) const |
obtain Entity from EntitySeed. | |
GridImp & | asImp () |
Barton-Nackman trick. | |
const GridImp & | asImp () const |
Barton-Nackman trick. | |
Detailed Description
class Dune::UGGrid< dim >
Front-end for the grid manager of the finite element toolbox UG3.
This is the implementation of the grid interface using the UG3 grid management system. It is best described in this paper. To our knowledge, the original code is not available anymore, but the relevant parts have been forked into the Dune module dune-uggrid, available from .
UGGrid provides conforming grids in two and three space dimensions. The grids can be mixed, i.e. 2d grids can contain triangles and quadrilaterals and 3d grids can contain tetrahedra and hexahedra and also pyramids and prisms. The grid refinement rules are very flexible. Local adaptive red/green refinement is the default, but a special method in the UGGrid class allows you to directly access a number of anisotropic refinements rules. Last but not least, the UG grid manager is completely parallelized, and you can use boundaries parametrized by either analytical expressions or high-resolution piecewise linear surfaces.
In your Dune application, you can instantiate objects of the type UGGrid<2> or UGGrid<3>. You can have more than one, if you choose. It is even possible to have 2d and 3d grids at the same time, even though the original UG system never intended to support this!
See the documentation for the factory class GridFactory<UGGrid<dimworld> > to learn how to create UGGrid objects.
Member Enumeration Documentation
◆ ClosureType
enum Dune::UGGrid::ClosureType |
◆ RefinementType
enum Dune::UGGrid::RefinementType |
Member Function Documentation
◆ comm()
|
inline |
the communication
◆ getChildrenOfSubface()
void Dune::UGGrid< dim >::getChildrenOfSubface | ( | const typename Traits::template Codim< 0 >::Entity & | e, |
int | elementSide, | ||
int | maxl, | ||
std::vector< typename Traits::template Codim< 0 >::Entity > & | childElements, | ||
std::vector< unsigned char > & | childElementSides | ||
) | const |
Rudimentary substitute for a hierarchic iterator on faces.
- Parameters
-
e,elementSide Grid face specified by an element and one of its sides maxl The finest level that should be traversed by the iterator [out] childElements For each subface: element index, elementSide, and level [out] childElementSides Indices for transformation because Dune numbers the faces of several elements differently than UG
◆ getMark()
|
inlineinherited |
◆ globalRefine()
void Dune::UGGrid< dim >::globalRefine | ( | int | n | ) |
Does uniform refinement.
- Parameters
-
n Number of uniform refinement steps
◆ loadBalance() [1/4]
bool Dune::UGGrid< dim >::loadBalance | ( | const std::vector< Rank > & | targetProcessors, |
unsigned int | fromLevel | ||
) |
Distribute this grid over a distributed machine.
- Parameters
-
[in] targetProcessors For each leaf element the rank of the process the element shall be sent to [in] fromLevel The lowest level that gets redistributed (set to 0 when in doubt)
This method allows to (re-)distribute the grid controlled by an external grid repartitioning library. You need to get that library to assign a target rank to each interior element in the leaf grid. With this information in a std::vector, call this method, and UG will do the actual repartitioning. Each leaf element will be sent to the assigned target rank. For all other elements we look at where there children are being sent to. The parent is then sent to where most of its children are (Familienzusammenfuehrung).
The size of the input array targetProcessors is expected to be equal to the number of elements in the 'all'-partition, i.e., the number Interior elements plus the number of Ghost elements. To get the array entry corresponding to an Interior element, a MultipleCodimMultipleGeomTypeMapper with layout class MCMGElementLayout is used. Array entries corresponding to Ghost elements are ignored.
In some cases you may also want to leave the lowest levels on one process, to have them all together for multigrid coarse grid corrections. In that case, use the fromLevel parameter with a value other than zero, to redistribute only elements above a certain level.
The fromLevel argument is also needed to allow the compiler to distinguish this method from the loadBalance method with a single template DataHandle argument.
- Note
- In theory you can assign a target rank to any element on any level, and UG will magically transfer the element to that rank and make everything come out right. This is not supported by the UGGrid interface, because I didn't see a use case for it. If you do need it please ask on the Dune mailing list.
- Returns
- true
◆ loadBalance() [2/4]
|
inline |
Distributes the grid over the processes of a parallel machine, and sends data along with it.
- Parameters
-
[in] targetProcessors For each leaf element the rank of the process the element shall be sent to [in] fromLevel The lowest level that gets redistributed (set to 0 when in doubt) [in,out] dataHandle A data handle object that does the gathering and scattering of data
- Template Parameters
-
DataHandle works like the data handle for the communicate methods.
- Returns
- true
References Dune::GridDefaultImplementation< dim, dim, double, UGGridFamily< dim > >::leafGridView(), and Dune::GridDefaultImplementation< dim, dim, double, UGGridFamily< dim > >::loadBalance().
◆ loadBalance() [3/4]
|
inline |
Distributes the grid and some data over the available nodes in a distributed machine.
- Template Parameters
-
DataHandle works like the data handle for the communicate methods.
- Returns
- True, if grid has changed, false otherwise
References Dune::GridDefaultImplementation< dim, dim, double, UGGridFamily< dim > >::leafGridView(), and Dune::GridDefaultImplementation< dim, dim, double, UGGridFamily< dim > >::loadBalance().
◆ loadBalance() [4/4]
bool Dune::UGGrid< dim >::loadBalance | ( | int | minlevel = 0 | ) |
Distributes this grid over the available nodes in a distributed machine.
- Bug:
- The return value is always 'true'
- Parameters
-
minlevel The coarsest grid level that gets distributed
◆ loadState()
void Dune::UGGrid< dim >::loadState | ( | const std::string & | filename | ) |
Read entire grid hierarchy from disk.
Test implementation – not working!
◆ mark() [1/4]
bool Dune::UGGrid< dim >::mark | ( | const typename Traits::template Codim< 0 >::Entity & | e, |
typename UG_NS< dim >::RefinementRule | rule, | ||
int | side = 0 |
||
) |
Mark method accepting a UG refinement rule.
- Parameters
-
e element to be marked for refinement rule One of the UG refinement rules side If rule==UG::D2::BLUE (one quadrilateral is split into two rectangles) you can choose the orientation of the cut by setting side==0 or side==1
The available values for RefinementRule are: (see the RefinementRule enum in ug/gm/gm.h)
2D
- NO_REFINEMENT
- COPY
- RED
- BLUE
- COARSE
- BISECTION_1
- BISECTION_2_Q
- BISECTION_2_T1
- BISECTION_2_T2
- BISECTION_3
3D
- NO_REFINEMENT
- COPY
- RED
- COARSE
- TETRA_RED_HEX
- PRISM_BISECT_1_2
- PRISM_QUADSECT
- PRISM_BISECT_HEX0
- PRISM_BISECT_HEX1
- PRISM_BISECT_HEX2
- PRISM_ROTATE_LEFT
- PRISM_ROTATE_RGHT
- PRISM_QUADSECT_HEXPRI0
- PRISM_RED_HEX
- PRISM_BISECT_0_1
- PRISM_BISECT_0_2
- PRISM_BISECT_0_3
- HEX_BISECT_0_1
- HEX_BISECT_0_2
- HEX_BISECT_0_3
- HEX_TRISECT_0
- HEX_TRISECT_5
- HEX_QUADSECT_0
- HEX_QUADSECT_1
- HEX_QUADSECT_2
- HEX_BISECT_HEXPRI0
- HEX_BISECT_HEXPRI1
◆ mark() [2/4]
|
inlineinherited |
◆ mark() [3/4]
|
inlineinherited |
Marks an entity to be refined/coarsened in a subsequent adapt.
- Parameters
-
[in] refCount Number of subdivisions that should be applied. Negative value means coarsening. [in] e Entity to Entity that should be refined
- Returns
- true if Entity was marked, false otherwise.
- Note
- default implementation is: return false; for grids with no adaptation.
- for the grid programmer: this method is implemented as a template method, because the Entity type is not defined when the class is instantiated You won't need this trick in the implementation. In your implementation you should use it as typename Traits::template Codim<0>::Entity & e ).bool mark(int refCount, const typename Traits::template Codim< 0 >::Entity &e)Mark element for refinement.
◆ mark() [4/4]
bool Dune::UGGrid< dim >::mark | ( | int | refCount, |
const typename Traits::template Codim< 0 >::Entity & | e | ||
) |
Mark element for refinement.
- Parameters
-
refCount - 1: mark for red refinement
- -1: mark for coarsening
- 0: delete a possible refinement mark
e Element to be marked
- Returns
- true, if element was marked
- false, if nothing changed
◆ maxLevel()
int Dune::UGGrid< dim >::maxLevel | ( | ) | const |
Return maximum level defined in this grid. Levels are numbered 0 ... maxlevel with 0 the coarsest level.
Referenced by Dune::UGGrid< dim >::levelIndexSet().
◆ saveState()
void Dune::UGGrid< dim >::saveState | ( | const std::string & | filename | ) | const |
Save entire grid hierarchy to disk.
Test implementation – not working!
◆ setPosition()
void Dune::UGGrid< dim >::setPosition | ( | const typename Traits::template Codim< dim >::Entity & | e, |
const FieldVector< double, dim > & | pos | ||
) |
Sets a vertex to a new position.
Changing a vertex' position changes its position on all grid levels!
The documentation for this class was generated from the following file:
- dune/grid/uggrid.hh