Dune Namespace Reference

The namespace Dune is the main namespace for all Dune code. More...

Namespaces

namespace  Capabilities
 
namespace  cpgrid
 
namespace  GeometryHelpers
 

Classes

class  CartesianIndexMapper
 Interface class to access the logical Cartesian grid as used in industry standard simulator decks. More...
 
class  CartesianIndexMapper< CpGrid >
 
class  CartesianIndexMapper< PolyhedralGrid< dim, dimworld, coord_t > >
 
class  CpGrid
 [ provides Dune::Grid ] More...
 
struct  CpGridFamily
 
struct  CpGridTraits
 
struct  DGFGridFactory< CpGrid >
 
struct  DGFGridFactory< PolyhedralGrid< dim, dimworld, coord_t > >
 
struct  DGFGridInfo< CpGrid >
 
struct  DGFGridInfo< PolyhedralGrid< dim, dimworld > >
 
class  GridFactory< PolyhedralGrid< dim, dimworld, coord_t > >
 
class  NonBlockingExchangeImplementation
 
struct  OrderByFirst
 
class  PersistentContainer< CpGrid, Data >
 
class  PersistentContainer< PolyhedralGrid< dim, dimworld >, Data >
 
class  Point2PointCommunicator
 Point-2-Point communicator for exchange messages between processes. More...
 
class  PolyhedralGrid
 identical grid wrapper More...
 
struct  PolyhedralGridBasicGeometry
 
class  PolyhedralGridEntity
 
class  PolyhedralGridEntity< 0, dim, Grid >
 construct a null entity More...
 
class  PolyhedralGridEntityBasic
 
class  PolyhedralGridEntityPointer
 
class  PolyhedralGridEntitySeed
 
struct  PolyhedralGridFamily
 
class  PolyhedralGridGeometry
 
class  PolyhedralGridIdSet
 
class  PolyhedralGridIndexSet
 
class  PolyhedralGridIntersection
 
class  PolyhedralGridIntersectionIterator
 
class  PolyhedralGridIterator
 
class  PolyhedralGridLocalGeometry
 
class  PolyhedralGridView
 
struct  PolyhedralGridViewTraits
 
class  SimpleMessageBuffer
 
class  SubGridPart
 A class to represent a part of a grid, similar to a GridView. More...
 
struct  SubGridPartTraits
 

Enumerations

enum  EdgeWeightMethod { uniformEdgeWgt =0 , defaultTransEdgeWgt =1 , logTransEdgeWgt =2 }
 enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Metis' graph partitioner. More...
 
enum  PartitionMethod { simple =0 , zoltan =1 , metis =2 , zoltanGoG =3 }
 enum for choosing methods for partitioning a graph. More...
 

Functions

void partition (const CpGrid &grid, const std::array< int, 3 > &initial_split, int &num_part, std::vector< int > &cell_part, bool recursive=false, bool ensureConnectivity=true)
 
void addOverlapLayer (const CpGrid &grid, const std::vector< int > &cell_part, std::vector< std::set< int > > &cell_overlap, int mypart, int overlapLayers, bool all=false, int level=-1)
 Adds a layer of overlap cells to a partitioning. More...
 
int addOverlapLayer (const CpGrid &grid, const std::vector< int > &cell_part, std::vector< std::tuple< int, int, char > > &exportList, std::vector< std::tuple< int, int, char, int > > &importList, const Communication< Dune::MPIHelper::MPICommunicator > &cc, bool addCornerCells, const double *trans, int layers=1, int level=-1)
 Adds a layer of overlap cells to a partitioning. More...
 
template<typename T >
FieldVector< T, 3 > cross (const FieldVector< T, 3 > &a, const FieldVector< T, 3 > &b)
 
template<class Vector >
Vector::field_type inner (const Vector &a, const Vector &b)
 
template<typename T , template< typename, int > class Point>
determinantOf (const Point< T, 2 > *a)
 
