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;
74 template <
int>
class Entity;
75 template<
int,
int>
class Geometry;
76 class HierarchicIterator;
77 class IntersectionIterator;
78 template<
int, PartitionIteratorType>
class Iterator;
79 class LevelGlobalIdSet;
82 class IntersectionIterator;
90 const std::vector<std::array<int,3>>&,
91 const std::vector<std::array<int,3>>&,
92 const std::vector<std::array<int,3>>&,
93 const std::vector<std::string>&);
95void testCase(
const std::string&,
97 const std::vector<std::array<int,3>>&,
98 const std::vector<std::array<int,3>>&,
99 const std::vector<std::array<int,3>>&,
100 const std::vector<std::string>&,
104 const std::vector<std::array<int,3>>&,
105 const std::vector<std::array<int,3>>&);
110 const std::array<int, 3>&,
114 const std::vector<std::array<int,3>>&,
115 const std::vector<std::array<int,3>>&,
116 const std::vector<std::array<int,3>>&,
117 const std::vector<std::string>&);
119void refinePatch_and_check(
const std::array<int,3>&,
120 const std::array<int,3>&,
121 const std::array<int,3>&);
177 template <PartitionIteratorType pitype>
189 template <PartitionIteratorType pitype>
195 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> >
LeafGridView;
202 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>>
LeafGridView;
215 using CollectiveCommunication = cpgrid::CpGridDataTraits::CollectiveCommunication;
237 :
public GridDefaultImplementation<3, 3, double, CpGridFamily>
244 template<
typename Gr
id,
typename Gr
idView>
friend class Opm::LookUpData;
249 const std::vector<std::array<int,3>>&,
250 const std::vector<std::array<int,3>>&,
251 const std::vector<std::array<int,3>>&,
252 const std::vector<std::string>&);
253 friend void ::testCase(
const std::string&,
255 const std::vector<std::array<int,3>>&,
256 const std::vector<std::array<int,3>>&,
257 const std::vector<std::array<int,3>>&,
258 const std::vector<std::string>&,
261 const std::vector<std::array<int,3>>&,
262 const std::vector<std::array<int,3>>&);
266 const std::array<int,3>&,
270 const std::vector<std::array<int,3>>&,
271 const std::vector<std::array<int,3>>&,
272 const std::vector<std::array<int,3>>&,
273 const std::vector<std::string>&);
275 void ::refinePatch_and_check(
const std::array<int,3>&,
276 const std::array<int,3>&,
277 const std::array<int,3>&);
333 Opm::EclipseState* ecl_state,
334 bool periodic_extension,
bool turn_normals,
bool clip_z,
357 Opm::EclipseState* ecl_state,
358 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false);
383 const std::array<double, 3>& cellsize,
384 const std::array<int, 3>& shift = {0,0,0});
401 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>&
chooseData()
const;
410 void getIJK(
const int c, std::array<int,3>& ijk)
const;
431 std::string
name()
const;
438 typename Traits::template Codim<codim>::LevelIterator
lbegin (
int level)
const;
441 typename Traits::template Codim<codim>::LevelIterator
lend (
int level)
const;
444 template<
int codim, PartitionIteratorType PiType>
445 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lbegin (
int level)
const;
447 template<
int codim, PartitionIteratorType PiType>
448 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lend (
int level)
const;
452 typename Traits::template Codim<codim>::LeafIterator
leafbegin()
const;
455 typename Traits::template Codim<codim>::LeafIterator
leafend()
const;
458 template<
int codim, PartitionIteratorType PiType>
459 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafbegin()
const;
461 template<
int codim, PartitionIteratorType PiType>
462 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafend()
const;
465 int size (
int level,
int codim)
const;
468 int size (
int codim)
const;
471 int size (
int level, GeometryType type)
const;
474 int size (GeometryType type)
const;
491 const std::vector<Dune::GeometryType>& geomTypes(
const int)
const;
510 void addLgrUpdateLeafView(
const std::array<int,3>& cells_per_dim,
const std::array<int,3>& startIJK,
511 const std::array<int,3>& endIJK,
const std::string& lgr_name);
530 const std::vector<std::array<int,3>>& startIJK_vec,
531 const std::vector<std::array<int,3>>& endIJK_vec,
532 const std::vector<std::string>& lgr_name_vec);
535 const std::map<std::string,int>& getLgrNameToLevel()
const;
543 std::array<double,3> getEclCentroid(
const int& idx)
const;
570 bool mark(int refCount, const typename Traits::template Codim<0>::EntityPointer & e)
572 return hostgrid_->mark(refCount, getHostEntity<0>(*e));
579 int getMark(const typename Traits::template Codim<0>::EntityPointer & e) const
581 return hostgrid_->getMark(getHostEntity<0>(*e));
586 return hostgrid_->preAdapt();
593 return hostgrid_->adapt();
598 return hostgrid_->postAdapt();
601 end of refinement section */
619 void setZoltanParams(
const std::map<std::string,std::string>& params);
632 return get<0>(scatterGrid(
defaultTransEdgeWgt,
false,
nullptr,
false,
nullptr,
true, overlapLayers, useZoltan ));
656 std::pair<bool,std::vector<std::pair<std::string,bool>>>
658 const double* transmissibilities =
nullptr,
659 int overlapLayers=1,
bool useZoltan=
true)
661 return scatterGrid(
defaultTransEdgeWgt,
false, wells,
false, transmissibilities,
false, overlapLayers, useZoltan);
689 std::pair<bool,std::vector<std::pair<std::string,bool>>>
691 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
692 bool addCornerCells=
false,
int overlapLayers=1,
693 bool useZoltan =
true)
695 return scatterGrid(method, ownersFirst, wells,
false, transmissibilities, addCornerCells, overlapLayers, useZoltan);
717 template<
class DataHandle>
718 std::pair<bool, std::vector<std::pair<std::string,bool> > >
720 const std::vector<cpgrid::OpmWellType> * wells,
721 const double* transmissibilities =
nullptr,
722 int overlapLayers=1,
bool useZoltan =
true)
724 auto ret =
loadBalance(wells, transmissibilities, overlapLayers, useZoltan);
761 template<
class DataHandle>
762 std::pair<bool, std::vector<std::pair<std::string,bool> > >
764 const std::vector<cpgrid::OpmWellType> * wells,
765 bool serialPartitioning,
766 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
767 bool addCornerCells=
false,
int overlapLayers=1,
bool useZoltan =
true,
768 double zoltanImbalanceTol = 1.1,
769 bool allowDistributedWells =
false)
771 auto ret = scatterGrid(method, ownersFirst, wells, serialPartitioning, transmissibilities,
772 addCornerCells, overlapLayers, useZoltan, zoltanImbalanceTol, allowDistributedWells);
797 template<
class DataHandle>
798 std::pair<bool, std::vector<std::pair<std::string,bool> > >
800 const std::vector<cpgrid::OpmWellType> * wells,
801 bool ownersFirst=
false,
802 bool addCornerCells=
false,
int overlapLayers=1)
808 addCornerCells, overlapLayers,
false,
827 template<
class DataHandle>
829 decltype(data.fixedSize(0,0)) overlapLayers=1,
bool useZoltan =
true)
852 bool loadBalance(
const std::vector<int>& parts,
bool ownersFirst=
false,
853 bool addCornerCells=
false,
int overlapLayers=1)
859 addCornerCells, overlapLayers,
false,
876 template<
class DataHandle>
877 bool loadBalance(DataHandle& data,
const std::vector<int>& parts,
bool ownersFirst=
false,
878 bool addCornerCells=
false,
int overlapLayers=1)
880 bool ret =
loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
895 const double* transmissibilities,
897 const double zoltanImbalanceTol)
const;
906 template<
class DataHandle>
907 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int )
const
919 template<
class DataHandle>
920 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir)
const;
936 typedef Dune::FieldVector<double, 3> Vector;
939 const std::vector<double>& zcornData()
const;
965 int cellFace(
int cell,
int local_index)
const;
981 int faceCell(
int face,
int local_index)
const;
991 int numFaceVertices(
int face)
const;
997 int faceVertex(
int face,
int local_index)
const;
1006 const Vector faceAreaNormalEcl(
int face)
const;
1040 :
public RandomAccessIteratorFacade<CentroidIterator<codim>,
1041 FieldVector<double, 3>,
1042 const FieldVector<double, 3>&, int>
1046 typedef typename std::vector<
cpgrid::Geometry<3-codim, 3> >::const_iterator
1054 const FieldVector<double,3>& dereference()
const
1056 return iter_->center();
1062 const FieldVector<double,3>& elementAt(
int n)
1064 return iter_[n]->center();
1066 void advance(
int n){
1073 int distanceTo(
const CentroidIterator& o)
1077 bool equals(
const CentroidIterator& o)
const{
1092 int boundaryId(
int face)
const;
1100 template<
class Cell2FacesRowIterator>
1102 faceTag(
const Cell2FacesRowIterator& cell_face)
const;
1124 template<
class DataHandle>
1133 template<
class DataHandle>
1182 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1184 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1187 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1192 const CommunicationType& cellCommunication()
const;
1194 ParallelIndexSet& getCellIndexSet();
1196 RemoteIndices& getCellRemoteIndices();
1198 const ParallelIndexSet& getCellIndexSet()
const;
1200 const RemoteIndices& getCellRemoteIndices()
const;
1232 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1235 const std::vector<cpgrid::OpmWellType> * wells,
1236 bool serialPartitioning,
1237 const double* transmissibilities,
1238 bool addCornerCells,
1240 bool useZoltan =
true,
1241 double zoltanImbalanceTol = 1.1,
1242 bool allowDistributedWells =
true,
1243 const std::vector<int>& input_cell_part = {});
1249 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1253 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1255 std::map<std::string,int> lgr_names_ = {{
"GLOBAL", 0}};
1261 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1267 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1271 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1277 std::map<std::string,std::string> zoltanParams;
1283#include <opm/grid/cpgrid/Entity.hpp>
1284#include <opm/grid/cpgrid/Iterators.hpp>
1285#include <opm/grid/cpgrid/CpGridData.hpp>
1291 namespace Capabilities
1295 struct hasEntity<CpGrid, 0>
1297 static const bool v =
true;
1302 struct hasEntity<CpGrid, 3>
1304 static const bool v =
true;
1308 struct canCommunicate<CpGrid,0>
1310 static const bool v =
true;
1314 struct canCommunicate<CpGrid,3>
1316 static const bool v =
true;
1321 struct hasBackupRestoreFacilities<CpGrid>
1323 static const bool v =
false;
1328 template<
class DataHandle>
1331 current_view_data_->
communicate(data, iftype, dir);
1335 template<
class DataHandle>
1339 if(distributed_data_.empty())
1340 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balanced grid!");
1349 template<
class DataHandle>
1353 if(distributed_data_.empty())
1354 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balance grid!");
1355 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1363 template<
class Cell2FacesRowIterator>
1376 const int cell = cell_face.getCellIndex();
1377 const int face = *cell_face;
1378 assert (0 <= cell); assert (cell <
numCells());
1379 assert (0 <= face); assert (face <
numFaces());
1384 const F2C& f2c = current_view_data_->face_to_cell_[f];
1385 const face_tag tag = current_view_data_->face_tag_[f];
1387 assert ((f2c.size() == 1) || (f2c.size() == 2));
1389 int inside_cell = 0;
1391 if ( f2c.size() == 2 )
1393 if ( f2c[1].index() == cell )
1398 const bool normal_is_in = ! f2c[inside_cell].orientation();
1403 return normal_is_in ? 0 : 1;
1406 return normal_is_in ? 2 : 3;
1410 return normal_is_in ? 4 : 5;
1415 OPM_THROW(std::logic_error,
"Unhandled face tag. This should never happen!");
1424#include <opm/grid/cpgrid/PersistentContainer.hpp>
1425#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1426#include "cpgrid/Intersection.hpp"
1427#include "cpgrid/Geometry.hpp"
1428#include "cpgrid/Indexsets.hpp"
An iterator over the centroids of the geometry of the entities.
Definition CpGrid.hpp:1043
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition CpGrid.hpp:1050
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition CpGrid.hpp:1047
[ provides Dune::Grid ]
Definition CpGrid.hpp:238
std::string name() const
Get the grid name.
Definition CpGrid.cpp:636
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:690
void switchToGlobalView()
Switch to the global view.
Definition CpGrid.cpp:1361
void readSintefLegacyFormat(const std::string &grid_prefix)
Read the Sintef legacy grid format ('topogeom').
Definition CpGrid.cpp:1403
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition CpGrid.cpp:976
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:1475
CpGridFamily GridFamily
Family typedef, why is this not defined by Grid<>?
Definition CpGrid.hpp:288
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:1350
void globalRefine(int)
global refinement
Definition CpGrid.cpp:933
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition CpGrid.cpp:916
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
Definition CpGrid.cpp:1295
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition CpGrid.cpp:911
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition CpGrid.hpp:1365
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGrid.cpp:631
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:1448
unsigned int overlapSize(int) const
Size of the overlap on the leaf level.
Definition CpGrid.cpp:954
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition CpGrid.hpp:1137
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:1280
int numCellFaces() const
Get the sum of all faces attached to all cells.
Definition CpGrid.cpp:1073
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition CpGrid.cpp:614
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:799
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Definition CpGrid.cpp:1290
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition CpGrid.cpp:641
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition CpGrid.cpp:1001
double faceArea(int face) const
Get the area of a face.
Definition CpGrid.cpp:1270
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
Definition CpGrid.cpp:1275
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition CpGrid.cpp:921
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:887
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:1083
bool loadBalance(int overlapLayers=1, bool useZoltan=true)
Distributes this grid over the available nodes in a distributed machine.
Definition CpGrid.hpp:629
int numVertices() const
Get The number of vertices.
Definition CpGrid.cpp:1022
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:621
int cellFace(int cell, int local_index) const
Get a specific face of a cell.
Definition CpGrid.cpp:1032
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGrid.cpp:626
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:657
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:1482
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:763
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:828
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1042
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:1413
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition CpGrid.cpp:960
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
Definition CpGrid.cpp:944
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:719
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition CpGrid.cpp:928
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & chooseData() const
Returns either data_ or distributed_data_(if non empty).
Definition CpGrid.cpp:604
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:1356
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:1305
void switchToDistributedView()
Switch to the distributed view.
Definition CpGrid.cpp:1366
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
Definition CpGrid.cpp:1300
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:877
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:1037
int numCells() const
Get the number of cells.
Definition CpGrid.cpp:1012
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition CpGrid.cpp:1351
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:852
int numFaces() const
Get the number of faces.
Definition CpGrid.cpp:1017
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
Definition CpGrid.cpp:1265
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
Definition CpGrid.cpp:1123
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition CpGrid.hpp:907
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:1285
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition CpGrid.hpp:1336
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:131
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition CpGridData.hpp:863
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 field properties of leaf grid view elements via CartesianIndex (carte...
Definition LookUpData.hh:173
LookUpData class - To search field properties of leaf grid view elements via element/elementIndex.
Definition LookUpData.hh:63
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:225
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:179
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition CpGrid.hpp:181
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition CpGrid.hpp:183
Traits associated with a specific codim.
Definition CpGrid.hpp:155
cpgrid::Entity< cd > Entity
The type of the entity.
Definition CpGrid.hpp:164
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity.
Definition CpGrid.hpp:158
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition CpGrid.hpp:161
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition CpGrid.hpp:170
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition CpGrid.hpp:167
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition CpGrid.hpp:173
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:191
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:195
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:193
Definition CpGrid.hpp:135
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition CpGrid.hpp:205
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition CpGrid.hpp:144
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:202
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:200
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition CpGrid.hpp:209
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition CpGrid.hpp:146
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition CpGrid.hpp:207
GlobalIdSet LocalIdSet
The type of the local id set.
Definition CpGrid.hpp:211
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition CpGrid.hpp:214
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition CpGrid.hpp:149
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition CpGrid.hpp:142
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition CpGrid.hpp:140
CpGrid Grid
The type that implements the grid.
Definition CpGrid.hpp:137
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56