38#ifndef OPM_CPGRID_HEADER
39#define OPM_CPGRID_HEADER
42#include <opm/grid/utility/platform_dependent/disable_warnings.h>
44#include <dune/grid/common/grid.hh>
45#include <opm/grid/cpgrid/CpGridDataTraits.hpp>
46#include <opm/grid/cpgrid/OrientedEntityTable.hpp>
48#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
49#include "common/GridEnums.hpp"
50#include <opm/grid/utility/OpmWellType.hpp>
71 template <
int>
class Entity;
72 template<
int,
int>
class Geometry;
73 class HierarchicIterator;
74 class IntersectionIterator;
75 template<
int, PartitionIteratorType>
class Iterator;
76 class LevelGlobalIdSet;
79 class IntersectionIterator;
87 const std::array<int, 3>&,
91 const std::array<int,3>&,
92 const std::array<int,3>&,
93 const std::array<int,3>&);
95void refinePatch_and_check(
const std::array<int,3>&,
96 const std::array<int,3>&,
97 const std::array<int,3>&);
153 template <PartitionIteratorType pitype>
165 template <PartitionIteratorType pitype>
171 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> >
LeafGridView;
178 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>>
LeafGridView;
191 using CollectiveCommunication = cpgrid::CpGridDataTraits::CollectiveCommunication;
213 :
public GridDefaultImplementation<3, 3, double, CpGridFamily>
224 const std::array<int,3>&,
228 const std::array<int,3>&,
229 const std::array<int,3>&,
230 const std::array<int,3>&);
232 void ::refinePatch_and_check(
const std::array<int,3>&,
233 const std::array<int,3>&,
234 const std::array<int,3>&);
290 Opm::EclipseState* ecl_state,
291 bool periodic_extension,
bool turn_normals,
bool clip_z,
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false);
340 const std::array<double, 3>& cellsize,
341 const std::array<int, 3>& shift = {0,0,0});
364 void getIJK(
const int c, std::array<int,3>& ijk)
const;
385 std::string
name()
const;
392 typename Traits::template Codim<codim>::LevelIterator
lbegin (
int level)
const;
395 typename Traits::template Codim<codim>::LevelIterator
lend (
int level)
const;
398 template<
int codim, PartitionIteratorType PiType>
399 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lbegin (
int level)
const;
401 template<
int codim, PartitionIteratorType PiType>
402 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lend (
int level)
const;
406 typename Traits::template Codim<codim>::LeafIterator
leafbegin()
const;
409 typename Traits::template Codim<codim>::LeafIterator
leafend()
const;
412 template<
int codim, PartitionIteratorType PiType>
413 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafbegin()
const;
415 template<
int codim, PartitionIteratorType PiType>
416 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafend()
const;
419 int size (
int level,
int codim)
const;
422 int size (
int codim)
const;
425 int size (
int level, GeometryType type)
const;
428 int size (GeometryType type)
const;
445 const std::vector<Dune::GeometryType>& geomTypes(
const int)
const;
463 void addLgrUpdateLeafView(
const std::array<int,3>& cells_per_dim,
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK);
526 void setZoltanParams(
const std::map<std::string,std::string>& params);
539 return get<0>(scatterGrid(
defaultTransEdgeWgt,
false,
nullptr,
false,
nullptr,
true, overlapLayers, useZoltan ));
563 std::pair<bool,std::vector<std::pair<std::string,bool>>>
565 const double* transmissibilities =
nullptr,
566 int overlapLayers=1,
bool useZoltan=
true)
568 return scatterGrid(
defaultTransEdgeWgt,
false, wells,
false, transmissibilities,
false, overlapLayers, useZoltan);
596 std::pair<bool,std::vector<std::pair<std::string,bool>>>
598 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
599 bool addCornerCells=
false,
int overlapLayers=1,
600 bool useZoltan =
true)
602 return scatterGrid(method, ownersFirst, wells,
false, transmissibilities, addCornerCells, overlapLayers, useZoltan);
624 template<
class DataHandle>
625 std::pair<bool, std::vector<std::pair<std::string,bool> > >
627 const std::vector<cpgrid::OpmWellType> * wells,
628 const double* transmissibilities =
nullptr,
629 int overlapLayers=1,
bool useZoltan =
true)
631 auto ret =
loadBalance(wells, transmissibilities, overlapLayers, useZoltan);
668 template<
class DataHandle>
669 std::pair<bool, std::vector<std::pair<std::string,bool> > >
671 const std::vector<cpgrid::OpmWellType> * wells,
672 bool serialPartitioning,
673 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
674 bool addCornerCells=
false,
int overlapLayers=1,
bool useZoltan =
true,
675 double zoltanImbalanceTol = 1.1,
676 bool allowDistributedWells =
false)
678 auto ret = scatterGrid(method, ownersFirst, wells, serialPartitioning, transmissibilities,
679 addCornerCells, overlapLayers, useZoltan, zoltanImbalanceTol, allowDistributedWells);
704 template<
class DataHandle>
705 std::pair<bool, std::vector<std::pair<std::string,bool> > >
707 const std::vector<cpgrid::OpmWellType> * wells,
708 bool ownersFirst=
false,
709 bool addCornerCells=
false,
int overlapLayers=1)
715 addCornerCells, overlapLayers,
false,
734 template<
class DataHandle>
736 decltype(data.fixedSize(0,0)) overlapLayers=1,
bool useZoltan =
true)
759 bool loadBalance(
const std::vector<int>& parts,
bool ownersFirst=
false,
760 bool addCornerCells=
false,
int overlapLayers=1)
766 addCornerCells, overlapLayers,
false,
783 template<
class DataHandle>
784 bool loadBalance(DataHandle& data,
const std::vector<int>& parts,
bool ownersFirst=
false,
785 bool addCornerCells=
false,
int overlapLayers=1)
787 bool ret =
loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
802 const double* transmissibilities,
804 const double zoltanImbalanceTol)
const;
813 template<
class DataHandle>
814 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int )
const
826 template<
class DataHandle>
827 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir)
const;
843 typedef Dune::FieldVector<double, 3> Vector;
846 const std::vector<double>& zcornData()
const;
872 int cellFace(
int cell,
int local_index)
const;
888 int faceCell(
int face,
int local_index)
const;
898 int numFaceVertices(
int face)
const;
904 int faceVertex(
int face,
int local_index)
const;
911 const Vector faceCenterEcl(
int cell_index,
int face)
const;
913 const Vector faceAreaNormalEcl(
int face)
const;
947 :
public RandomAccessIteratorFacade<CentroidIterator<codim>,
948 FieldVector<double, 3>,
949 const FieldVector<double, 3>&, int>
961 const FieldVector<double,3>& dereference()
const
963 return iter_->center();
969 const FieldVector<double,3>& elementAt(
int n)
971 return iter_[n]->center();
999 int boundaryId(
int face)
const;
1007 template<
class Cell2FacesRowIterator>
1009 faceTag(
const Cell2FacesRowIterator& cell_face)
const;
1031 template<
class DataHandle>
1040 template<
class DataHandle>
1089 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1091 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1094 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1099 const CommunicationType& cellCommunication()
const;
1101 ParallelIndexSet& getCellIndexSet();
1103 RemoteIndices& getCellRemoteIndices();
1105 const ParallelIndexSet& getCellIndexSet()
const;
1107 const RemoteIndices& getCellRemoteIndices()
const;
1139 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1142 const std::vector<cpgrid::OpmWellType> * wells,
1143 bool serialPartitioning,
1144 const double* transmissibilities,
1145 bool addCornerCells,
1147 bool useZoltan =
true,
1148 double zoltanImbalanceTol = 1.1,
1149 bool allowDistributedWells =
true,
1150 const std::vector<int>& input_cell_part = {});
1156 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1160 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1166 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1172 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1176 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1182 std::map<std::string,std::string> zoltanParams;
1188#include <opm/grid/cpgrid/Entity.hpp>
1189#include <opm/grid/cpgrid/Iterators.hpp>
1190#include <opm/grid/cpgrid/CpGridData.hpp>
1196 namespace Capabilities
1202 static const bool v =
true;
1209 static const bool v =
true;
1215 static const bool v =
true;
1221 static const bool v =
true;
1228 static const bool v =
false;
1233 template<
class DataHandle>
1236 current_view_data_->
communicate(data, iftype, dir);
1240 template<
class DataHandle>
1244 if(distributed_data_.empty())
1245 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balanced grid!");
1254 template<
class DataHandle>
1258 if(distributed_data_.empty())
1259 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balance grid!");
1260 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1268 template<
class Cell2FacesRowIterator>
1281 const int cell = cell_face.getCellIndex();
1282 const int face = *cell_face;
1283 assert (0 <= cell); assert (cell <
numCells());
1284 assert (0 <= face); assert (face <
numFaces());
1289 const F2C& f2c = current_view_data_->face_to_cell_[f];
1290 const face_tag tag = current_view_data_->face_tag_[f];
1292 assert ((f2c.size() == 1) || (f2c.size() == 2));
1294 int inside_cell = 0;
1296 if ( f2c.size() == 2 )
1298 if ( f2c[1].index() == cell )
1303 const bool normal_is_in = ! f2c[inside_cell].orientation();
1308 return normal_is_in ? 0 : 1;
1311 return normal_is_in ? 2 : 3;
1315 return normal_is_in ? 4 : 5;
1320 OPM_THROW(std::logic_error,
"Unhandled face tag. This should never happen!");
1329#include <opm/grid/cpgrid/PersistentContainer.hpp>
1330#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1331#include "cpgrid/Intersection.hpp"
1332#include "cpgrid/Geometry.hpp"
1333#include "cpgrid/Indexsets.hpp"
An iterator over the centroids of the geometry of the entities.
Definition: CpGrid.hpp:950
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition: CpGrid.hpp:957
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition: CpGrid.hpp:954
[ provides Dune::Grid ]
Definition: CpGrid.hpp:214
std::string name() const
Get the grid name.
Definition: CpGrid.cpp:584
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, bool useZoltan=true)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:597
void switchToGlobalView()
Switch to the global view.
Definition: CpGrid.cpp:1251
void readSintefLegacyFormat(const std::string &grid_prefix)
Read the Sintef legacy grid format ('topogeom').
Definition: CpGrid.cpp:1293
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition: CpGrid.cpp:923
CpGridFamily GridFamily
Family typedef, why is this not defined by Grid<>?
Definition: CpGrid.hpp:245
Traits::template Codim< codim >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
void gatherData(DataHandle &handle) const
Moves data from the distributed view to the global (all data on process) view.
Definition: CpGrid.hpp:1255
void globalRefine(int)
global refinement
Definition: CpGrid.cpp:881
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: CpGrid.cpp:864
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
Definition: CpGrid.cpp:1185
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: CpGrid.cpp:859
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition: CpGrid.hpp:1270
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition: CpGrid.cpp:579
void 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.
Definition: CpGrid.cpp:496
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
void processEclipseFormat(const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false)
Read the Eclipse grid format ('grdecl').
Definition: CpGrid.cpp:1338
unsigned int overlapSize(int) const
Size of the overlap on the leaf level.
Definition: CpGrid.cpp:901
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition: CpGrid.hpp:1044
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this level
const Vector & faceNormal(int face) const
Get the unit normal of a face.
Definition: CpGrid.cpp:1170
int numCellFaces() const
Get the sum of all faces attached to all cells.
Definition: CpGrid.cpp:1020
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition: CpGrid.cpp:564
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, const std::vector< int > &parts, const std::vector< cpgrid::OpmWellType > *wells, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:706
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Definition: CpGrid.cpp:1180
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition: CpGrid.cpp:589
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition: CpGrid.cpp:948
double faceArea(int face) const
Get the area of a face.
Definition: CpGrid.cpp:1160
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
Definition: CpGrid.cpp:1165
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: CpGrid.cpp:869
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition: CpGrid.cpp:835
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the global grid.
Definition: CpGrid.cpp:559
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
Definition: CpGrid.cpp:1030
bool loadBalance(int overlapLayers=1, bool useZoltan=true)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:536
int numVertices() const
Get The number of vertices.
Definition: CpGrid.cpp:969
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: CpGrid.cpp:569
int cellFace(int cell, int local_index) const
Get a specific face of a cell.
Definition: CpGrid.cpp:979
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition: CpGrid.cpp:574
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(const std::vector< cpgrid::OpmWellType > *wells, const double *transmissibilities=nullptr, int overlapLayers=1, bool useZoltan=true)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:564
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, EdgeWeightMethod method, const std::vector< cpgrid::OpmWellType > *wells, bool serialPartitioning, const double *transmissibilities=nullptr, bool ownersFirst=false, bool addCornerCells=false, int overlapLayers=1, bool useZoltan=true, double zoltanImbalanceTol=1.1, bool allowDistributedWells=false)
Distributes this grid over the available nodes in a distributed machine.
Definition: CpGrid.hpp:670
std::vector< int > zoltanPartitionWithoutScatter(const std::vector< cpgrid::OpmWellType > *wells, const double *transmissibilities, const int numParts, const double zoltanImbalanceTol) const
Partitions the grid using Zoltan without decomposing and distributing it among processes.
Definition: CpGrid.cpp:149
bool loadBalance(DataHandle &data, decltype(data.fixedSize(0, 0)) overlapLayers=1, bool useZoltan=true)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:735
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition: CpGrid.cpp:989
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
void writeSintefLegacyFormat(const std::string &grid_prefix) const
Write the Sintef legacy grid format ('topogeom').
Definition: CpGrid.cpp:1303
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition: CpGrid.cpp:907
std::pair< bool, std::vector< std::pair< std::string, bool > > > loadBalance(DataHandle &data, const std::vector< cpgrid::OpmWellType > *wells, const double *transmissibilities=nullptr, int overlapLayers=1, bool useZoltan=true)
Distributes this grid and data over the available nodes in a distributed machine.
Definition: CpGrid.hpp:626
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: CpGrid.cpp:876
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
const InterfaceMap & pointScatterGatherInterface() const
Get an interface for gathering/scattering data attached to points with communication.
Definition: CpGrid.cpp:1246
CpGrid()
Default constructor.
Definition: CpGrid.cpp:122
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level and PartitionIteratorType.
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition: CpGrid.cpp:1195
void switchToDistributedView()
Switch to the distributed view.
Definition: CpGrid.cpp:1256
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
Definition: CpGrid.cpp:1190
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.
Definition: CpGrid.hpp:784
const cpgrid::OrientedEntityTable< 0, 1 >::row_type cellFaceRow(int cell) const
Get a list of indices identifying all faces of a cell.
Definition: CpGrid.cpp:984
int numCells() const
Get the number of cells.
Definition: CpGrid.cpp:959
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition: CpGrid.cpp:1241
void addLgrUpdateLeafView(const std::array< int, 3 > &cells_per_dim, const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK)
Create a grid out of a coarse one and a refinement(LGR) of a selected block-shaped patch of cells fro...
Definition: CpGrid.cpp:1365
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &) const
given an EntitySeed (or EntityPointer) return an entity object
Definition: CpGrid.cpp:892
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.
Definition: CpGrid.hpp:759
int numFaces() const
Get the number of faces.
Definition: CpGrid.cpp:964
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
Definition: CpGrid.cpp:1155
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
Definition: CpGrid.cpp:1035
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition: CpGrid.hpp:814
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
double cellVolume(int cell) const
Get the volume of the cell.
Definition: CpGrid.cpp:1175
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition: CpGrid.hpp:1241
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:122
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:738
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
This class encapsulates geometry for vertices, intersections, and cells.
Definition: Geometry.hpp:74
The global id set for Dune.
Definition: Indexsets.hpp:304
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:108
Definition: Indexsets.hpp:55
Definition: Intersection.hpp:288
Definition: Intersection.hpp:66
Iterator intended to be used as LeafIterator and LevelIterator (no difference due to no adaptivity) f...
Definition: Iterators.hpp:60
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
Copyright 2019 Equinor AS.
Definition: CartesianIndexMapper.hpp:10
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's graph...
Definition: GridEnums.hpp:34
@ defaultTransEdgeWgt
Use the transmissibilities as edge weights.
Definition: GridEnums.hpp:38
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:29
Low-level corner-point processing routines and supporting data structures.
face_tag
Connection taxonomy.
Definition: preprocess.h:66
@ K_FACE
Connection topologically normal to I-J plane.
Definition: preprocess.h:69
@ J_FACE
Connection topologically normal to I-K plane.
Definition: preprocess.h:68
@ NNC_FACE
Arbitrary non-neighbouring connection.
Definition: preprocess.h:70
@ I_FACE
Connection topologically normal to J-K plane.
Definition: preprocess.h:67
Definition: CpGrid.hpp:201
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:155
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition: CpGrid.hpp:157
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition: CpGrid.hpp:159
Traits associated with a specific codim.
Definition: CpGrid.hpp:131
cpgrid::Entity< cd > Entity
The type of the entity.
Definition: CpGrid.hpp:140
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity.
Definition: CpGrid.hpp:134
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition: CpGrid.hpp:137
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition: CpGrid.hpp:146
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition: CpGrid.hpp:143
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition: CpGrid.hpp:149
Traits associated with a specific grid partition type.
Definition: CpGrid.hpp:167
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:171
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:169
Definition: CpGrid.hpp:111
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition: CpGrid.hpp:181
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition: CpGrid.hpp:120
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition: CpGrid.hpp:178
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition: CpGrid.hpp:176
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition: CpGrid.hpp:185
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition: CpGrid.hpp:122
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition: CpGrid.hpp:183
GlobalIdSet LocalIdSet
The type of the local id set.
Definition: CpGrid.hpp:187
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition: CpGrid.hpp:190
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition: CpGrid.hpp:125
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition: CpGrid.hpp:118
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition: CpGrid.hpp:116
CpGrid Grid
The type that implements the grid.
Definition: CpGrid.hpp:113
Raw corner-point specification of a particular geological model.
Definition: preprocess.h:56