template<typename T , template< typename, int > class Point>
determinantOf (const Point< T, 3 > *a)
 
template<typename T , template< typename, int > class Point, int Dim>
simplex_volume (const Point< T, Dim > *a)
 
template<typename T , template< typename, int > class Point>
area (const Point< T, 2 > *c)
 
template<typename T , template< typename, int > class Point>
area (const Point< T, 3 > *c)
 
template<typename T , template< typename, int > class Point>
volume (const Point< T, 3 > *c)
 Computes the volume of a 3D simplex (embedded i 3D space). More...
 
template<typename T , template< typename, int > class Point>
signed_area (const Point< T, 3 > *c, const Point< T, 3 > &normal)
 
template<int dim>
cpgrid::Entity< dim > createEntity (const CpGrid &, int, bool)
 
template<typename K , int n>
FieldVector< K, n > operator- (const FieldVector< K, n > &v)
 

Detailed Description

The namespace Dune is the main namespace for all Dune code.

Enumeration Type Documentation

◆ EdgeWeightMethod

enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Metis' graph partitioner.

uniform methods means all edges have weight 1. defaultTrans uses transmissibility as weights. logTrans uses the logarithm of the transmissibility. The uniform and logTrans edge-weighting methods produce partitioning results with lower edge-cut, fewer overlap/ghost cells and less communication overhead than when using defaultTrans. However, the impact on parallel linear solver performance is negative.

Enumerator
uniformEdgeWgt 

All edge have a uniform weight of 1.

defaultTransEdgeWgt 

Use the transmissibilities as edge weights.

logTransEdgeWgt 

Use the log of the transmissibilities as edge weights.

◆ PartitionMethod

enum for choosing methods for partitioning a graph.

Enumerator
simple 

Use simple approach based on rectangular partitioning the underlying cartesian grid.

zoltan 

Use Zoltan for partitioning.

metis 

Use METIS for partitioning.

zoltanGoG 

use Zoltan on GraphOfGrid for partitioning

Function Documentation

◆ addOverlapLayer() [1/2]

void Dune::addOverlapLayer ( const CpGrid grid,
const std::vector< int > &  cell_part,
std::vector< std::set< int > > &  cell_overlap,
int  mypart,
int  overlapLayers,
bool  all = false,
int  level = -1 
)

Adds a layer of overlap cells to a partitioning.

Parameters
[in]gridThe grid that is partitioned.
[in]cell_parta vector containing each cells partition number.
[out]cell_overlapa vector of sets that contains for each cell all the partition numbers that it is an overlap cell of.
[in]mypartThe partition number of the processor.
[in]allWhether to compute the overlap for all partions or just the one associated by mypart.
[in]levelIndicating level grid that is partitioned. Default -1 for leaf grid view.

◆ addOverlapLayer() [2/2]

int Dune::addOverlapLayer ( const CpGrid grid,
const std::vector< int > &  cell_part,
std::vector< std::tuple< int, int, char > > &  exportList,
std::vector< std::tuple< int, int, char, int > > &  importList,
const Communication< Dune::MPIHelper::MPICommunicator > &  cc,
bool  addCornerCells,
const double *  trans,
int  layers = 1,
int  level = -1 
)

Adds a layer of overlap cells to a partitioning.

Parameters
[in]gridThe grid that is partitioned.
[in]cell_parta vector containing each cells partition number.
[in,out]exportListList indices to export, each entry is a tuple of global index, process rank (to export to), attribute on remote.
[in,out]importListList indices to import, each entry is a tuple of global index, process rank (to import from), attribute here, local index here
[in]ccThe communication object
[in]addCornerCellsSwitch for adding corner cells to overlap layer.
[in]transThe transmissibilities on cell faces. When trans[i]==0, no overlap is added.
[in]layerNumber of overlap layers
[in]levelIndicating level grid that is partitioned. Default -1 for leaf grid view.

◆ area() [1/2]

