Loading [MathJax]/extensions/tex2jax.js

dune-mmesh (1.4)

Dune::MMeshImplicitGridFactory< Grid > Class Template Reference

specialization of the implicit GridFactory for MMesh More...

#include <dune/mmesh/grid/implicitgridfactory.hh>

Public Types

typedef Grid::ctype ctype
 type of (scalar) coordinates
 
typedef Grid::HostGridType HostGrid
 type of the hostgrid
 
typedef FieldVector< ctype, dimensionworldWorldVector
 type of vector for world coordinates
 
typedef FieldMatrix< ctype, dimensionworld, dimensionworldWorldMatrix
 type of matrix from world coordinates to world coordinates
 
typedef Dune::BoundarySegment< dimension, dimensionworldBoundarySegment
 type of a Dune boundary segment
 
typedef Grid::IdType IdType
 type of an id
 
typedef std::unordered_map< IdType, std::size_t > BoundarySegments
 type of the boundary segment id map
 
typedef std::unordered_map< IdType, std::size_t > InterfaceSegments
 type of the interface segment set
 

Public Member Functions

 MMeshImplicitGridFactory ()
 
void insertElement (const GeometryType &type, const std::vector< unsigned int > &v)
 insert an element into the macro grid More...
 
template<int d = dimension>
std::enable_if_t< d==2, bool > isElement (const std::vector< unsigned int > &v, Element_handle &fh) const
 Returns if there is a face with the given vertices in the triangulation 2. More...
 
template<int d = dimension>
std::enable_if_t< d==3, bool > isElement (const std::vector< unsigned int > &v, Element_handle &fh) const
 Returns if there is a cell with the given vertices in the triangulation 3. More...
 
void insertBoundarySegment (const std::vector< unsigned int > &vertices)
 insert boundary segment More...
 
void insertVertex (const WorldVector &pos)
 Insert a vertex into the macro grid. More...
 
void insertInterface (const std::vector< unsigned int > &vertices, const std::size_t marker=1)
 insert an interface into the macro grid More...
 
unsigned int insertionIndex (const typename Codim< 0 >::Entity &entity) const
 return index of inserted vertex within the macro grid More...
 
unsigned int insertionIndex (const typename Codim< dimension >::Entity &entity) const
 return insertion index of entity More...
 
unsigned int insertionIndex (const typename Grid::LeafIntersection &intersection) const
 return insertion index of boundary intersection More...
 
const BoundarySegmentsboundarySegments () const
 returns the boundary segment to index map
 
const BoundaryIds & boundaryIds () const
 returns the boundary segment index to boundary id map
 
void addBoundaryId (std::size_t boundarySegmentIndex, std::size_t boundaryId)
 add a boundary id
 
std::unique_ptr< Grid > createGrid ()
 finalize grid creation and hand over the grid More...
 

Static Public Attributes

static const int dimension = Grid::dimension
 dimension of the grid
 
static const int dimensionworld = Grid::dimensionworld
 dimension of the world
 
static const bool supportsBoundaryIds = true
 are boundary ids supported by this factory?
 
static const bool supportPeriodicity = false
 the factory is not able to create periodic meshes
 

Detailed Description

template<class Grid>
class Dune::MMeshImplicitGridFactory< Grid >

specialization of the implicit GridFactory for MMesh

The implicit grid factory for MMesh

Constructor & Destructor Documentation

◆ MMeshImplicitGridFactory()

template<class Grid >
Dune::MMeshImplicitGridFactory< Grid >::MMeshImplicitGridFactory ( )
inline

default constructor

Member Function Documentation

◆ createGrid()

template<class Grid >
std::unique_ptr< Grid > Dune::MMeshImplicitGridFactory< Grid >::createGrid ( )
inline

finalize grid creation and hand over the grid

This version of createGrid is original to the MMesh grid factroy, allowing to specity a grid name.

Returns
a pointer to the newly created grid

◆ insertBoundarySegment()

template<class Grid >
void Dune::MMeshImplicitGridFactory< Grid >::insertBoundarySegment ( const std::vector< unsigned int > &  vertices)
inline

insert boundary segment

Parameters
[in]verticesVertices

References Dune::MMeshImplicitGridFactory< Grid >::dimension.

◆ insertElement()

template<class Grid >
void Dune::MMeshImplicitGridFactory< Grid >::insertElement ( const GeometryType &  type,
const std::vector< unsigned int > &  v 
)
inline

insert an element into the macro grid

Parameters
[in]typeGeometryType of the new element
[in]vindices of the element vertices (starting with 0)
Note
The implicit grid factory ignores the insertion of specific elements. The CGAL backend constructs the elements by inserting vertices implicitly.

References Dune::MMeshImplicitGridFactory< Grid >::isElement().

◆ insertInterface()

template<class Grid >
void Dune::MMeshImplicitGridFactory< Grid >::insertInterface ( const std::vector< unsigned int > &  vertices,
const std::size_t  marker = 1 
)
inline

insert an interface into the macro grid

Parameters
[in]verticesindices of the interface vertices (starting with 0)
[in]markermarker value of the interface segment (default 1)

References Dune::MMeshImplicitGridFactory< Grid >::dimension.

◆ insertionIndex() [1/3]

template<class Grid >
unsigned int Dune::MMeshImplicitGridFactory< Grid >::insertionIndex ( const typename Codim< 0 >::Entity &  entity) const
inline

return index of inserted vertex within the macro grid

Parameters
[in]posposition of the vertex (in world coordinates)

◆ insertionIndex() [2/3]

template<class Grid >
unsigned int Dune::MMeshImplicitGridFactory< Grid >::insertionIndex ( const typename Codim< dimension >::Entity &  entity) const
inline

return insertion index of entity

Parameters
[in]entityEntity

◆ insertionIndex() [3/3]

template<class Grid >
unsigned int Dune::MMeshImplicitGridFactory< Grid >::insertionIndex ( const typename Grid::LeafIntersection &  intersection) const
inline

return insertion index of boundary intersection

Parameters
[in]intersectionLeaf intersection

◆ insertVertex()

template<class Grid >
void Dune::MMeshImplicitGridFactory< Grid >::insertVertex ( const WorldVector pos)
inline

Insert a vertex into the macro grid.

Parameters
[in]posposition of the vertex (in world coordinates)
Note
This method assumes that the vertices are inserted consecutively with respect to their index.

◆ isElement() [1/2]

template<class Grid >
template<int d = dimension>
std::enable_if_t< d==2, bool > Dune::MMeshImplicitGridFactory< Grid >::isElement ( const std::vector< unsigned int > &  v,
Element_handle &  fh 
) const
inline

Returns if there is a face with the given vertices in the triangulation 2.

Parameters
[in]vindices of the element vertices

References Dune::MMeshImplicitGridFactory< Grid >::dimension.

Referenced by Dune::MMeshImplicitGridFactory< Grid >::insertElement().

◆ isElement() [2/2]

template<class Grid >
template<int d = dimension>
std::enable_if_t< d==3, bool > Dune::MMeshImplicitGridFactory< Grid >::isElement ( const std::vector< unsigned int > &  v,
Element_handle &  fh 
) const
inline

Returns if there is a cell with the given vertices in the triangulation 3.

Parameters
[in]vindices of the element vertices

References Dune::MMeshImplicitGridFactory< Grid >::dimension.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 7, 22:57, 2025)