Dune Core Modules (2.3.1)
[ provides Dune::Grid ] More...
#include <dune/grid/yaspgrid.hh>
Classes | |
struct | YGridLevel |
A single grid level within a YaspGrid. More... | |
Public Types | |
typedef yaspgrid_ctype | ctype |
Type used for coordinates. | |
typedef FieldVector< int, dim > | iTupel |
define types used for arguments | |
typedef ReservedVector< YGridLevel, 32 >::const_iterator | YGridLevelIterator |
Iterator over the grid levels. | |
typedef YaspGridFamily< dim > | GridFamily |
the GridFamily of this grid | |
typedef SubYGrid< dim, ctype >::TransformingSubIterator | TSI |
shorthand for some data types | |
Exported constants | |
enum | |
A constant that exports the template parameter dim. | |
enum | |
A constant that exports the template parameter dimworld. | |
Exported types | |
typedef Partition< All_Partition >::LevelGridView | LevelGridView |
View types for All_Partition. | |
typedef Partition< All_Partition >::LeafGridView | LeafGridView |
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. | |
typedef GridFamily::Traits::CollectiveCommunication | CollectiveCommunication |
A type that is a model of Dune::CollectiveCommunication. It provides a portable way for collective communication on the set of processes used by the grid. | |
Public Member Functions | |
const Torus< dim > & | torus () const |
return reference to torus | |
YGridLevelIterator | begin () const |
return iterator pointing to coarsest level | |
YGridLevelIterator | begin (int i) const |
return iterator pointing to given level | |
YGridLevelIterator | end () const |
return iterator pointing to one past the finest level | |
void | MultiYGridSetup (fTupel L, iTupel s, std::bitset< dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) |
The constructor of the old MultiYGrid class. | |
void | MultiYGridSetup (fTupel L, Dune::array< int, dim > s, std::bitset< dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) |
The constructor of the old MultiYGrid class. | |
YaspGrid (Dune::MPIHelper::MPICommunicator comm, Dune::FieldVector< ctype, dim > L, Dune::FieldVector< int, dim > s, Dune::FieldVector< bool, dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) | |
YaspGrid (Dune::FieldVector< ctype, dim > L, Dune::FieldVector< int, dim > s, Dune::FieldVector< bool, dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) | |
YaspGrid (Dune::MPIHelper::MPICommunicator comm, Dune::FieldVector< ctype, dim > L, Dune::array< int, dim > s, std::bitset< dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) | |
YaspGrid (Dune::FieldVector< ctype, dim > L, Dune::array< int, dim > s, std::bitset< dim > periodic, int overlap, const YLoadBalance< dim > *lb=defaultLoadbalancer()) | |
YaspGrid (Dune::FieldVector< ctype, dim > L, Dune::array< int, dim > elements) | |
int | maxLevel () const |
void | globalRefine (int refCount) |
refine the grid refCount times. What about overlap? | |
void | refineOptions (bool keepPhysicalOverlap) |
set options for refinement More... | |
bool | mark (int refCount, const typename Traits::template Codim< 0 >::Entity &e) |
Marks an entity to be refined/coarsened in a subsequent adapt. More... | |
int | getMark (const typename Traits::template Codim< 0 >::Entity &e) const |
returns adaptation mark for given entity More... | |
bool | adapt () |
map adapt to global refine | |
bool | preAdapt () |
returns true, if the grid will be coarsened | |
void | postAdapt () |
clean up some markers | |
template<int cd, PartitionIteratorType pitype> | |
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator | lbegin (int level) const |
one past the end on this level | |
template<int cd, PartitionIteratorType pitype> | |
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator | lend (int level) const |
Iterator to one past the last entity of given codim on level for partition type. | |
template<int cd> | |
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator | lbegin (int level) const |
version without second template parameter for convenience | |
template<int cd> | |
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator | lend (int level) const |
version without second template parameter for convenience | |
template<int cd, PartitionIteratorType pitype> | |
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator | leafbegin () const |
return LeafIterator which points to the first entity in maxLevel | |
template<int cd, PartitionIteratorType pitype> | |
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator | leafend () const |
return LeafIterator which points behind the last entity in maxLevel | |
template<int cd> | |
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator | leafbegin () const |
return LeafIterator which points to the first entity in maxLevel | |
template<int cd> | |
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator | leafend () const |
return LeafIterator which points behind the last entity in maxLevel | |
int | overlapSize (int level, int codim) const |
return size (= distance in graph) of overlap region | |
int | overlapSize (int codim) const |
return size (= distance in graph) of overlap region | |
int | ghostSize (int level, int codim) const |
return size (= distance in graph) of ghost region | |
int | ghostSize (int codim) const |
return size (= distance in graph) of ghost region | |
int | size (int level, int codim) const |
number of entities per level and codim in this process | |
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 |
returns the number of boundary segments within the macro grid | |
template<class DataHandleImp , class DataType > | |
void | communicate (CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir, int level) const |
template<class DataHandleImp , class DataType > | |
void | communicate (CommDataHandleIF< DataHandleImp, DataType > &data, InterfaceType iftype, CommunicationDirection dir) const |
template<class DataHandle , int codim> | |
void | communicateCodim (DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int level) const |
const CollectiveCommunication< MPI_Comm > & | comm () const |
return a collective communication object | |
Traits::template Partition< pitype >::LevelGridView | DUNE_DEPRECATED_MSG ("The method levelView has been renamed to levelGridView.") levelView(int level) const |
View for a grid level. | |
Traits::template Partition< pitype >::LeafGridView | DUNE_DEPRECATED_MSG ("The method levelView has been renamed to leafGridView.") leafView() const |
View for the leaf grid. | |
Traits::template Partition< All_Partition >::LevelGridView | levelView (int level) const DUNE_DEPRECATED_MSG("The method levelView has been renamed to levelGridView.") |
View for a grid level for All_Partition. | |
Traits::template Partition< All_Partition >::LeafGridView | leafView () const DUNE_DEPRECATED_MSG("The method leafView has been renamed to leafGridView.") |
View for the leaf grid for All_Partition. | |
Traits::template Partition< pitype >::LevelGridView | levelGridView (int level) const |
View for a grid level. | |
Traits::template Partition< All_Partition >::LevelGridView | levelGridView (int level) const |
View for a grid level for All_Partition. | |
Traits::template Partition< pitype >::LeafGridView | leafGridView () const |
View for the leaf grid. | |
Traits::template Partition< All_Partition >::LeafGridView | leafGridView () const |
View for the leaf grid for All_Partition. | |
bool | mark (int refCount, const typename Traits ::template Codim< 0 >::Entity &e) |
Marks an entity to be refined/coarsened in a subsequent adapt. More... | |
void | communicate (CommDataHandleIF< DataHandleImp, DataTypeImp > &data, InterfaceType iftype, CommunicationDirection dir, int level) const |
void | communicate (CommDataHandleIF< DataHandleImp, DataTypeImp > &data, InterfaceType iftype, CommunicationDirection dir) const |
bool | loadBalance () |
default implementation of load balance does nothing and returns false | |
bool | loadBalance (DataHandle &data) |
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... | |
Protected Member Functions | |
YGridLevel | makelevel (int level, fTupel L, iTupel s, std::bitset< dim > periodic, iTupel o_interior, iTupel s_interior, int overlap) |
Make a new YGridLevel structure. More... | |
void | intersections (const SubYGrid< dim, ctype > &sendgrid, const SubYGrid< dim, ctype > &recvgrid, const iTupel &size, std::deque< Intersection > &sendlist, std::deque< Intersection > &recvlist) |
Construct list of intersections with neighboring processors. More... | |
Static Protected Member Functions | |
static ReturnImplementationType< InterfaceType >::ImplementationType & | getRealImplementation (InterfaceType &i) |
return real implementation of interface class | |
Parallel data distribution and communication | |
Codim< EntitySeed::codimension >::EntityPointer | entityPointer (const EntitySeed &seed) const |
obtain EntityPointer from EntitySeed. | |
GridImp & | asImp () |
Barton-Nackman trick. | |
const GridImp & | asImp () const |
Barton-Nackman trick. | |
Detailed Description
class Dune::YaspGrid< dim >
[ provides Dune::Grid ]
Provides a distributed structured cube mesh.
YaspGrid stands for yet another structured parallel grid. It implements the dune grid interface for structured grids with codim 0 and dim, with arbitrary overlap (including zero), periodic boundaries and fast implementation allowing on-the-fly computations.
- Template Parameters
-
dim The dimension of the grid and its surrounding world
- History:
- started on July 31, 2004 by PB based on abstractions developed in summer 2003
Constructor & Destructor Documentation
◆ YaspGrid() [1/5]
|
inline |
Constructor
- Parameters
-
comm MPI communicator where this mesh is distributed to L extension of the domain s number of cells on coarse mesh in each direction periodic tells if direction is periodic or not overlap size of overlap on coarsest grid (same in all directions) lb pointer to an overloaded YLoadBalance instance
- Deprecated:
- Will be removed after dune-grid 2.3. Use the corresponding constructor taking array<int> and std::bitset instead.
References Dune::YaspGrid< dim >::MultiYGridSetup().
◆ YaspGrid() [2/5]
|
inline |
Constructor for a sequential YaspGrid
Sequential here means that the whole grid is living on one process even if your program is running in parallel.
- See also
- YaspGrid(Dune::MPIHelper::MPICommunicator, Dune::FieldVector<ctype, dim>, Dune::FieldVector<int, dim>, Dune::FieldVector<bool, dim>, int) for constructing one parallel grid decomposed between the processors.
- Parameters
-
L extension of the domain s number of cells on coarse mesh in each direction periodic tells if direction is periodic or not overlap size of overlap on coarsest grid (same in all directions) lb pointer to an overloaded YLoadBalance instance
- Deprecated:
- Will be removed after dune-grid 2.3. Use the corresponding constructor taking array<int> and std::bitset instead.
References Dune::YaspGrid< dim >::MultiYGridSetup().
◆ YaspGrid() [3/5]
|
inline |
Constructor
- Parameters
-
comm MPI communicator where this mesh is distributed to L extension of the domain s number of cells on coarse mesh in each direction periodic tells if direction is periodic or not overlap size of overlap on coarsest grid (same in all directions) lb pointer to an overloaded YLoadBalance instance
References Dune::YaspGrid< dim >::MultiYGridSetup().
◆ YaspGrid() [4/5]
|
inline |
Constructor for a sequential YaspGrid
Sequential here means that the whole grid is living on one process even if your program is running in parallel.
- See also
- YaspGrid(Dune::MPIHelper::MPICommunicator, Dune::FieldVector<ctype, dim>, Dune::FieldVector<int, dim>, Dune::FieldVector<bool, dim>, int) for constructing one parallel grid decomposed between the processors.
- Parameters
-
L extension of the domain s number of cells on coarse mesh in each direction periodic tells if direction is periodic or not overlap size of overlap on coarsest grid (same in all directions) lb pointer to an overloaded YLoadBalance instance
References Dune::YaspGrid< dim >::MultiYGridSetup().
◆ YaspGrid() [5/5]
|
inline |
Constructor for a sequential YaspGrid without periodicity
Sequential here means that the whole grid is living on one process even if your program is running in parallel.
- See also
- YaspGrid(Dune::MPIHelper::MPICommunicator, Dune::FieldVector<ctype, dim>, Dune::FieldVector<int, dim>, Dune::FieldVector<bool, dim>, int) for constructing one parallel grid decomposed between the processors.
- Parameters
-
L extension of the domain (lower left is always (0,...,0) elements number of cells on coarse mesh in each direction
References Dune::DenseVector< V >::begin(), and Dune::YaspGrid< dim >::makelevel().
Member Function Documentation
◆ communicate() [1/4]
|
inline |
The new communication interface
communicate objects for all codims on the leaf grid
References Dune::YaspGrid< dim >::maxLevel().
◆ communicate() [2/4]
|
inline |
The new communication interface
communicate objects for all codims on a given level
◆ communicate() [3/4]
|
inlineinherited |
dummy communicate, doing nothing
◆ communicate() [4/4]
|
inlineinherited |
dummy communicate, doing nothing
◆ communicateCodim()
|
inline |
The new communication interface
communicate objects for one codim
References Dune::All_All_Interface, Dune::BackwardCommunication, Dune::YaspGrid< dim >::begin(), Dune::InteriorBorder_All_Interface, Dune::InteriorBorder_InteriorBorder_Interface, Dune::Overlap_All_Interface, Dune::Overlap_OverlapFront_Interface, and Dune::YaspGrid< dim >::torus().
◆ getMark() [1/2]
|
inlineinherited |
◆ getMark() [2/2]
|
inline |
returns adaptation mark for given entity
- Parameters
-
[in] e Entity for which adaptation mark should be determined
- Returns
- int adaptation mark, here the default value 0 is returned
References Dune::YaspGrid< dim >::maxLevel().
◆ intersections()
|
inlineprotected |
Construct list of intersections with neighboring processors.
- Parameters
-
recvgrid the grid stored in this processor sendgrid the subgrid to be sent to neighboring processors size needed to shift local grid in periodic case
- Returns
- two lists: Intersections to be sent and Intersections to be received
- Note
- sendgrid/recvgrid may be SubYGrids. Since intersection method is virtual it should work properly
References Dune::YaspGrid< dim >::size().
Referenced by Dune::YaspGrid< dim >::makelevel().
◆ makelevel()
|
inlineprotected |
Make a new YGridLevel structure.
- Parameters
-
L size of the whole domain in each direction s number of cells in each direction periodic indicate periodicity for each direction o_interior origin of interior (non-overlapping) cell decomposition s_interior size of interior cell decomposition overlap to be used on this grid level
References Dune::YaspGrid< dim >::intersections(), and Dune::YaspGrid< dim >::YGridLevel::level_.
Referenced by Dune::YaspGrid< dim >::globalRefine(), Dune::YaspGrid< dim >::MultiYGridSetup(), and Dune::YaspGrid< dim >::YaspGrid().
◆ mark() [1/3]
|
inlineinherited |
◆ mark() [2/3]
|
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)Marks an entity to be refined/coarsened in a subsequent adapt.Definition: yaspgrid.hh:949
◆ mark() [3/3]
|
inline |
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
- On yaspgrid marking one element will mark all other elements of the level as well
- If refCount is lower than refCount of a previous mark-call, nothing is changed
References Dune::YaspGrid< dim >::maxLevel().
◆ maxLevel()
|
inline |
Return maximum level defined in this grid. Levels are numbered 0 ... maxlevel with 0 the coarsest level.
Referenced by Dune::YaspGrid< dim >::begin(), Dune::YaspGrid< dim >::communicate(), Dune::YaspGrid< dim >::end(), Dune::YaspGrid< dim >::getMark(), Dune::YaspGrid< dim >::globalRefine(), Dune::YaspGrid< dim >::leafbegin(), Dune::YaspGrid< dim >::leafend(), Dune::YaspGrid< dim >::mark(), Dune::operator<<(), Dune::YaspGrid< dim >::overlapSize(), and Dune::YaspGrid< dim >::size().
◆ refineOptions()
|
inline |
set options for refinement
- Parameters
-
keepPhysicalOverlap [true] keep the physical size of the overlap, [false] keep the number of cells in the overlap. Default is [true].
The documentation for this class was generated from the following file:
- dune/grid/yaspgrid.hh