template<typename T , template< typename, int > class Point>
T Dune::area ( const Point< T, 2 > *  c)
inline

Computes the area of a 2-dimensional triangle. Input is an array of corner points. Same function also exists for 3-dimensional triangles.

References simplex_volume().

Referenced by Dune::GeometryHelpers::polygonArea(), Dune::GeometryHelpers::polygonCentroid(), Dune::GeometryHelpers::polygonNormal(), and Dune::cpgrid::Geometry< 3, cdim >::refineCellifiedPatch().

◆ area() [2/2]

template<typename T , template< typename, int > class Point>
T Dune::area ( const Point< T, 3 > *  c)
inline

Computes the area of a 3-dimensional triangle. Input is an array of corner points. Same function also exists for 2-dimensional triangles.

References cross().

◆ createEntity()

template<int dim>
cpgrid::Entity< dim > Dune::createEntity ( const CpGrid ,
int  ,
bool   
)

◆ cross()

template<typename T >
FieldVector< T, 3 > Dune::cross ( const FieldVector< T, 3 > &  a,
const FieldVector< T, 3 > &  b 
)
Template Parameters

param

Returns

References b().

Referenced by area(), Dune::GeometryHelpers::polygonNormal(), and signed_area().

◆ determinantOf() [1/2]

template<typename T , template< typename, int > class Point>
T Dune::determinantOf ( const Point< T, 2 > *  a)
inline

Calculates the determinant of a 2 x 2 matrix, represented in memory as an array of two-dimensional points. Same function also exists for 3 x 3 matrices.

Referenced by simplex_volume().

◆ determinantOf() [2/2]

template<typename T , template< typename, int > class Point>
T Dune::determinantOf ( const Point< T, 3 > *  a)
inline

Calculates the determinant of a 3 x 3 matrix, represented in memory as an array of three-dimensional points. Same function also exists for 2 x 2 matrices.

◆ inner()

template<class Vector >
Vector::field_type Dune::inner ( const Vector &  a,
const Vector &  b 
)
Template Parameters

param

Returns

References b().

Referenced by signed_area().

◆ operator-()

template<typename K , int n>
FieldVector< K, n > Dune::operator- ( const FieldVector< K, n > &  v)

◆ partition()

void Dune::partition ( const CpGrid grid,
const std::array< int, 3 > &  initial_split,
int &  num_part,
std::vector< int > &  cell_part,
bool  recursive = false,
bool  ensureConnectivity = true 
)

Partition a CpGrid based on (ijk) coordinates, with splitting to ensure that each partition is connected.

Parameters
[in]gridthe grid to partition
[in]initial_splitthe number of parts in which to partition the grid, in each cardinal direction. Their product is the expected number of partitions produced.
[out]num_partthe resulting number of partitions. This may be lower than expected, because of inactive cells, or higher than expected, because of splits to ensure connectedness.
[out]cell_parta vector containing, for each cell, its partition number

◆ signed_area()

template<typename T , template< typename, int > class Point>
T Dune::signed_area ( const Point< T, 3 > *  c,
const Point< T, 3 > &  normal 
)

Computes the signed area of a triangle embedded in 3D space. Input is an array of corner points and a normal to determine the sign.

References cross(), and inner().

◆ simplex_volume()

template<typename T , template< typename, int > class Point, int Dim>
T Dune::simplex_volume ( const Point< T, Dim > *  a)
inline

Computes the volume of a simplex consisting of (Dim+1) vertices embedded in Euclidean space of dimension (Dim)

References determinantOf().

Referenced by area(), Dune::GeometryHelpers::polygonCellCentroid(), Dune::GeometryHelpers::polygonCellVolume(), Dune::cpgrid::Geometry< 3, cdim >::refineCellifiedPatch(), and volume().

◆ volume()

template<typename T , template< typename, int > class Point>
T Dune::volume ( const Point< T, 3 > *  c)
inline