[ provides Dune::Grid ]
More...
#include <CpGrid.hpp>
|
| CpGrid () |
| Default constructor. More...
|
|
| CpGrid (MPIHelper::MPICommunicator comm) |
|
void | processEclipseFormat (const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false) |
|
|
A cornerpoint grid can be seen as a degenerated and distorted cartesian grid. Therefore it provides mappings from cells to the underlying cartesian index.
|
void | createCartesian (const std::array< int, 3 > &dims, const std::array< double, 3 > &cellsize, const std::array< int, 3 > &shift={0, 0, 0}) |
|
const std::array< int, 3 > & | logicalCartesianSize () const |
|
const std::vector< int > & | globalCell () const |
|
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & | currentData () const |
| Returns either data_ or distributed_data_(if non empty). More...
|
|
std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & | currentData () |
| Returns either data_ or distributed_data_(if non empty). More...
|
|
void | getIJK (const int c, std::array< int, 3 > &ijk) const |
| Extract Cartesian index triplet (i,j,k) of an active cell. More...
|
|
bool | uniqueBoundaryIds () const |
|
void | setUniqueBoundaryIds (bool uids) |
|
|
It provides additional methods not in the DUNE interface but needed by OPM. Vertices, faces, and cells are not represented as entities but identified by indices.
|
typedef Dune::FieldVector< double, 3 > | Vector |
|
const std::vector< double > & | zcornData () const |
|
int | numCells (int level=-1) const |
| Get the number of cells. More...
|
|
int | numFaces (int level=-1) const |
| Get the number of faces. More...
|
|
int | numVertices () const |
| Get The number of vertices. More...
|
|
int | numCellFaces (int cell, int level=-1) const |
| Get the number of faces of a cell. More...
|
|
int | cellFace (int cell, int local_index, int level=-1) const |
| Get a specific face of a cell. More...
|
|
const cpgrid::OrientedEntityTable< 0, 1 >::row_type | cellFaceRow (int cell) const |
| Get a list of indices identifying all faces of a cell. More...
|
|
int | faceCell (int face, int local_index, int level=-1) const |
| Get the index identifying a cell attached to a face. More...
|
|
int | numCellFaces () const |
| Get the sum of all faces attached to all cells. More...
|
|
int | numFaceVertices (int face) const |
|
int | faceVertex (int face, int local_index) const |
| Get the index identifying a vertex of a face. More...
|
|
double | cellCenterDepth (int cell_index) const |
| Get vertical position of cell center ("zcorn" average). More...
|
|
const Vector | faceCenterEcl (int cell_index, int face, const Dune::cpgrid::Intersection &intersection) const |
|
const Vector | faceAreaNormalEcl (int face) const |
|
const Vector & | vertexPosition (int vertex) const |
| Get the Position of a vertex. More...
|
|
double | faceArea (int face) const |
| Get the area of a face. More...
|
|
const Vector & | faceCentroid (int face) const |
| Get the coordinates of the center of a face. More...
|
|
const Vector & | faceNormal (int face) const |
| Get the unit normal of a face. More...
|
|
double | cellVolume (int cell) const |
| Get the volume of the cell. More...
|
|
const Vector & | cellCentroid (int cell) const |
| Get the coordinates of the center of a cell. More...
|
|
CentroidIterator< 0 > | beginCellCentroids () const |
| Get an iterator over the cell centroids positioned at the first one. More...
|
|
CentroidIterator< 1 > | beginFaceCentroids () const |
| Get an iterator over the face centroids positioned at the first one. More...
|
|
int | boundaryId (int face) const |
|
template<class Cell2FacesRowIterator > |
int | faceTag (const Cell2FacesRowIterator &cell_face) const |
| Get the cartesian tag associated with a face tag. More...
|
|
|
std::string | name () const |
| Get the grid name. More...
|
|
int | maxLevel () const |
| Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview), 0 = the coarsest level. More...
|
|
template<int codim> |
Traits::template Codim< codim >::LevelIterator | lbegin (int level) const |
| Iterator to first entity of given codim on level. More...
|
|
template<int codim> |
Traits::template Codim< codim >::LevelIterator | lend (int level) const |
| one past the end on this level More...
|
|
template<int codim, PartitionIteratorType PiType> |
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator | lbegin (int level) const |
| Iterator to first entity of given codim on level and PartitionIteratorType. More...
|
|
template<int codim, PartitionIteratorType PiType> |
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator | lend (int level) const |
| one past the end on this level More...
|
|
template<int codim> |
Traits::template Codim< codim >::LeafIterator | leafbegin () const |
| Iterator to first leaf entity of given codim. More...
|
|
template<int codim> |
Traits::template Codim< codim >::LeafIterator | leafend () const |
| one past the end of the sequence of leaf entities More...
|
|
template<int codim, PartitionIteratorType PiType> |
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator | leafbegin () const |
| Iterator to first leaf entity of given codim and PartitionIteratorType. More...
|
|
template<int codim, PartitionIteratorType PiType> |
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator | leafend () const |
| one past the end of the sequence of leaf entities More...
|
|
int | size (int level, int codim) const |
| Number of grid entities per level and codim. More...
|
|
int | size (int codim) const |
| number of leaf entities per codim in this process More...
|
|
int | size (int level, GeometryType type) const |
| number of entities per level and geometry type in this process More...
|
|
int | size (GeometryType type) const |
| number of leaf entities per geometry type in this process More...
|
|
const Traits::GlobalIdSet & | globalIdSet () const |
| Access to the GlobalIdSet. More...
|
|
const Traits::LocalIdSet & | localIdSet () const |
| Access to the LocalIdSet. More...
|
|
const Traits::LevelIndexSet & | levelIndexSet (int level) const |
| Access to the LevelIndexSets. More...
|
|
const Traits::LeafIndexSet & | leafIndexSet () const |
| Access to the LeafIndexSet. More...
|
|
void | globalRefine (int refCount) |
| Refine the grid refCount times using the default refinement rule. This behaves like marking all elements for refinement and then calling preAdapt, adapt and postAdapt. The state after globalRefine is comparable to the state after postAdapt. More...
|
|
const std::vector< Dune::GeometryType > & | geomTypes (const int) const |
|
template<int codim> |
cpgrid::Entity< codim > | entity (const cpgrid::Entity< codim > &seed) const |
| given an EntitySeed (or EntityPointer) return an entity object More...
|
|
void | addLgrsUpdateLeafView (const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, const std::vector< std::string > &lgr_name_vec) |
| Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint patches of cells from that coarse grid. More...
|
|
const std::map< std::string, int > & | getLgrNameToLevel () const |
|
std::array< double, 3 > | getEclCentroid (const int &idx) const |
|
std::array< double, 3 > | getEclCentroid (const cpgrid::Entity< 0 > &elem) const |
|
Dune::cpgrid::Intersection | getParentIntersectionFromLgrBoundaryFace (const Dune::cpgrid::Intersection &intersection) const |
|
bool | mark (int refCount, const cpgrid::Entity< 0 > &element) |
| Mark entity for refinement (or coarsening). More...
|
|
int | getMark (const cpgrid::Entity< 0 > &element) const |
| Return refinement mark for entity. More...
|
|
bool | preAdapt () |
| Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called after elements have been marked for refinement. More...
|
|
bool | adapt () |
| Triggers the grid refinement process. Returns true if the grid has changed, false otherwise. More...
|
|
bool | adapt (const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &assignRefinedLevel, const std::vector< std::string > &lgr_name_vec, const std::vector< std::array< int, 3 > > &startIJK_vec=std::vector< std::array< int, 3 > >{}, const std::vector< std::array< int, 3 > > &endIJK_vec=std::vector< std::array< int, 3 > >{}) |
| Triggers the grid refinement process, allowing to select diffrent refined level grids. More...
|
|
void | postAdapt () |
| Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'. More...
|
|
void | syncDistributedGlobalCellIds () |
| Synchronizes cell global ids across processes after load balancing. More...
|
|
std::vector< std::unordered_map< std::size_t, std::size_t > > | mapLocalCartesianIndexSetsToLeafIndexSet () const |
| Compute for each level grid, a map from the global_cell_[ cell index in level grid ] to the leaf index of the equivalent cell on the leaf grid view. Notice that cells that vanished and do not appear on the leaf grid view will not be considered. global_cell_[ cell index in level grid ] coincide with (local) Cartesian Index. More...
|
|
std::vector< std::array< int, 2 > > | mapLeafIndexSetToLocalCartesianIndexSets () const |
| Reverse map: from leaf index cell to { level, local/level Cartesian index of the cell }. More...
|
|
unsigned int | overlapSize (int) const |
| Size of the overlap on the leaf level. More...
|
|
unsigned int | ghostSize (int) const |
| Size of the ghost cell layer on the leaf level. More...
|
|
unsigned int | overlapSize (int, int) const |
| Size of the overlap on a given level. More...
|
|
unsigned int | ghostSize (int, int) const |
| Size of the ghost cell layer on a given level. More...
|
|
unsigned int | numBoundarySegments () const |
| returns the number of boundary segments within the macro grid More...
|
|
void | setPartitioningParams (const std::map< std::string, std::string > ¶ms) |
|
bool | loadBalance (int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan, double imbalanceTol=1.1, int level=-1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
bool | loadBalanceSerial (int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan, int edgeWeightMethod=Dune::EdgeWeightMethod::defaultTransEdgeWgt, double imbalanceTol=1.1, int level=-1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
std::pair< bool, std::vector< std::pair< std::string, bool > > > | loadBalance (const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, int level=-1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
std::pair< bool, std::vector< std::pair< std::string, bool > > > | loadBalance (EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, int level=-1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
template<class DataHandle > |
std::pair< bool, std::vector< std::pair< std::string, bool > > > | loadBalance (DataHandle &data, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, const double *transmissibilities=nullptr, int overlapLayers=1, int partitionMethod=1, int level=-1) |
| Distributes this grid and data over the available nodes in a distributed machine. More...
|
|
template<class DataHandle > |
std::pair< bool, std::vector< std::pair< std::string, bool > > > | loadBalance (DataHandle &data, EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, bool serialPartitioning, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltanGoG, double imbalanceTol=1.1, bool allowDistributedWells=false) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
template<class DataHandle > |
std::pair< bool, std::vector< std::pair< std::string, bool > > > | loadBalance (DataHandle &data, const std::vector< int > &parts, const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections={}, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
template<class DataHandle > |
bool | loadBalance (DataHandle &data, decltype(data.fixedSize(0, 0)) overlapLayers=1, int partitionMethod=Dune::PartitionMethod::zoltan) |
| Distributes this grid and data over the available nodes in a distributed machine. More...
|
|
bool | loadBalance (const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1) |
| Distributes this grid over the available nodes in a distributed machine. More...
|
|
template<class DataHandle > |
bool | loadBalance (DataHandle &data, const std::vector< int > &parts, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1) |
| Distributes this grid and data over the available nodes in a distributed machine. More...
|
|
std::vector< int > | zoltanPartitionWithoutScatter (const std::vector< cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const double *transmissibilities, const int numParts, const double imbalanceTol) const |
| Partitions the grid using Zoltan without decomposing and distributing it among processes. More...
|
|
template<class DataHandle > |
void | communicate (DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const |
| communicate objects for all codims on a given level More...
|
|
template<class DataHandle > |
void | communicate (DataHandle &data, InterfaceType iftype, CommunicationDirection dir) const |
| communicate objects for all codims on a given level. More...
|
|
const CpGridTraits::Communication & | comm () const |
| Get the collective communication object. More...
|
|
int | replaceLgr1CornerIdxByLgr2CornerIdx (const std::array< int, 3 > &cells_per_dim_lgr1, int cornerIdxLgr1, const std::array< int, 3 > &cells_per_dim_lgr2) const |
| A refined corner appears in two single-cell-refinements. Given the corner index in the first single-cell-refinement, compute the corner index in the neighboring single-cell-refinement. More...
|
|
int | replaceLgr1CornerIdxByLgr2CornerIdx (const std::array< int, 3 > &cells_per_dim_lgr1, int cornerIdxLgr1, int elemLgr1, int parentFaceLastAppearanceIdx, const std::array< int, 3 > &cells_per_dim_lgr2) const |
| A new refined corner lays on an edge and appears in at least two single-cell-refinements. Given the corner index in one single-cell-refinement, compute the corner index in a neighboring single-cell-refinement. More...
|
|
int | replaceLgr1FaceIdxByLgr2FaceIdx (const std::array< int, 3 > &cells_per_dim_lgr1, int faceIdxInLgr1, const std::shared_ptr< cpgrid::CpGridData > &elemLgr1_ptr, const std::array< int, 3 > &cells_per_dim_lgr2) const |
| A new refined face lays on the boudndary of a single-cell-refinement appears in at most two single-cell-refinements. Given the face index in one single-cell-refinement, compute the face index in a neighboring single-cell-refinement. More...
|
|
◆ CommunicationType
The type of the owner-overlap-copy communication.
◆ GridFamily
Family typedef, why is this not defined by Grid<>?
◆ InterfaceMap
The type of the map describing communication interfaces.
◆ ParallelIndexSet
The type of the parallel index set.
◆ RemoteIndices
The type of the remote indices information.
◆ Vector
◆ CpGrid() [1/2]
◆ CpGrid() [2/2]
Dune::CpGrid::CpGrid |
( |
MPIHelper::MPICommunicator |
comm | ) |
|
|
explicit |
◆ adapt() [1/2]
bool Dune::CpGrid::adapt |
( |
| ) |
|
Triggers the grid refinement process. Returns true if the grid has changed, false otherwise.
◆ adapt() [2/2]
bool Dune::CpGrid::adapt |
( |
const std::vector< std::array< int, 3 > > & |
cells_per_dim_vec, |
|
|
const std::vector< int > & |
assignRefinedLevel, |
|
|
const std::vector< std::string > & |
lgr_name_vec, |
|
|
const std::vector< std::array< int, 3 > > & |
startIJK_vec = std::vector< std::array< int, 3 > >{} , |
|
|
const std::vector< std::array< int, 3 > > & |
endIJK_vec = std::vector< std::array< int, 3 > >{} |
|
) |
| |
Triggers the grid refinement process, allowing to select diffrent refined level grids.
- Parameters
-
[in] | cells_per_dim_vec | For each set of marked elements for refinement, that will belong to a same refined level grid, number of (refined) cells in each direction that each parent cell should be refined to. |
[in] | assignRefinedLevel | Vector with size equal to total amount of cells of the starting grid where the marked elements belong. In each entry, the refined level grid where the refined entities of the (parent) marked element should belong is stored. |
[in] | lgr_name_vector | Each refined level grid name, e.g. {"LGR1", "LGR2"}. |
[in] | startIJK_vec | Default empty vector. When isCARFIN, the starting ijk Cartesian index of each block of cells to be refined. |
[in] | endIJK_vec | Default empty vector. When isCARFIN, the final ijk Cartesian index of each block of cells to be refined. |
◆ addLgrsUpdateLeafView()
void Dune::CpGrid::addLgrsUpdateLeafView |
( |
const std::vector< std::array< int, 3 > > & |
cells_per_dim_vec, |
|
|
const std::vector< std::array< int, 3 > > & |
startIJK_vec, |
|
|
const std::vector< std::array< int, 3 > > & |
endIJK_vec, |
|
|
const std::vector< std::string > & |
lgr_name_vec |
|
) |
| |
Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint patches of cells from that coarse grid.
Level0 refers to the coarse grid, assumed to be this-> data_[0]. Level1 and level2 refer to the LGRs (stored in this->data_[1] data_[2]). LeafView (stored in this-> data_[3]) is built with the level0-entities which weren't involded in the refinenment, together with the new born entities created in level1 and level2. Old-corners and old-faces (from coarse grid) lying on the boundary of the patches, get replaced by new-born-equivalent corners and new-born-faces.
- Parameters
-
[in] | cells_per_dim_vec | Vector of Number of (refined) cells in each direction that each parent cell should be refined to. |
[in] | startIJK_vec | Vector of Cartesian triplet indices where each patch starts. |
[in] | endIJK_vec | Vector of Cartesian triplet indices where each patch ends. Last cell part of each patch(lgr) will be {endIJK_vec[<patch-number>][0]-1, ..., endIJK_vec[<patch-number>][2]-1}. |
[in] | lgr_name_vec | Names (std::string) for the LGRs/levels. |
◆ beginCellCentroids()
Get an iterator over the cell centroids positioned at the first one.
◆ beginFaceCentroids()
Get an iterator over the face centroids positioned at the first one.
◆ boundaryId()
int Dune::CpGrid::boundaryId |
( |
int |
face | ) |
const |
◆ cellCenterDepth()
double Dune::CpGrid::cellCenterDepth |
( |
int |
cell_index | ) |
const |
Get vertical position of cell center ("zcorn" average).
cell_index The index of the specific cell.
◆ cellCentroid()
const Vector & Dune::CpGrid::cellCentroid |
( |
int |
cell | ) |
const |
Get the coordinates of the center of a cell.
- Parameters
-
cell | The index identifying the face. |
◆ cellCommunication()
Get the owner-overlap-copy communication for cells.
Suitable e.g. for parallel linear algebra used by CCFV
◆ cellFace()
int Dune::CpGrid::cellFace |
( |
int |
cell, |
|
|
int |
local_index, |
|
|
int |
level = -1 |
|
) |
| const |
Get a specific face of a cell.
- Parameters
-
cell | The index identifying the cell. |
local_index | The local index (in [0,numFaces(cell))) of the face in this cell. |
level | Integer representing the level grid to be considered. Default leaf grid view (current_view_data_) set to -1. |
- Returns
- The index identifying the face.
◆ cellFaceRow()
Get a list of indices identifying all faces of a cell.
- Parameters
-
cell | The index identifying the cell. |
◆ cellScatterGatherInterface()
const InterfaceMap & Dune::CpGrid::cellScatterGatherInterface |
( |
| ) |
const |
Get an interface for gathering/scattering data attached to cells with communication.
Scattering means sending data from the indices of the global grid on process 0 to the distributed grid on all ranks independent of the grid. Gathering is the other way around. The interface can be used with VariableSizeCommunicator and a custom index based data handle to scatter (forward direction of the communicator) and gather data (backward direction of the communicator). Here is a small example that prints the received values when scattering: struct Handle{
typedef int DataType;
const std::vector<int>& vals;
bool fixedsize() { return true; }
size_t size(std::size_t) { return 1; }
void gather(auto& B buf, size_t i)[ buf.write(vals[i]); }
void scatter(auto& B buf, size_t i, std::size_t) {
int val;
buf.read(val);
cout<<i<<": "<<val<<" "; }
};
Handle handle;
handle.vals.resize(grid.size(0), -1);
Dune::VariableSizeCommunicator<> comm(grid.comm(),
grid.cellScatterGatherInterface());
const CpGridTraits::Communication & comm() const Get the collective communication object.
int size(int level, int codim) const Number of grid entities per level and codim.
Referenced by scatterData().
◆ cellVolume()
double Dune::CpGrid::cellVolume |
( |
int |
cell | ) |
const |
Get the volume of the cell.
- Parameters
-
cell | The index identifying the cell. |
◆ comm()
Get the collective communication object.
◆ communicate() [1/2]
template<class DataHandle >
void Dune::CpGrid::communicate |
( |
DataHandle & |
data, |
|
|
InterfaceType |
iftype, |
|
|
CommunicationDirection |
dir |
|
) |
| const |
communicate objects for all codims on a given level.
The new communication interface. - Template Parameters
-
DataHandle | The type of the data handle describing the data. |
- Parameters
-
data | The data handle describing the data. Has to adhere to the Dune::DataHandleIF interface. |
iftype | The interface to use for the communication. |
dir | The direction of the communication along the interface (forward or backward). |
level | discarded as CpGrid is not adaptive. |
References Dune::cpgrid::CpGridData::communicate(), and data.
◆ communicate() [2/2]
template<class DataHandle >
void Dune::CpGrid::communicate |
( |
DataHandle & |
data, |
|
|
InterfaceType |
iftype, |
|
|
CommunicationDirection |
dir, |
|
|
int |
|
|
) |
| const |
|
inline |
communicate objects for all codims on a given level
The new communication interface. - Parameters
-
data | The data handle describing the data. Has to adhere to the Dune::DataHandleIF interface. |
iftype | The interface to use for the communication. |
dir | The direction of the communication along the interface (forward or backward). |
level | discarded as CpGrid is not adaptive. |
References communicate(), and data.
Referenced by communicate().
◆ createCartesian()
void Dune::CpGrid::createCartesian |
( |
const std::array< int, 3 > & |
dims, |
|
|
const std::array< double, 3 > & |
cellsize, |
|
|
const std::array< int, 3 > & |
shift = {0, 0, 0} |
|
) |
| |
Create a cartesian grid. - Parameters
-
dims | the number of cells in each cartesian direction. |
cellsize | the size of each cell in each dimension. |
shift | The origin of the grid, i.e. the corner of the cell with index (0,0,0) where the left, bottom, and top face of that cell intersect, is at the coordinate origin per default. This parameter shifts that corner to lie at (shift[0]*cellsize[0], ..., shift[2]*cellsize[2]). |
◆ currentData() [1/2]
Returns either data_ or distributed_data_(if non empty).
◆ currentData() [2/2]
◆ entity()
given an EntitySeed (or EntityPointer) return an entity object
◆ faceArea()
double Dune::CpGrid::faceArea |
( |
int |
face | ) |
const |
Get the area of a face.
- Parameters
-
cell | The index identifying the face. |
◆ faceAreaNormalEcl()
const Vector Dune::CpGrid::faceAreaNormalEcl |
( |
int |
face | ) |
const |
◆ faceCell()
int Dune::CpGrid::faceCell |
( |
int |
face, |
|
|
int |
local_index, |
|
|
int |
level = -1 |
|
) |
| const |
Get the index identifying a cell attached to a face.
Note that a face here is always oriented. If there are two neighboring cells then the orientation will be from local_index 0 to local_index 1 - Parameters
-
face | The index identifying the face. |
local_index | The local_index of the cell. |
level | Integer representing the level grid to be considered. Default leaf grid view (current_view_data_) set to -1. |
- Returns
- The index identifying a cell or -1 if there is no such cell due the face being part of the grid boundary or the cell being stored on another process.
Referenced by Dune::cpgrid::FaceCellsContainerProxy::operator()(), and Dune::cpgrid::FaceCellsProxy::operator[]().
◆ faceCenterEcl()
◆ faceCentroid()
const Vector & Dune::CpGrid::faceCentroid |
( |
int |
face | ) |
const |
Get the coordinates of the center of a face.
- Parameters
-
cell | The index identifying the face. |
◆ faceNormal()
const Vector & Dune::CpGrid::faceNormal |
( |
int |
face | ) |
const |
Get the unit normal of a face.
- Parameters
-
cell | The index identifying the face. |
- See also
- faceCell
◆ faceTag()
template<class Cell2FacesRowIterator >
int Dune::CpGrid::faceTag |
( |
const Cell2FacesRowIterator & |
cell_face | ) |
const |
Get the cartesian tag associated with a face tag.
The tag tells us in which direction the face would point in the underlying cartesian grid. - Parameters
-
An | iterator that points to the face and was obtained by iterating over Opm::UgGridHelpers::cell2Faces(grid). |
References I_FACE, J_FACE, K_FACE, NNC_FACE, numCells(), and numFaces().
◆ faceVertex()
int Dune::CpGrid::faceVertex |
( |
int |
face, |
|
|
int |
local_index |
|
) |
| const |
Get the index identifying a vertex of a face.
- Parameters
-
cell | The index identifying the face. |
local_index | The local_index (in [0,numFaceVertices(vertex) - 1]]) of the vertex. |
◆ gatherData()
template<class DataHandle >
void Dune::CpGrid::gatherData |
( |
DataHandle & |
handle | ) |
const |
Moves data from the distributed view to the global (all data on process) view.
- Template Parameters
-
DataHandle | The type of the data handle describing the data and responsible for gathering and scattering the data. |
- Parameters
-
handle | The data handle describing the data and responsible for gathering and scattering the data. |
◆ geomTypes()
const std::vector< Dune::GeometryType > & Dune::CpGrid::geomTypes |
( |
const int |
| ) |
const |
◆ getCellIndexSet() [1/2]
◆ getCellIndexSet() [2/2]
◆ getCellRemoteIndices() [1/2]
◆ getCellRemoteIndices() [2/2]
const RemoteIndices & Dune::CpGrid::getCellRemoteIndices |
( |
| ) |
const |
◆ getEclCentroid() [1/2]
std::array< double, 3 > Dune::CpGrid::getEclCentroid |
( |
const cpgrid::Entity< 0 > & |
elem | ) |
const |
◆ getEclCentroid() [2/2]
std::array< double, 3 > Dune::CpGrid::getEclCentroid |
( |
const int & |
idx | ) |
const |
◆ getIJK()
void Dune::CpGrid::getIJK |
( |
const int |
c, |
|
|
std::array< int, 3 > & |
ijk |
|
) |
| const |
◆ getLgrNameToLevel()
const std::map< std::string, int > & Dune::CpGrid::getLgrNameToLevel |
( |
| ) |
const |
◆ getMark()
int Dune::CpGrid::getMark |
( |
const cpgrid::Entity< 0 > & |
element | ) |
const |
Return refinement mark for entity.
- Returns
- refinement mark (1,0,-1) Currently, only 1 (refinement), or 0 (doing nothing).
◆ getParentIntersectionFromLgrBoundaryFace()
◆ ghostSize() [1/2]
unsigned int Dune::CpGrid::ghostSize |
( |
int |
| ) |
const |
Size of the ghost cell layer on the leaf level.
◆ ghostSize() [2/2]
unsigned int Dune::CpGrid::ghostSize |
( |
int |
, |
|
|
int |
|
|
) |
| const |
Size of the ghost cell layer on a given level.
◆ globalCell()
const std::vector< int > & Dune::CpGrid::globalCell |
( |
| ) |
const |
◆ globalIdSet()
const Traits::GlobalIdSet & Dune::CpGrid::globalIdSet |
( |
| ) |
const |
Access to the GlobalIdSet.
◆ globalRefine()
void Dune::CpGrid::globalRefine |
( |
int |
refCount | ) |
|
Refine the grid refCount times using the default refinement rule. This behaves like marking all elements for refinement and then calling preAdapt, adapt and postAdapt. The state after globalRefine is comparable to the state after postAdapt.
◆ lbegin() [1/2]
template<int codim>
Traits::template Codim< codim >::LevelIterator Dune::CpGrid::lbegin |
( |
int |
level | ) |
const |
Iterator to first entity of given codim on level.
◆ lbegin() [2/2]
template<int codim, PartitionIteratorType PiType>
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator Dune::CpGrid::lbegin |
( |
int |
level | ) |
const |
Iterator to first entity of given codim on level and PartitionIteratorType.
◆ leafbegin() [1/2]
template<int codim>
Traits::template Codim< codim >::LeafIterator Dune::CpGrid::leafbegin |
( |
| ) |
const |
Iterator to first leaf entity of given codim.
◆ leafbegin() [2/2]
template<int codim, PartitionIteratorType PiType>
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator Dune::CpGrid::leafbegin |
( |
| ) |
const |
Iterator to first leaf entity of given codim and PartitionIteratorType.
◆ leafend() [1/2]
template<int codim>
Traits::template Codim< codim >::LeafIterator Dune::CpGrid::leafend |
( |
| ) |
const |
one past the end of the sequence of leaf entities
◆ leafend() [2/2]
template<int codim, PartitionIteratorType PiType>
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator Dune::CpGrid::leafend |
( |
| ) |
const |
one past the end of the sequence of leaf entities
◆ leafIndexSet()
const Traits::LeafIndexSet & Dune::CpGrid::leafIndexSet |
( |
| ) |
const |
Access to the LeafIndexSet.
◆ lend() [1/2]
template<int codim>
Traits::template Codim< codim >::LevelIterator Dune::CpGrid::lend |
( |
int |
level | ) |
const |
one past the end on this level
◆ lend() [2/2]
template<int codim, PartitionIteratorType PiType>
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator Dune::CpGrid::lend |
( |
int |
level | ) |
const |
one past the end on this level
◆ levelIndexSet()
const Traits::LevelIndexSet & Dune::CpGrid::levelIndexSet |
( |
int |
level | ) |
const |
Access to the LevelIndexSets.
◆ loadBalance() [1/9]
std::pair< bool, std::vector< std::pair< std::string, bool > > > Dune::CpGrid::loadBalance |
( |
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections = {} , |
|
|
const double * |
transmissibilities = nullptr , |
|
|
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltanGoG , |
|
|
int |
level = -1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
This will construct the corresponding graph to the grid and use the transmissibilities specified as weights associated with its edges. The graph will be passed to the load balancer. - Parameters
-
wells | The wells of the eclipse If null wells will be neglected. If this is not null then complete well information of of the last scheduler step of the eclipse state will be used to make sure that all the possible completion cells of each well are stored on one process. This done by adding an edge with a very high edge weight for all possible pairs of cells in the completion set of a well. |
possibleFutureConnections | An optional unordered_map<string, set<array<int,3>>> containing possible future connections that might be opened during an ACTIONX. The fist entry is the name of the well and the second entry is a set containing the cartesian coordinates of the grid cells that get perforated of a possible future connection. The possible future connections are handed over to the grid partitioner to make sure these will be no the same partition when partitioning the grid. |
transmissibilities | The transmissibilities used as the edge weights. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
level | Level grid to be distributed. Integer between 0,..., maxLevel(). Defualt value set to -1, representing the leaf grid view. |
- Warning
- Throw if level>0. Currently, for CpGrid with LGRs, only distributing level zero grid is supported.
-
May only be called once.
- Returns
- A pair consisting of a boolean indicating whether loadbalancing actually happened and a vector containing a pair of name and a boolean, indicating whether this well has perforated cells local to the process, for all wells (sorted by name)
◆ loadBalance() [2/9]
bool Dune::CpGrid::loadBalance |
( |
const std::vector< int > & |
parts, |
|
|
bool |
ownersFirst = false , |
|
|
bool |
addCornerCells = false , |
|
|
int |
overlapLayers = 1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
- Parameters
-
parts | The partitioning information. For a cell with local index i the entry parts[i] is the partion number. Partition numbers need to start with zero and need to be consectutive also parts.size()==grid.leafGridView().size() and the ranks communicator need to be able to map all parts. Needs to valid at rank 0. Number of parts cannot exceed the number of ranks. Parts need to numbered consecutively starting from zero. |
ownersFirst | Order owner cells before copy/overlap cells. |
addCornerCells | Add corner cells to the overlap layer. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
- Warning
- May only be called once.
References Dune::defaultTransEdgeWgt, and Dune::simple.
◆ loadBalance() [3/9]
template<class DataHandle >
std::pair< bool, std::vector< std::pair< std::string, bool > > > Dune::CpGrid::loadBalance |
( |
DataHandle & |
data, |
|
|
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections = {} , |
|
|
const double * |
transmissibilities = nullptr , |
|
|
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = 1 , |
|
|
int |
level = -1 |
|
) |
| |
|
inline |
Distributes this grid and data over the available nodes in a distributed machine.
- Parameters
-
data | A data handle describing how to distribute attached data. |
wells | The wells of the eclipse Default: null If this is not null then complete well information of of the last scheduler step of the eclipse state will be used to make sure that all the possible completion cells of each well are stored on one process. This done by adding an edge with a very high edge weight for all possible pairs of cells in the completion set of a well. |
possibleFutureConnections | An optional unordered_map<string, set<array<int,3>>> containing possible future connections that might be opened during an ACTIONX. The fist entry is the name of the well and the second entry is a set containing the cartesian coordinates of the grid cells that get perforated of a possible future connection. The possible future connections are handed over to the grid partitioner to make sure these will be no the same partition when partitioning the grid. |
transmissibilities | The transmissibilities used to calculate the edge weights. |
overlapLayers | The number of layers of overlap cells to be added (default: 1) |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
level | Level grid to be distributed. Integer between 0,..., maxLevel(). Defualt value set to -1, representing the leaf grid view. |
- Warning
- Throw if level>0. Currently, for CpGrid with LGRs, only distributing level zero grid is supported.
- Template Parameters
-
DataHandle | The type implementing DUNE's DataHandle interface. |
- Warning
- May only be called once.
- Returns
- A pair consisting of a boolean indicating whether loadbalancing actually happened and a vector containing a pair of name and a boolean, indicating whether this well has perforated cells local to the process, for all wells (sorted by name)
◆ loadBalance() [4/9]
template<class DataHandle >
bool Dune::CpGrid::loadBalance |
( |
DataHandle & |
data, |
|
|
const std::vector< int > & |
parts, |
|
|
bool |
ownersFirst = false , |
|
|
bool |
addCornerCells = false , |
|
|
int |
overlapLayers = 1 |
|
) |
| |
|
inline |
Distributes this grid and data over the available nodes in a distributed machine.
- Parameters
-
data | A data handle describing how to distribute attached data. |
parts | The partitioning information. For a cell with local index i the entry parts[i] is the partion number. Partition numbers need to start with zero and need to be consectutive also parts.size()==grid.leafGridView().size() and the ranks communicator need to be able to map all parts. Needs to valid at rank 0. Number of parts cannot exceed the number of ranks. Parts need to numbered consecutively starting from zero. |
ownersFirst | Order owner cells before copy/overlap cells. |
addCornerCells | Add corner cells to the overlap layer. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
- Warning
- May only be called once.
References data, loadBalance(), and scatterData().
◆ loadBalance() [5/9]
template<class DataHandle >
std::pair< bool, std::vector< std::pair< std::string, bool > > > Dune::CpGrid::loadBalance |
( |
DataHandle & |
data, |
|
|
const std::vector< int > & |
parts, |
|
|
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections = {} , |
|
|
bool |
ownersFirst = false , |
|
|
bool |
addCornerCells = false , |
|
|
int |
overlapLayers = 1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
- Parameters
-
data | A data handle describing how to distribute attached data. |
parts | The partitioning information. For a cell with local index i the entry parts[i] is the partion number. Partition numbers need to start with zero and need to be consectutive also parts.size()==grid.leafGridView().size() and the ranks communicator need to be able to map all parts. Needs to valid at rank 0. Number of parts cannot exceed the number of ranks. Parts need to numbered consecutively starting from zero. |
ownersFirst | Order owner cells before copy/overlap cells. |
addCornerCells | Add corner cells to the overlap layer. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
- Template Parameters
-
DataHandle | The type implementing DUNE's DataHandle interface. |
- Warning
- May only be called once.
- Returns
- A pair consisting of a boolean indicating whether loadbalancing actually happened and a vector containing a pair of name and a boolean, indicating whether this well has perforated cells local to the process, for all wells (sorted by name)
◆ loadBalance() [6/9]
template<class DataHandle >
bool Dune::CpGrid::loadBalance |
( |
DataHandle & |
data, |
|
|
decltype(data.fixedSize(0, 0)) |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltan |
|
) |
| |
|
inline |
Distributes this grid and data over the available nodes in a distributed machine.
- Parameters
-
data | A data handle describing how to distribute attached data. |
overlapLayers | The number of layers of overlap cells to be added (default: 1) |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
- Template Parameters
-
DataHandle | The type implementing DUNE's DataHandle interface. |
- Warning
- May only be called once.
References data, loadBalance(), and scatterData().
◆ loadBalance() [7/9]
template<class DataHandle >
std::pair< bool, std::vector< std::pair< std::string, bool > > > Dune::CpGrid::loadBalance |
( |
DataHandle & |
data, |
|
|
EdgeWeightMethod |
method, |
|
|
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections, |
|
|
bool |
serialPartitioning, |
|
|
const double * |
transmissibilities = nullptr , |
|
|
bool |
ownersFirst = false , |
|
|
bool |
addCornerCells = false , |
|
|
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltanGoG , |
|
|
double |
imbalanceTol = 1.1 , |
|
|
bool |
allowDistributedWells = false |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
This will construct the corresponding graph to the grid and use the transmissibilities specified to calculate the weights associated with its edges. The graph will be passed to the load balancer. - Parameters
-
data | A data handle describing how to distribute attached data. |
method | The edge-weighting method to be used on the graph partitioner. |
wells | The information about all possible wells. If null then the wells will be neglected. Otherwise the wells will be used to make sure that all the possible completion cells of each well are stored on one process. This is done by adding an edge with a very high edge weight for all possible pairs of cells in the completion set of a well. |
possibleFutureConnections | An optional unordered_map<string, set<array<int,3>>> containing possible future connections that might be opened during an ACTIONX. The fist entry is the name of the well and the second entry is a set containing the cartesian coordinates of the grid cells that get perforated of a possible future connection. The possible future connections are handed over to the grid partitioner to make sure these will be no the same partition when partitioning the grid. |
serialPartitioning | If true, the partitioning will be done on a single process. |
transmissibilities | The transmissibilities used to calculate the edge weights. |
ownersFirst | Order owner cells before copy/overlap cells. |
addCornerCells | Add corner cells to the overlap layer. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
imbalanceTol | Set the imbalance tolerance used by the partitioner |
allowDistributedWells | Allow the perforation of a well to be distributed to the interior region of multiple processes. |
- Template Parameters
-
DataHandle | The type implementing DUNE's DataHandle interface. |
- Warning
- May only be called once.
- Returns
- A pair consisting of a boolean indicating whether loadbalancing actually happened and a vector containing a pair of name and a boolean, indicating whether this well has perforated cells local to the process, for all wells (sorted by name)
References data, and scatterData().
◆ loadBalance() [8/9]
std::pair< bool, std::vector< std::pair< std::string, bool > > > Dune::CpGrid::loadBalance |
( |
EdgeWeightMethod |
method, |
|
|
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections = {} , |
|
|
const double * |
transmissibilities = nullptr , |
|
|
bool |
ownersFirst = false , |
|
|
bool |
addCornerCells = false , |
|
|
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltanGoG , |
|
|
double |
imbalanceTol = 1.1 , |
|
|
int |
level = -1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
This will construct the corresponding graph to the grid and use the transmissibilities specified to calculate the weights associated with its edges. The graph will be passed to the load balancer. - Parameters
-
method | The edge-weighting method to be used on the graph partitioner. |
wells | The wells of the eclipse If null wells will be neglected. If this is not null then complete well information of of the last scheduler step of the eclipse state will be used to make sure that all the possible completion cells of each well are stored on one process. This done by adding an edge with a very high edge weight for all possible pairs of cells in the completion set of a well. |
possibleFutureConnections | An optional unordered_map<string, set<array<int,3>>> containing possible future connections that might be opened during an ACTIONX. The fist entry is the name of the well and the second entry is a set containing the cartesian coordinates of the grid cells that get perforated of a possible future connection. The possible future connections are handed over to the grid partitioner to make sure these will be no the same partition when partitioning the grid. |
transmissibilities | The transmissibilities used to calculate the edge weights. |
ownersFirst | Order owner cells before copy/overlap cells. |
addCornerCells | Add corner cells to the overlap layer. |
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
imbalanceTol | |
level | Level grid to be distributed. Integer between 0,..., maxLevel(). Defualt value set to -1, representing the leaf grid view. |
- Warning
- Throw if level>0. Currently, for CpGrid with LGRs, only distributing level zero grid is supported.
-
May only be called once.
- Returns
- A pair consisting of a boolean indicating whether loadbalancing actually happened and a vector containing a pair of name and a boolean, indicating whether this well has perforated cells local to the process, for all wells (sorted by name)
◆ loadBalance() [9/9]
bool Dune::CpGrid::loadBalance |
( |
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltan , |
|
|
double |
imbalanceTol = 1.1 , |
|
|
int |
level = -1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
- Parameters
-
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
imbalanceTol | |
level | Level grid to be distributed. Integer between 0,..., maxLevel(). Defualt value set to -1, representing the leaf grid view. |
- Warning
- Throw if level>0. Currently, for CpGrid with LGRs, only distributing level zero grid is supported.
-
May only be called once.
References Dune::defaultTransEdgeWgt.
Referenced by loadBalance().
◆ loadBalanceSerial()
bool Dune::CpGrid::loadBalanceSerial |
( |
int |
overlapLayers = 1 , |
|
|
int |
partitionMethod = Dune::PartitionMethod::zoltan , |
|
|
int |
edgeWeightMethod = Dune::EdgeWeightMethod::defaultTransEdgeWgt , |
|
|
double |
imbalanceTol = 1.1 , |
|
|
int |
level = -1 |
|
) |
| |
|
inline |
Distributes this grid over the available nodes in a distributed machine.
- Parameters
-
overlapLayers | The number of layers of cells of the overlap region (default: 1). |
partitionMethod | The method used to partition the grid, one of Dune::PartitionMethod |
edgeWeightMethod | The edge-weighting method to be used on the graph partitioner. |
imbalanceTol | |
level | Level grid to be distributed. Integer between 0,..., maxLevel(). Defualt value set to -1, representing the leaf grid view. |
- Warning
- Throw if level>0. Currently, for CpGrid with LGRs, only distributing level zero grid is supported.
-
May only be called once.
◆ localIdSet()
const Traits::LocalIdSet & Dune::CpGrid::localIdSet |
( |
| ) |
const |
Access to the LocalIdSet.
◆ logicalCartesianSize()
const std::array< int, 3 > & Dune::CpGrid::logicalCartesianSize |
( |
| ) |
const |
The logical cartesian size of the global grid. This function is not part of the Dune grid interface, and should be used with caution.
◆ mapLeafIndexSetToLocalCartesianIndexSets()
std::vector< std::array< int, 2 > > Dune::CpGrid::mapLeafIndexSetToLocalCartesianIndexSets |
( |
| ) |
const |
Reverse map: from leaf index cell to { level, local/level Cartesian index of the cell }.
◆ mapLocalCartesianIndexSetsToLeafIndexSet()
std::vector< std::unordered_map< std::size_t, std::size_t > > Dune::CpGrid::mapLocalCartesianIndexSetsToLeafIndexSet |
( |
| ) |
const |
Compute for each level grid, a map from the global_cell_[ cell index in level grid ] to the leaf index of the equivalent cell on the leaf grid view. Notice that cells that vanished and do not appear on the leaf grid view will not be considered. global_cell_[ cell index in level grid ] coincide with (local) Cartesian Index.
◆ mark()
bool Dune::CpGrid::mark |
( |
int |
refCount, |
|
|
const cpgrid::Entity< 0 > & |
element |
|
) |
| |
Mark entity for refinement (or coarsening).
------------— Adaptivity (begin) ------------—
Refinement on CpGrid is partially supported for Cartesian grids, with the keyword CARFIN. Nested refinement is not supported yet, so the so-called "host grid" is defined by default equal to the GLOBAL Grid (level zero). Therefore, we mark elements (from the GLOBAL grid) for refinement. This only works for entities of codim 0. In distributed grids, element markings are synchronized across processes. If an element is marked by one process, all other processes that share that element will also receive and apply the same mark.
- Parameters
-
[in] | refCount | To mark the element for
- refinement, refCount == 1
- doing nothing, refCount == 0
- coarsening, refCount == -1 (not applicable yet)
|
[in] | element | Entity<0>. Currently, an element from the GLOBAL grid (level zero). |
- Returns
- true, if marking was succesfull. false, if marking was not possible.
◆ maxLevel()
int Dune::CpGrid::maxLevel |
( |
| ) |
const |
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview), 0 = the coarsest level.
◆ name()
std::string Dune::CpGrid::name |
( |
| ) |
const |
Get the grid name.
It's the same as the class name. What did you expect, something funny?
◆ numBoundarySegments()
unsigned int Dune::CpGrid::numBoundarySegments |
( |
| ) |
const |
returns the number of boundary segments within the macro grid
◆ numCellFaces() [1/2]
int Dune::CpGrid::numCellFaces |
( |
| ) |
const |
Get the sum of all faces attached to all cells.
Each face identified by a unique index is counted as often as there are neigboring cells attached to it. - See also
- numCellFaces(int)const
◆ numCellFaces() [2/2]
int Dune::CpGrid::numCellFaces |
( |
int |
cell, |
|
|
int |
level = -1 |
|
) |
| const |
Get the number of faces of a cell.
Due to faults, and collapsing vertices (along pillars) this number is quite arbitrary. Its lower bound is 4, but there is no upper bound. - Parameters
-
cell | the index identifying the cell. |
level | Integer representing the level grid to be considered. Default leaf grid view (current_view_data_) set to -1. |
◆ numCells()
int Dune::CpGrid::numCells |
( |
int |
level = -1 | ) |
const |
Get the number of cells.
- Parameters
-
level | Integer representing the level grid to be considered. Default leaf grid view (current_view_data_) set to -1. |
Referenced by faceTag().
◆ numFaces()
int Dune::CpGrid::numFaces |
( |
int |
level = -1 | ) |
const |
Get the number of faces.
- Parameters
-
level | Integer representing the level grid to be considered. Default leaf grid view (current_view_data_) set to -1. |
Referenced by faceTag().
◆ numFaceVertices()
int Dune::CpGrid::numFaceVertices |
( |
int |
face | ) |
const |
◆ numVertices()
int Dune::CpGrid::numVertices |
( |
| ) |
const |
Get The number of vertices.
◆ overlapSize() [1/2]
unsigned int Dune::CpGrid::overlapSize |
( |
int |
| ) |
const |
Size of the overlap on the leaf level.
◆ overlapSize() [2/2]
unsigned int Dune::CpGrid::overlapSize |
( |
int |
, |
|
|
int |
|
|
) |
| const |
Size of the overlap on a given level.
◆ pointScatterGatherInterface()
const InterfaceMap & Dune::CpGrid::pointScatterGatherInterface |
( |
| ) |
const |
◆ postAdapt()
void Dune::CpGrid::postAdapt |
( |
| ) |
|
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
◆ preAdapt()
bool Dune::CpGrid::preAdapt |
( |
| ) |
|
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called after elements have been marked for refinement.
◆ processEclipseFormat()
void Dune::CpGrid::processEclipseFormat |
( |
const grdecl & |
input_data, |
|
|
bool |
remove_ij_boundary, |
|
|
bool |
turn_normals = false |
|
) |
| |
Read the Eclipse grid format ('grdecl'). - Parameters
-
input_data | the data in grdecl format, declared in preprocess.h. |
remove_ij_boundary | if true, will remove (i, j) boundaries. Used internally. |
◆ replaceLgr1CornerIdxByLgr2CornerIdx() [1/2]
int Dune::CpGrid::replaceLgr1CornerIdxByLgr2CornerIdx |
( |
const std::array< int, 3 > & |
cells_per_dim_lgr1, |
|
|
int |
cornerIdxLgr1, |
|
|
const std::array< int, 3 > & |
cells_per_dim_lgr2 |
|
) |
| const |
|
protected |
A refined corner appears in two single-cell-refinements. Given the corner index in the first single-cell-refinement, compute the corner index in the neighboring single-cell-refinement.
- Parameters
-
[in] | cells_per_dim_lgr1 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr1 single-cell-refinement. |
[in] | cornerIdxInLgr1 | Corner index in the elemLgr1 single-cell-refinement. |
[in] | cells_per_dim_lgr2 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr2 single-cell-refinement. |
◆ replaceLgr1CornerIdxByLgr2CornerIdx() [2/2]
int Dune::CpGrid::replaceLgr1CornerIdxByLgr2CornerIdx |
( |
const std::array< int, 3 > & |
cells_per_dim_lgr1, |
|
|
int |
cornerIdxLgr1, |
|
|
int |
elemLgr1, |
|
|
int |
parentFaceLastAppearanceIdx, |
|
|
const std::array< int, 3 > & |
cells_per_dim_lgr2 |
|
) |
| const |
|
protected |
A new refined corner lays on an edge and appears in at least two single-cell-refinements. Given the corner index in one single-cell-refinement, compute the corner index in a neighboring single-cell-refinement.
- Parameters
-
[in] | cells_per_dim_lgr1 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr1 single-cell-refinement. |
[in] | cornerIdxInLgr1 | Corner index in the elemLgr1 single-cell-refinement. |
[in] | parentFaceLastAppearanceIdx | Parent face index where the refined corner appears for last time. |
[in] | cells_per_dim_lgr2 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr2 single-cell-refinement. |
◆ replaceLgr1FaceIdxByLgr2FaceIdx()
int Dune::CpGrid::replaceLgr1FaceIdxByLgr2FaceIdx |
( |
const std::array< int, 3 > & |
cells_per_dim_lgr1, |
|
|
int |
faceIdxInLgr1, |
|
|
const std::shared_ptr< cpgrid::CpGridData > & |
elemLgr1_ptr, |
|
|
const std::array< int, 3 > & |
cells_per_dim_lgr2 |
|
) |
| const |
|
protected |
A new refined face lays on the boudndary of a single-cell-refinement appears in at most two single-cell-refinements. Given the face index in one single-cell-refinement, compute the face index in a neighboring single-cell-refinement.
- Parameters
-
[in] | cells_per_dim_lgr1 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr1 single-cell-refinement. |
[in] | faceIdxInLgr1 | Face index in the elemLgr1 single-cell-refinement. |
[in] | elemLgr1_ptr | Pointer to the elemLgr1 single-cell-refinement grid. |
[in] | cells_per_dim_lgr2 | Total children cells in each direction (x-,y-, and z-direction) of the elemLgr2 single-cell-refinement. |
◆ scatterData()
template<class DataHandle >
void Dune::CpGrid::scatterData |
( |
DataHandle & |
handle | ) |
const |
Moves data from the global (all data on process) view to the distributed view.
This method does not do communication but assumes that the global grid is present on every process and simply copies data to the distributed view. - Template Parameters
-
DataHandle | The type of the data handle describing the data and responsible for gathering and scattering the data. |
- Parameters
-
handle | The data handle describing the data and responsible for gathering and scattering the data. |
References cellScatterGatherInterface(), and pointScatterGatherInterface().
Referenced by loadBalance().
◆ setPartitioningParams()
void Dune::CpGrid::setPartitioningParams |
( |
const std::map< std::string, std::string > & |
params | ) |
|
◆ setUniqueBoundaryIds()
void Dune::CpGrid::setUniqueBoundaryIds |
( |
bool |
uids | ) |
|
Set whether we want to have unique boundary ids. - Parameters
-
uids | if true, each boundary intersection will have a unique boundary id. |
◆ size() [1/4]
int Dune::CpGrid::size |
( |
GeometryType |
type | ) |
const |
number of leaf entities per geometry type in this process
◆ size() [2/4]
int Dune::CpGrid::size |
( |
int |
codim | ) |
const |
number of leaf entities per codim in this process
◆ size() [3/4]
int Dune::CpGrid::size |
( |
int |
level, |
|
|
GeometryType |
type |
|
) |
| const |
number of entities per level and geometry type in this process
◆ size() [4/4]
int Dune::CpGrid::size |
( |
int |
level, |
|
|
int |
codim |
|
) |
| const |
Number of grid entities per level and codim.
◆ sortedNumAquiferCells()
const std::vector< int > & Dune::CpGrid::sortedNumAquiferCells |
( |
| ) |
const |
Get sorted active cell indices of numerical aquifer.
◆ switchToDistributedView()
void Dune::CpGrid::switchToDistributedView |
( |
| ) |
|
Switch to the distributed view.
◆ switchToGlobalView()
void Dune::CpGrid::switchToGlobalView |
( |
| ) |
|
Switch to the global view.
◆ syncDistributedGlobalCellIds()
void Dune::CpGrid::syncDistributedGlobalCellIds |
( |
| ) |
|
Synchronizes cell global ids across processes after load balancing.
LGRs (Local Grid Refinements) can be added either in the undistributed view first and then in the distributed view, or vice versa. This method ensures consistency by rewriting the global cell ids in the distributed view using the corresponding ids from the undistributed view.
◆ uniqueBoundaryIds()
bool Dune::CpGrid::uniqueBoundaryIds |
( |
| ) |
const |
Is the grid currently using unique boundary ids? - Returns
- true if each boundary intersection has a unique id false if we use the (default) 1-6 ids for i- i+ j- j+ k- k+ boundaries.
◆ vertexPosition()
const Vector & Dune::CpGrid::vertexPosition |
( |
int |
vertex | ) |
const |
Get the Position of a vertex.
- Parameters
-
cell | The index identifying the cell. |
- Returns
- The coordinates of the vertex.
◆ zcornData()
const std::vector< double > & Dune::CpGrid::zcornData |
( |
| ) |
const |
◆ zoltanPartitionWithoutScatter()
std::vector< int > Dune::CpGrid::zoltanPartitionWithoutScatter |
( |
const std::vector< cpgrid::OpmWellType > * |
wells, |
|
|
const std::unordered_map< std::string, std::set< int > > & |
possibleFutureConnections, |
|
|
const double * |
transmissibilities, |
|
|
const int |
numParts, |
|
|
const double |
imbalanceTol |
|
) |
| const |
Partitions the grid using Zoltan without decomposing and distributing it among processes.
- Parameters
-
wells | The wells of the eclipse. |
possibleFutureConnections | An optional unordered_map<string, set<array<int,3>>> containing possible future connections that might be opened during an ACTIONX. The fist entry is the name of the well and the second entry is a set containing the cartesian coordinates of the grid cells that get perforated of a possible future connection. The possible future connections are handed over to the grid partitioner to make sure these will be no the same partition when partitioning the grid. |
transmissibilities | The transmissibilities used to calculate the edge weights. |
numParts | Number of parts in the partition. |
- Returns
- An array with the domain index for each cell.
◆ cpgrid::CpGridData
◆ cpgrid::Entity< 0 >
◆ cpgrid::Entity< 1 >
◆ cpgrid::Entity< 2 >
◆ cpgrid::Entity< 3 >
◆ createEntity
The documentation for this class was generated from the following file:
|