39#ifndef OPM_CPGRID_HEADER
40#define OPM_CPGRID_HEADER
43#include <opm/grid/utility/platform_dependent/disable_warnings.h>
45#include <dune/grid/common/grid.hh>
46#include <opm/grid/cpgrid/CpGridDataTraits.hpp>
47#include <opm/grid/cpgrid/OrientedEntityTable.hpp>
49#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
50#include "common/GridEnums.hpp"
51#include <opm/grid/utility/OpmWellType.hpp>
62template<
typename Gr
id,
typename Gr
idView>
class LookUpData;
63template<
typename Gr
id,
typename Gr
idView>
class LookUpCartesianData;
70 template<
typename Gr
id,
typename Gr
idView>
class LookUpData;
76 template <
int>
class Entity;
77 template<
int,
int>
class Geometry;
78 class HierarchicIterator;
79 class IntersectionIterator;
80 template<
int, PartitionIteratorType>
class Iterator;
81 class LevelGlobalIdSet;
84 class IntersectionIterator;
92 const std::vector<std::array<int,3>>&,
93 const std::vector<std::array<int,3>>&);
98 const std::array<int, 3>&,
102 const std::vector<std::array<int,3>>&,
103 const std::vector<std::array<int,3>>&,
104 const std::vector<std::array<int,3>>&,
105 const std::vector<std::string>&);
107void refinePatch_and_check(
const std::array<int,3>&,
108 const std::array<int,3>&,
109 const std::array<int,3>&);
165 template <PartitionIteratorType pitype>
177 template <PartitionIteratorType pitype>
183 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> >
LeafGridView;
190 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>>
LeafGridView;
203 using CollectiveCommunication = cpgrid::CpGridDataTraits::CollectiveCommunication;
225 :
public GridDefaultImplementation<3, 3, double, CpGridFamily>
232 template<
typename Gr
id,
typename Gr
idView>
friend class Opm::LookUpData;
237 const std::vector<std::array<int,3>>&,
238 const std::vector<std::array<int,3>>&);
242 const std::array<int,3>&,
246 const std::vector<std::array<int,3>>&,
247 const std::vector<std::array<int,3>>&,
248 const std::vector<std::array<int,3>>&,
249 const std::vector<std::string>&);
251 void ::refinePatch_and_check(
const std::array<int,3>&,
252 const std::array<int,3>&,
253 const std::array<int,3>&);
309 Opm::EclipseState* ecl_state,
310 bool periodic_extension,
bool turn_normals,
bool clip_z,
333 Opm::EclipseState* ecl_state,
334 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false);
359 const std::array<double, 3>& cellsize,
360 const std::array<int, 3>& shift = {0,0,0});
383 void getIJK(
const int c, std::array<int,3>& ijk)
const;
404 std::string
name()
const;
411 typename Traits::template Codim<codim>::LevelIterator
lbegin (
int level)
const;
414 typename Traits::template Codim<codim>::LevelIterator
lend (
int level)
const;
417 template<
int codim, PartitionIteratorType PiType>
418 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lbegin (
int level)
const;
420 template<
int codim, PartitionIteratorType PiType>
421 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lend (
int level)
const;
425 typename Traits::template Codim<codim>::LeafIterator
leafbegin()
const;
428 typename Traits::template Codim<codim>::LeafIterator
leafend()
const;
431 template<
int codim, PartitionIteratorType PiType>
432 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafbegin()
const;
434 template<
int codim, PartitionIteratorType PiType>
435 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafend()
const;
438 int size (
int level,
int codim)
const;
441 int size (
int codim)
const;
444 int size (
int level, GeometryType type)
const;
447 int size (GeometryType type)
const;
464 const std::vector<Dune::GeometryType>& geomTypes(
const int)
const;
483 void addLgrUpdateLeafView(
const std::array<int,3>& cells_per_dim,
const std::array<int,3>& startIJK,
484 const std::array<int,3>& endIJK,
const std::string& lgr_name);
503 const std::vector<std::array<int,3>>& startIJK_vec,
504 const std::vector<std::array<int,3>>& endIJK_vec,
505 const std::vector<std::string>& lgr_name_vec);
508 const std::map<std::string,int>& getLgrNameToLevel()
const;
516 std::array<double,3> getEclCentroid(
const int& idx)
const;
540 bool mark(int refCount, const typename Traits::template Codim<0>::EntityPointer & e)
542 return hostgrid_->mark(refCount, getHostEntity<0>(*e));
549 int getMark(const typename Traits::template Codim<0>::EntityPointer & e) const
551 return hostgrid_->getMark(getHostEntity<0>(*e));
556 return hostgrid_->preAdapt();
563 return hostgrid_->adapt();
568 return hostgrid_->postAdapt();
571 end of refinement section */
589 void setZoltanParams(
const std::map<std::string,std::string>& params);
602 return get<0>(scatterGrid(
defaultTransEdgeWgt,
false,
nullptr,
false,
nullptr,
true, overlapLayers, useZoltan ));
626 std::pair<bool,std::vector<std::pair<std::string,bool>>>
628 const double* transmissibilities =
nullptr,
629 int overlapLayers=1,
bool useZoltan=
true)
631 return scatterGrid(
defaultTransEdgeWgt,
false, wells,
false, transmissibilities,
false, overlapLayers, useZoltan);
659 std::pair<bool,std::vector<std::pair<std::string,bool>>>
661 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
662 bool addCornerCells=
false,
int overlapLayers=1,
663 bool useZoltan =
true)
665 return scatterGrid(method, ownersFirst, wells,
false, transmissibilities, addCornerCells, overlapLayers, useZoltan);
687 template<
class DataHandle>
688 std::pair<bool, std::vector<std::pair<std::string,bool> > >
690 const std::vector<cpgrid::OpmWellType> * wells,
691 const double* transmissibilities =
nullptr,
692 int overlapLayers=1,
bool useZoltan =
true)
694 auto ret =
loadBalance(wells, transmissibilities, overlapLayers, useZoltan);
731 template<
class DataHandle>
732 std::pair<bool, std::vector<std::pair<std::string,bool> > >
734 const std::vector<cpgrid::OpmWellType> * wells,
735 bool serialPartitioning,
736 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
737 bool addCornerCells=
false,
int overlapLayers=1,
bool useZoltan =
true,
738 double zoltanImbalanceTol = 1.1,
739 bool allowDistributedWells =
false)
741 auto ret = scatterGrid(method, ownersFirst, wells, serialPartitioning, transmissibilities,
742 addCornerCells, overlapLayers, useZoltan, zoltanImbalanceTol, allowDistributedWells);
767 template<
class DataHandle>
768 std::pair<bool, std::vector<std::pair<std::string,bool> > >
770 const std::vector<cpgrid::OpmWellType> * wells,
771 bool ownersFirst=
false,
772 bool addCornerCells=
false,
int overlapLayers=1)
778 addCornerCells, overlapLayers,
false,
797 template<
class DataHandle>
799 decltype(data.fixedSize(0,0)) overlapLayers=1,
bool useZoltan =
true)
822 bool loadBalance(
const std::vector<int>& parts,
bool ownersFirst=
false,
823 bool addCornerCells=
false,
int overlapLayers=1)
829 addCornerCells, overlapLayers,
false,
846 template<
class DataHandle>
847 bool loadBalance(DataHandle& data,
const std::vector<int>& parts,
bool ownersFirst=
false,
848 bool addCornerCells=
false,
int overlapLayers=1)
850 bool ret =
loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
865 const double* transmissibilities,
867 const double zoltanImbalanceTol)
const;
876 template<
class DataHandle>
877 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int )
const
889 template<
class DataHandle>
890 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir)
const;
906 typedef Dune::FieldVector<double, 3> Vector;
909 const std::vector<double>& zcornData()
const;
935 int cellFace(
int cell,
int local_index)
const;
951 int faceCell(
int face,
int local_index)
const;
961 int numFaceVertices(
int face)
const;
967 int faceVertex(
int face,
int local_index)
const;
974 const Vector faceCenterEcl(
int cell_index,
int face)
const;
976 const Vector faceAreaNormalEcl(
int face)
const;
1010 :
public RandomAccessIteratorFacade<CentroidIterator<codim>,
1011 FieldVector<double, 3>,
1012 const FieldVector<double, 3>&, int>
1016 typedef typename std::vector<
cpgrid::Geometry<3-codim, 3> >::const_iterator
1024 const FieldVector<double,3>& dereference()
const
1026 return iter_->center();
1032 const FieldVector<double,3>& elementAt(
int n)
1034 return iter_[n]->center();
1036 void advance(
int n){
1043 int distanceTo(
const CentroidIterator& o)
1047 bool equals(
const CentroidIterator& o)
const{
1062 int boundaryId(
int face)
const;
1070 template<
class Cell2FacesRowIterator>
1072 faceTag(
const Cell2FacesRowIterator& cell_face)
const;
1094 template<
class DataHandle>
1103 template<
class DataHandle>
1152 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1154 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1157 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1162 const CommunicationType& cellCommunication()
const;
1164 ParallelIndexSet& getCellIndexSet();
1166 RemoteIndices& getCellRemoteIndices();
1168 const ParallelIndexSet& getCellIndexSet()
const;
1170 const RemoteIndices& getCellRemoteIndices()
const;
1202 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1205 const std::vector<cpgrid::OpmWellType> * wells,
1206 bool serialPartitioning,
1207 const double* transmissibilities,
1208 bool addCornerCells,
1210 bool useZoltan =
true,
1211 double zoltanImbalanceTol = 1.1,
1212 bool allowDistributedWells =
true,
1213 const std::vector<int>& input_cell_part = {});
1219 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1223 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1225 std::map<std::string,int> lgr_names_ = {{
"GLOBAL", 0}};
1231 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1237 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1241 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1247 std::map<std::string,std::string> zoltanParams;
1253#include <opm/grid/cpgrid/Entity.hpp>
1254#include <opm/grid/cpgrid/Iterators.hpp>
1255#include <opm/grid/cpgrid/CpGridData.hpp>
1261 namespace Capabilities
1265 struct hasEntity<CpGrid, 0>
1267 static const bool v =
true;
1272 struct hasEntity<CpGrid, 3>
1274 static const bool v =
true;
1278 struct canCommunicate<CpGrid,0>
1280 static const bool v =
true;
1284 struct canCommunicate<CpGrid,3>
1286 static const bool v =
true;
1291 struct hasBackupRestoreFacilities<CpGrid>
1293 static const bool v =
false;
1298 template<
class DataHandle>
1301 current_view_data_->
communicate(data, iftype, dir);
1305 template<
class DataHandle>
1309 if(distributed_data_.empty())
1310 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balanced grid!");
1319 template<
class DataHandle>
1323 if(distributed_data_.empty())
1324 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balance grid!");
1325 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1333 template<
class Cell2FacesRowIterator>
1346 const int cell = cell_face.getCellIndex();
1347 const int face = *cell_face;
1348 assert (0 <= cell); assert (cell <
numCells());
1349 assert (0 <= face); assert (face <
numFaces());
1354 const F2C& f2c = current_view_data_->face_to_cell_[f];
1355 const face_tag tag = current_view_data_->face_tag_[f];
1357 assert ((f2c.size() == 1) || (f2c.size() == 2));
1359 int inside_cell = 0;
1361 if ( f2c.size() == 2 )
1363 if ( f2c[1].index() == cell )
1368 const bool normal_is_in = ! f2c[inside_cell].orientation();
1373 return normal_is_in ? 0 : 1;
1376 return normal_is_in ? 2 : 3;
1380 return normal_is_in ? 4 : 5;
1385 OPM_THROW(std::logic_error,
"Unhandled face tag. This should never happen!");
1394#include <opm/grid/cpgrid/PersistentContainer.hpp>
1395#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1396#include "cpgrid/Intersection.hpp"
1397#include "cpgrid/Geometry.hpp"
1398#include "cpgrid/Indexsets.hpp"
An iterator over the centroids of the geometry of the entities.
Definition CpGrid.hpp:1013
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition CpGrid.hpp:1020
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition CpGrid.hpp:1017
[ provides Dune::Grid ]
Definition CpGrid.hpp:226
std::string name() const
Get the grid name.
Definition CpGrid.cpp:633
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:660
void switchToGlobalView()
Switch to the global view.
Definition CpGrid.cpp:1308
void readSintefLegacyFormat(const std::string &grid_prefix)
Read the Sintef legacy grid format ('topogeom').
Definition CpGrid.cpp:1350
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition CpGrid.cpp:980
void addLgrUpdateLeafView(const std::array< int, 3 > &cells_per_dim, const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK, const std::string &lgr_name)
Create a grid out of a coarse one and a refinement(LGR) of a selected block-shaped patch of cells fro...
Definition CpGrid.cpp:1422
CpGridFamily GridFamily
Family typedef, why is this not defined by Grid<>?
Definition CpGrid.hpp:264
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:1320
void globalRefine(int)
global refinement
Definition CpGrid.cpp:937
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition CpGrid.cpp:913
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
Definition CpGrid.cpp:1242
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition CpGrid.cpp:908
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition CpGrid.hpp:1335
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGrid.cpp:628
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:534
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:1395
unsigned int overlapSize(int) const
Size of the overlap on the leaf level.
Definition CpGrid.cpp:958
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition CpGrid.hpp:1107
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:1227
int numCellFaces() const
Get the sum of all faces attached to all cells.
Definition CpGrid.cpp:1077
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition CpGrid.cpp:604
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:769
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Definition CpGrid.cpp:1237
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition CpGrid.cpp:638
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition CpGrid.cpp:1005
double faceArea(int face) const
Get the area of a face.
Definition CpGrid.cpp:1217
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
Definition CpGrid.cpp:1222
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition CpGrid.cpp:918
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:884
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the global grid.
Definition CpGrid.cpp:597
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
Definition CpGrid.cpp:1087
bool loadBalance(int overlapLayers=1, bool useZoltan=true)
Distributes this grid over the available nodes in a distributed machine.
Definition CpGrid.hpp:599
int numVertices() const
Get The number of vertices.
Definition CpGrid.cpp:1026
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:618
int cellFace(int cell, int local_index) const
Get a specific face of a cell.
Definition CpGrid.cpp:1036
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGrid.cpp:623
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:627
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...
Definition CpGrid.cpp:1429
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:733
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:187
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:798
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1046
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:1360
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition CpGrid.cpp:964
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
Definition CpGrid.cpp:948
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:689
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition CpGrid.cpp:932
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:1303
CpGrid()
Default constructor.
Definition CpGrid.cpp:160
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:1252
void switchToDistributedView()
Switch to the distributed view.
Definition CpGrid.cpp:1313
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
Definition CpGrid.cpp:1247
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:847
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:1041
int numCells() const
Get the number of cells.
Definition CpGrid.cpp:1016
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition CpGrid.cpp:1298
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:822
int numFaces() const
Get the number of faces.
Definition CpGrid.cpp:1021
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
Definition CpGrid.cpp:1212
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
Definition CpGrid.cpp:1092
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition CpGrid.hpp:877
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:1232
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition CpGrid.hpp:1306
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:135
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition CpGridData.hpp:832
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:75
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:278
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
LookUpCartesianData - To search data via CartesianIndex (cartesianMapper)
Definition LookUpData.hh:155
LookUpData class - To search data via element index.
Definition LookUpData.hh:64
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.cpp:68
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:213
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:167
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition CpGrid.hpp:169
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition CpGrid.hpp:171
Traits associated with a specific codim.
Definition CpGrid.hpp:143
cpgrid::Entity< cd > Entity
The type of the entity.
Definition CpGrid.hpp:152
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity.
Definition CpGrid.hpp:146
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition CpGrid.hpp:149
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition CpGrid.hpp:158
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition CpGrid.hpp:155
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition CpGrid.hpp:161
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:179
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:183
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:181
Definition CpGrid.hpp:123
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition CpGrid.hpp:193
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition CpGrid.hpp:132
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:190
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:188
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition CpGrid.hpp:197
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition CpGrid.hpp:134
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition CpGrid.hpp:195
GlobalIdSet LocalIdSet
The type of the local id set.
Definition CpGrid.hpp:199
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition CpGrid.hpp:202
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition CpGrid.hpp:137
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition CpGrid.hpp:130
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition CpGrid.hpp:128
CpGrid Grid
The type that implements the grid.
Definition CpGrid.hpp:125
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56