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>
47#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
49#include <opm/grid/common/GridEnums.hpp>
51#include <opm/grid/cpgrid/CpGridDataTraits.hpp>
52#include <opm/grid/cpgrid/DefaultGeometryPolicy.hpp>
53#include <opm/grid/cpgrid/OrientedEntityTable.hpp>
57#include <opm/grid/utility/OpmWellType.hpp>
75 template <
int>
class Entity;
79 template<
int, PartitionIteratorType>
class Iterator;
142 template <PartitionIteratorType pitype>
154 template <PartitionIteratorType pitype>
160 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid> >
LeafGridView;
167 typedef Dune::GridView<DefaultLeafGridViewTraits<CpGrid>>
LeafGridView;
180 using CollectiveCommunication = cpgrid::CpGridDataTraits::CollectiveCommunication;
202 :
public GridDefaultImplementation<3, 3, double, CpGridFamily>
227 explicit CpGrid(MPIHelper::MPICommunicator
comm);
266 std::vector<std::size_t>
268 Opm::EclipseState* ecl_state,
269 bool periodic_extension,
273 bool edge_conformal);
312 std::vector<std::size_t>
314 Opm::EclipseState* ecl_state,
315 bool periodic_extension,
316 bool turn_normals =
false,
318 bool edge_conformal =
false);
336 bool remove_ij_boundary,
337 bool turn_normals =
false,
338 bool edge_conformal =
false);
356 const std::array<double, 3>& cellsize,
357 const std::array<int, 3>& shift = {0,0,0});
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>&
currentData()
const;
377 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>&
currentData();
386 void getIJK(
const int c, std::array<int,3>& ijk)
const;
407 std::string
name()
const;
414 typename Traits::template Codim<codim>::LevelIterator
lbegin (
int level)
const;
417 typename Traits::template Codim<codim>::LevelIterator
lend (
int level)
const;
420 template<
int codim, PartitionIteratorType PiType>
421 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lbegin (
int level)
const;
423 template<
int codim, PartitionIteratorType PiType>
424 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lend (
int level)
const;
428 typename Traits::template Codim<codim>::LeafIterator
leafbegin()
const;
431 typename Traits::template Codim<codim>::LeafIterator
leafend()
const;
434 template<
int codim, PartitionIteratorType PiType>
435 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafbegin()
const;
437 template<
int codim, PartitionIteratorType PiType>
438 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafend()
const;
441 int size (
int level,
int codim)
const;
444 int size (
int codim)
const;
447 int size (
int level, GeometryType type)
const;
450 int size (GeometryType type)
const;
469 const std::vector<Dune::GeometryType>& geomTypes(
const int)
const;
494 const std::vector<std::array<int,3>>& startIJK_vec,
495 const std::vector<std::array<int,3>>& endIJK_vec,
496 const std::vector<std::string>& lgr_name_vec,
497 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
505 void autoRefine(
const std::array<int,3>& nxnynz);
508 const std::map<std::string,int>& getLgrNameToLevel()
const;
516 std::array<double,3> getEclCentroid(
const int& idx)
const;
580 const std::vector<int>& assignRefinedLevel,
581 const std::vector<std::string>& lgr_name_vec,
582 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
583 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
590 void updateCornerHistoryLevels(
const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
591 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
592 const std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
593 const int& corner_count,
594 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
595 const int& preAdaptMaxLevel,
596 const int& newLevels);
598 void globalIdsPartitionTypesLgrAndLeafGrids(
const std::vector<int>& assignRefinedLevel,
599 const std::vector<std::array<int,3>>& cells_per_dim_vec,
600 const std::vector<int>& lgr_with_at_least_one_active_cell);
608 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
630 void computeGlobalCellLgr(
const int& level,
const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
635 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
651 void markElemAssignLevelDetectActiveLgrs(
const std::vector<std::array<int,3>>& startIJK_vec,
652 const std::vector<std::array<int,3>>& endIJK_vec,
653 std::vector<int>& assignRefinedLevel,
654 std::vector<int>& lgr_with_at_least_one_active_cell);
657 void populateCellIndexSetRefinedGrid(
int level);
660 void populateCellIndexSetLeafGridView();
663 void populateLeafGlobalIdSet();
692 void setPartitioningParams(
const std::map<std::string,std::string>& params);
707 double imbalanceTol = 1.1,
741 double imbalanceTol = 1.1,
790 std::pair<bool,std::vector<std::pair<std::string,bool>>>
792 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
793 const double* transmissibilities =
nullptr,
798 false, transmissibilities,
false,
799 overlapLayers, partitionMethod, 1.1,
false,
838 std::pair<bool,std::vector<std::pair<std::string,bool>>>
840 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
841 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
842 bool addCornerCells=
false,
int overlapLayers=1,
844 double imbalanceTol = 1.1,
847 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections,
false,
848 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
880 template<
class DataHandle>
881 std::pair<bool, std::vector<std::pair<std::string,bool> > >
883 const std::vector<cpgrid::OpmWellType> * wells,
884 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
885 const double* transmissibilities =
nullptr,
886 int overlapLayers=1,
int partitionMethod = 1,
int level =-1)
888 auto ret =
loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
931 template<
class DataHandle>
932 std::pair<bool, std::vector<std::pair<std::string,bool> > >
934 const std::vector<cpgrid::OpmWellType> * wells,
935 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
936 bool serialPartitioning,
937 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
939 double imbalanceTol = 1.1,
940 bool allowDistributedWells =
false)
942 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
943 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
944 std::vector<int>{}, 0);
977 template<
class DataHandle>
978 std::pair<bool, std::vector<std::pair<std::string,bool> > >
980 const std::vector<cpgrid::OpmWellType> * wells,
981 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
982 bool ownersFirst=
false,
983 bool addCornerCells=
false,
int overlapLayers=1)
987 possibleFutureConnections,
1008 template<
class DataHandle>
1014 bool ret =
loadBalance(overlapLayers, partitionMethod);
1033 bool loadBalance(
const std::vector<int>& parts,
bool ownersFirst=
false,
1034 bool addCornerCells=
false,
int overlapLayers=1)
1058 template<
class DataHandle>
1059 bool loadBalance(DataHandle& data,
const std::vector<int>& parts,
bool ownersFirst=
false,
1060 bool addCornerCells=
false,
int overlapLayers=1)
1062 bool ret =
loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1085 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1086 const double* transmissibilities,
1088 const double imbalanceTol)
const;
1097 template<
class DataHandle>
1098 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int )
const
1110 template<
class DataHandle>
1111 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir)
const;
1127 typedef Dune::FieldVector<double, 3> Vector;
1130 const std::vector<double>& zcornData()
const;
1138 int numCells(
int level = -1)
const;
1144 int numFaces(
int level = -1)
const;
1166 int cellFace(
int cell,
int local_index,
int level = -1)
const;
1170 const cpgrid::OrientedEntityTable<0,1>::row_type
cellFaceRow(
int cell)
const;
1191 int faceCell(
int face,
int local_index,
int level = -1)
const;
1201 int numFaceVertices(
int face)
const;
1207 int faceVertex(
int face,
int local_index)
const;
1216 const Vector faceAreaNormalEcl(
int face)
const;
1250 :
public RandomAccessIteratorFacade<CentroidIterator<codim>,
1251 FieldVector<double, 3>,
1252 const FieldVector<double, 3>&, int>
1256 typedef typename std::vector<
cpgrid::Geometry<3-codim, 3> >::const_iterator
1264 const FieldVector<double,3>& dereference()
const
1266 return iter_->center();
1272 const FieldVector<double,3>& elementAt(
int n)
1274 return iter_[n]->center();
1276 void advance(
int n){
1302 int boundaryId(
int face)
const;
1310 template<
class Cell2FacesRowIterator>
1312 faceTag(
const Cell2FacesRowIterator& cell_face)
const;
1334 template<
class DataHandle>
1343 template<
class DataHandle>
1392 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1394 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1397 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1402 const CommunicationType& cellCommunication()
const;
1404 ParallelIndexSet& getCellIndexSet();
1406 RemoteIndices& getCellRemoteIndices();
1408 const ParallelIndexSet& getCellIndexSet()
const;
1410 const RemoteIndices& getCellRemoteIndices()
const;
1448 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1451 const std::vector<cpgrid::OpmWellType> * wells,
1452 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1453 bool serialPartitioning,
1454 const double* transmissibilities,
1455 bool addCornerCells,
1458 double imbalanceTol = 1.1,
1459 bool allowDistributedWells =
true,
1460 const std::vector<int>& input_cell_part = {},
1467 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1471 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1473 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1475 std::map<std::string,int> lgr_names_ = {{
"GLOBAL", 0}};
1481 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1487 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1491 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1497 std::map<std::string,std::string> partitioningParams;
1503#include <opm/grid/cpgrid/Entity.hpp>
1504#include <opm/grid/cpgrid/Iterators.hpp>
1505#include <opm/grid/cpgrid/CpGridData.hpp>
1511 namespace Capabilities
1517 static const bool v =
true;
1524 static const bool v =
true;
1530 static const bool v =
true;
1536 static const bool v =
true;
1543 static const bool v =
false;
1548 template<
class DataHandle>
1551 current_view_data_->
communicate(data, iftype, dir);
1555 template<
class DataHandle>
1559 if (distributed_data_.empty()) {
1560 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balanced grid!");
1562 distributed_data_[0]->scatterData(handle, data_[0].get(),
1563 distributed_data_[0].get(),
1570 template<
class DataHandle>
1574 if (distributed_data_.empty()) {
1575 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balance grid!");
1577 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1583 template<
class Cell2FacesRowIterator>
1596 const int cell = cell_face.getCellIndex();
1597 const int face = *cell_face;
1598 assert (0 <= cell); assert (cell <
numCells());
1599 assert (0 <= face); assert (face <
numFaces());
1601 typedef cpgrid::OrientedEntityTable<1,0>::row_type F2C;
1604 const F2C& f2c = current_view_data_->face_to_cell_[f];
1605 const face_tag tag = current_view_data_->face_tag_[f];
1607 assert ((f2c.size() == 1) || (f2c.size() == 2));
1609 int inside_cell = 0;
1611 if ( f2c.size() == 2 )
1613 if ( f2c[1].index() == cell )
1618 const bool normal_is_in = ! f2c[inside_cell].orientation();
1623 return normal_is_in ? 0 : 1;
1626 return normal_is_in ? 2 : 3;
1630 return normal_is_in ? 4 : 5;
1635 OPM_THROW(std::logic_error,
"Unhandled face tag. This should never happen!");
1644#include <opm/grid/cpgrid/PersistentContainer.hpp>
1645#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1646#include "cpgrid/Intersection.hpp"
1647#include "cpgrid/Geometry.hpp"
1648#include "cpgrid/Indexsets.hpp"
An iterator over the centroids of the geometry of the entities.
Definition CpGrid.hpp:1253
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition CpGrid.hpp:1260
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition CpGrid.hpp:1257
[ provides Dune::Grid ]
Definition CpGrid.hpp:203
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.
Definition CpGrid.hpp:979
std::string name() const
Get the grid name.
Definition CpGrid.cpp:761
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 inde...
Definition CpGrid.cpp:726
void switchToGlobalView()
Switch to the global view.
Definition CpGrid.cpp:1660
unsigned int numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition CpGrid.cpp:1111
CpGridFamily GridFamily
Family typedef, why is this not defined by Grid<>?
Definition CpGrid.hpp:218
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:1571
int numFaces(int level=-1) const
Get the number of faces.
Definition CpGrid.cpp:1153
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition CpGrid.cpp:1018
std::vector< std::array< int, 2 > > mapLeafIndexSetToLocalCartesianIndexSets() const
Reverse map: from leaf index cell to { level, local/level Cartesian index of the cell }...
Definition CpGrid.cpp:736
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData() const
Returns either data_ or distributed_data_(if non empty).
Definition CpGrid.cpp:670
CentroidIterator< 0 > beginCellCentroids() const
Get an iterator over the cell centroids positioned at the first one.
Definition CpGrid.cpp:1594
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition CpGrid.cpp:1013
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.
Definition CpGrid.hpp:933
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition CpGrid.hpp:1585
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGrid.cpp:756
bool refineAndUpdateGrid(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.
Definition CpGrid.cpp:1846
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.
Definition CpGrid.hpp:839
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:591
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
unsigned int overlapSize(int) const
Size of the overlap on the leaf level.
Definition CpGrid.cpp:1089
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition CpGrid.hpp:1347
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.
Definition CpGrid.hpp:882
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:1579
int numCellFaces() const
Get the sum of all faces attached to all cells.
Definition CpGrid.cpp:1221
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition CpGrid.cpp:680
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.
Definition CpGrid.hpp:738
const Vector & cellCentroid(int cell) const
Get the coordinates of the center of a cell.
Definition CpGrid.cpp:1589
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.
Definition CpGrid.hpp:705
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition CpGrid.cpp:766
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition CpGrid.cpp:1136
double faceArea(int face) const
Get the area of a face.
Definition CpGrid.cpp:1569
const Vector & faceCentroid(int face) const
Get the coordinates of the center of a face.
Definition CpGrid.cpp:1574
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition CpGrid.cpp:1023
void postAdapt()
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
Definition CpGrid.cpp:2385
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:989
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the global grid.
Definition CpGrid.cpp:663
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
Definition CpGrid.cpp:1231
int numVertices() const
Get The number of vertices.
Definition CpGrid.cpp:1159
Dune::cpgrid::Intersection getParentIntersectionFromLgrBoundaryFace(const Dune::cpgrid::Intersection &intersection) const
Definition CpGrid.cpp:1236
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:746
void processEclipseFormat(const grdecl &input_data, bool remove_ij_boundary, bool turn_normals=false, bool edge_conformal=false)
Read the Eclipse grid format ('grdecl').
Definition CpGrid.cpp:1743
void syncDistributedGlobalCellIds()
Synchronizes cell global ids across processes after load balancing.
Definition CpGrid.cpp:2530
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGrid.cpp:751
int faceCell(int face, int local_index, int level=-1) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1183
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
------------— Adaptivity (begin) ------------—
Definition CpGrid.cpp:1776
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, const std::vector< std::string > &lgr_parent_grid_name_vec=std::vector< std::string >{})
Create a grid out of a coarse one and (at most) 2 refinements(LGRs) of selected block-shaped disjoint...
Definition CpGrid.cpp:2593
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.
Definition CpGrid.hpp:1009
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition CpGrid.cpp:1095
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.
Definition CpGrid.hpp:791
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
Definition CpGrid.cpp:1079
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition CpGrid.cpp:1030
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:1655
CpGrid()
Default constructor.
Definition CpGrid.cpp:167
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:1604
bool adapt()
Triggers the grid refinement process.
Definition CpGrid.cpp:1812
void switchToDistributedView()
Switch to the distributed view.
Definition CpGrid.cpp:1666
CentroidIterator< 1 > beginFaceCentroids() const
Get an iterator over the face centroids positioned at the first one.
Definition CpGrid.cpp:1599
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
Definition CpGrid.cpp:1800
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
Definition CpGrid.cpp:1795
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:1059
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:1178
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition CpGrid.cpp:1650
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:1033
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.
Definition CpGrid.cpp:194
int numCells(int level=-1) const
Get the number of cells.
Definition CpGrid.cpp:1147
const Vector & vertexPosition(int vertex) const
Get the Position of a vertex.
Definition CpGrid.cpp:1564
double cellCenterDepth(int cell_index) const
Get vertical position of cell center ("zcorn" average).
Definition CpGrid.cpp:1422
void globalRefine(int refCount)
Refine the grid refCount times using the default refinement rule.
Definition CpGrid.cpp:1035
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition CpGrid.hpp:1098
int cellFace(int cell, int local_index, int level=-1) const
Get a specific face of a cell.
Definition CpGrid.cpp:1171
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:1584
void autoRefine(const std::array< int, 3 > &nxnynz)
Global refine the grid with different refinement factors in each direction.
Definition CpGrid.cpp:2750
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition CpGrid.hpp:1556
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:118
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:76
The global id set for Dune.
Definition Indexsets.hpp:487
Only needs to provide interface for doing nothing.
Definition Iterators.hpp:118
Definition Indexsets.hpp:201
Definition Indexsets.hpp:56
Definition Intersection.hpp:279
Definition Intersection.hpp:66
Iterator intended to be used as LeafIterator and LevelIterator (no difference due to no adaptivity) f...
Definition Iterators.hpp:60
Definition Indexsets.hpp:371
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
@ simple
Use simple approach based on rectangular partitioning the underlying cartesian grid.
Definition GridEnums.hpp:46
@ zoltanGoG
use Zoltan on GraphOfGrid for partitioning
Definition GridEnums.hpp:52
@ zoltan
Use Zoltan for partitioning.
Definition GridEnums.hpp:48
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
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:190
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:144
cpgrid::Iterator< cd, pitype > LevelIterator
The type of the iterator over the level entities of this codim on this partition.
Definition CpGrid.hpp:146
cpgrid::Iterator< cd, pitype > LeafIterator
The type of the iterator over the leaf entities of this codim on this partition.
Definition CpGrid.hpp:148
Traits associated with a specific codim.
Definition CpGrid.hpp:120
cpgrid::Entity< cd > Entity
The type of the entity.
Definition CpGrid.hpp:129
cpgrid::Geometry< 3-cd, 3 > Geometry
The type of the geometry associated with the entity.
Definition CpGrid.hpp:123
cpgrid::Geometry< 3-cd, 3 > LocalGeometry
The type of the local geometry associated with the entity.
Definition CpGrid.hpp:126
cpgrid::Iterator< cd, All_Partition > LeafIterator
The type of the iterator over all leaf entities of this codim.
Definition CpGrid.hpp:135
cpgrid::Iterator< cd, All_Partition > LevelIterator
The type of the iterator over all level entities of this codim.
Definition CpGrid.hpp:132
cpgrid::Entity< cd > EntitySeed
The type of the entity pointer for entities of this codim.
Definition CpGrid.hpp:138
Traits associated with a specific grid partition type.
Definition CpGrid.hpp:156
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:160
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:158
Definition CpGrid.hpp:100
cpgrid::IndexSet LevelIndexSet
The type of the level index set.
Definition CpGrid.hpp:170
cpgrid::IntersectionIterator LeafIntersectionIterator
The type of the intersection iterator at the leafs of the grid.
Definition CpGrid.hpp:109
Dune::GridView< DefaultLeafGridViewTraits< CpGrid > > LeafGridView
The type of the leaf grid view associated with this partition type.
Definition CpGrid.hpp:167
Dune::GridView< DefaultLevelGridViewTraits< CpGrid > > LevelGridView
The type of the level grid view associated with this partition type.
Definition CpGrid.hpp:165
cpgrid::GlobalIdSet GlobalIdSet
The type of the global id set.
Definition CpGrid.hpp:174
cpgrid::IntersectionIterator LevelIntersectionIterator
The type of the intersection iterator at the levels of the grid.
Definition CpGrid.hpp:111
cpgrid::IndexSet LeafIndexSet
The type of the leaf index set.
Definition CpGrid.hpp:172
GlobalIdSet LocalIdSet
The type of the local id set.
Definition CpGrid.hpp:176
cpgrid::CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition CpGrid.hpp:179
cpgrid::HierarchicIterator HierarchicIterator
The type of the hierarchic iterator.
Definition CpGrid.hpp:114
cpgrid::Intersection LevelIntersection
The type of the intersection at the levels of the grid.
Definition CpGrid.hpp:107
cpgrid::Intersection LeafIntersection
The type of the intersection at the leafs of the grid.
Definition CpGrid.hpp:105
CpGrid Grid
The type that implements the grid.
Definition CpGrid.hpp:102
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56