DUNE PDELab (2.8)

One-dimensional adaptive grid. More...

#include <dune/grid/onedgrid.hh>

Public Types

enum  RefinementType { LOCAL , COPY }
 The different forms of grid refinement supported by OneDGrid. More...
 
typedef OneDGridGeometry< 0, 1, OneDGrid >::ctype ctype
 The type used to store coordinates.
 
typedef OneDGridFamily GridFamily
 GridFamily of OneDGrid.
 
Exported constants
enum  
 A constant that exports the template parameter dim.
 
enum  
 A constant that exports the template parameter dimworld.
 
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.
 
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

 OneDGrid (const std::vector< ctype > &coords)
 Constructor with an explicit set of coordinates.
 
 OneDGrid (int numElements, const ctype &leftBoundary, const ctype &rightBoundary)
 Constructor for a uniform grid.
 
 ~OneDGrid ()
 Destructor.
 
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.
 
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 coarse grid boundary segments. More...
 
const Traits::GlobalIdSetglobalIdSet () const
 Get the set of global ids.
 
const Traits::LocalIdSetlocalIdSet () const
 Get the set of local ids.
 
const Traits::LevelIndexSetlevelIndexSet (int level) const
 Get an index set for the given level.
 
const Traits::LeafIndexSetleafIndexSet () const
 Get an index set for the leaf level.
 
bool mark (int refCount, const Traits::Codim< 0 >::Entity &e)
 Mark entity for refinement. More...
 
int getMark (const Traits::Codim< 0 >::Entity &e) const
 return current adaptation marker of given entity More...
 
bool preAdapt ()
 Does nothing except return true if some element has been marked for refinement.
 
bool adapt ()
 Triggers the grid refinement process.
 
void postAdapt ()
 Adaptation post-processing: Reset all adaptation state flags.
 
void setRefinementType (RefinementType type)
 Sets the type of grid refinement.
 
void globalRefine (int refCount)
 Does one uniform refinement step. More...
 
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.
 
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, i.e. here the default implementation returns 0. More...
 
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...
 

Static Public Member Functions

template<typename Seed >
static Traits::template Codim< Seed::codimension >::Entity entity (const Seed &seed)
 Create an Entity from an EntitySeed.
 

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

One-dimensional adaptive grid.

[ provides Dune::Grid ]

This implementation of the grid interface provides one-dimensional grids only. The OneDGrid can be nonuniform and provides local mesh refinement and coarsening.

Member Enumeration Documentation

◆ RefinementType

The different forms of grid refinement supported by OneDGrid.

Enumerator
LOCAL 

New level consists only of the refined elements.

COPY 

New level consists of the refined elements and the unrefined ones, too.

Member Function Documentation

◆ getMark() [1/3]

int Dune::OneDGrid::getMark ( const Traits::Codim< 0 >::Entity &  e) const

return current adaptation marker of given entity

Parameters
eEntity to the entity you want to mark
Returns
int current adaptation marker of entity pointer e

◆ getMark() [2/3]

int Dune::Grid< dim, dimworld, OneDGridGeometry< 0, 1, OneDGrid >::ctype , OneDGridFamily >::getMark ( const typename Codim< 0 >::Entity &  e) const
inlineinherited

returns adaptation mark for given entity

Parameters
[in]eEntity for which adaptation mark should be determined
Returns
int adaptation mark currently set for given Entity e

◆ getMark() [3/3]

int Dune::GridDefaultImplementation< dim, dimworld, OneDGridGeometry< 0, 1, OneDGrid >::ctype , OneDGridFamily >::getMark ( const typename Traits::template Codim< 0 >::Entity &  e) const
inlineinherited

returns adaptation mark for given entity, i.e. here the default implementation returns 0.

Parameters
[in]eEntity for which adaptation mark should be determined
Returns
int adaptation mark, here the default value 0 is returned

◆ globalRefine()

void Dune::OneDGrid::globalRefine ( int  refCount)

Does one uniform refinement step.

Parameters
refCountI don't know what this is good for. It doesn't actually do anything.

◆ mark() [1/3]

bool Dune::OneDGrid::mark ( int  refCount,
const Traits::Codim< 0 >::Entity &  e 
)

Mark entity for refinement.

Parameters
refCountif >0 mark for refinement, if <0 mark for coarsening
eEntity to the entity you want to mark
Returns
True, if marking was successful

◆ mark() [2/3]

bool Dune::Grid< dim, dimworld, OneDGridGeometry< 0, 1, OneDGrid >::ctype , OneDGridFamily >::mark ( int  refCount,
const typename Codim< 0 >::Entity &  e 
)
inlineinherited

Marks an entity to be refined/coarsened in a subsequent adapt.

Parameters
[in]refCountNumber of subdivisions that should be applied. Negative value means coarsening.
[in]eEntity that should be marked
Returns
true if Entity was marked, false otherwise.

◆ mark() [3/3]

bool Dune::GridDefaultImplementation< dim, dimworld, OneDGridGeometry< 0, 1, OneDGrid >::ctype , OneDGridFamily >::mark ( int  refCount,
const typename Traits ::template Codim< 0 >::Entity &  e 
)
inlineinherited

Marks an entity to be refined/coarsened in a subsequent adapt.

Parameters
[in]refCountNumber of subdivisions that should be applied. Negative value means coarsening.
[in]eEntity 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
    bool mark( int refCount,
    typename Traits::template Codim<0>::Entity & e ).
    bool mark(int refCount, const Traits::Codim< 0 >::Entity &e)
    Mark entity for refinement.
    Static tag representing a codimension.
    Definition: dimension.hh:22
    This template method will vanish due to the inheritance rules.

◆ maxLevel()

int Dune::OneDGrid::maxLevel ( ) const
inline

Return maximum level defined in this grid.

Levels are numbered 0 ... maxlevel with 0 the coarsest level.

◆ numBoundarySegments()

size_t Dune::OneDGrid::numBoundarySegments ( ) const
inline

Return the number of coarse grid boundary segments.

For this grid implementation, the return value is always 2, because only connected domains are supported, and then the coarse grid boundary consists of two points.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)