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();
392 void getIJK(
const int c, std::array<int,3>& ijk)
const;
413 std::string
name()
const;
420 typename Traits::template Codim<codim>::LevelIterator
lbegin (
int level)
const;
423 typename Traits::template Codim<codim>::LevelIterator
lend (
int level)
const;
426 template<
int codim, PartitionIteratorType PiType>
427 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lbegin (
int level)
const;
429 template<
int codim, PartitionIteratorType PiType>
430 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator
lend (
int level)
const;
434 typename Traits::template Codim<codim>::LeafIterator
leafbegin()
const;
437 typename Traits::template Codim<codim>::LeafIterator
leafend()
const;
440 template<
int codim, PartitionIteratorType PiType>
441 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafbegin()
const;
443 template<
int codim, PartitionIteratorType PiType>
444 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator
leafend()
const;
447 int size (
int level,
int codim)
const;
450 int size (
int codim)
const;
453 int size (
int level, GeometryType type)
const;
456 int size (GeometryType type)
const;
475 void globalRefine (
int refCount,
bool throwOnFailure =
false);
477 const std::vector<Dune::GeometryType>& geomTypes(
const int)
const;
502 const std::vector<std::array<int,3>>& startIJK_vec,
503 const std::vector<std::array<int,3>>& endIJK_vec,
504 const std::vector<std::string>& lgr_name_vec,
505 const std::vector<std::string>& lgr_parent_grid_name_vec = std::vector<std::string>{});
513 void autoRefine(
const std::array<int,3>& nxnynz);
516 const std::map<std::string,int>& getLgrNameToLevel()
const;
524 std::array<double,3> getEclCentroid(
const int& idx)
const;
594 const std::vector<std::array<int,3>>& cells_per_dim_vec,
595 const std::vector<int>& assignRefinedLevel,
596 const std::vector<std::string>& lgr_name_vec,
597 const std::vector<std::array<int,3>>& startIJK_vec = std::vector<std::array<int,3>>{},
598 const std::vector<std::array<int,3>>& endIJK_vec = std::vector<std::array<int,3>>{});
605 void updateCornerHistoryLevels(
const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
606 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
607 const std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
608 const int& corner_count,
609 const std::vector<std::array<int,2>>& preAdaptGrid_corner_history,
610 const int& preAdaptMaxLevel,
611 const int& newLevels);
613 void globalIdsPartitionTypesLgrAndLeafGrids(
const std::vector<int>& assignRefinedLevel,
614 const std::vector<std::array<int,3>>& cells_per_dim_vec,
615 const std::vector<int>& lgr_with_at_least_one_active_cell);
623 void getFirstChildGlobalIds([[maybe_unused]] std::vector<int>& parentToFirstChildGlobalIds);
645 void computeGlobalCellLgr(
const int& level,
const std::array<int,3>& startIJK, std::vector<int>& global_cell_lgr);
650 void computeGlobalCellLeafGridViewWithLgrs(std::vector<int>& global_cell_leaf);
666 void markElemAssignLevelDetectActiveLgrs(
const std::vector<std::array<int,3>>& startIJK_vec,
667 const std::vector<std::array<int,3>>& endIJK_vec,
668 std::vector<int>& assignRefinedLevel,
669 std::vector<int>& lgr_with_at_least_one_active_cell);
672 void populateCellIndexSetRefinedGrid(
int level);
675 void populateCellIndexSetLeafGridView();
678 void populateLeafGlobalIdSet();
707 void setPartitioningParams(
const std::map<std::string,std::string>& params);
722 double imbalanceTol = 1.1,
756 double imbalanceTol = 1.1,
805 std::pair<bool,std::vector<std::pair<std::string,bool>>>
807 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
808 const double* transmissibilities =
nullptr,
813 false, transmissibilities,
false,
814 overlapLayers, partitionMethod, 1.1,
false,
853 std::pair<bool,std::vector<std::pair<std::string,bool>>>
855 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
856 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
857 bool addCornerCells=
false,
int overlapLayers=1,
859 double imbalanceTol = 1.1,
862 return scatterGrid(method, ownersFirst, wells, possibleFutureConnections,
false,
863 transmissibilities, addCornerCells, overlapLayers, partitionMethod, imbalanceTol,
895 template<
class DataHandle>
896 std::pair<bool, std::vector<std::pair<std::string,bool> > >
898 const std::vector<cpgrid::OpmWellType> * wells,
899 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
900 const double* transmissibilities =
nullptr,
901 int overlapLayers=1,
int partitionMethod = 1,
int level =-1)
903 auto ret =
loadBalance(wells, possibleFutureConnections, transmissibilities, overlapLayers, partitionMethod, level);
946 template<
class DataHandle>
947 std::pair<bool, std::vector<std::pair<std::string,bool> > >
949 const std::vector<cpgrid::OpmWellType> * wells,
950 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
951 bool serialPartitioning,
952 const double* transmissibilities =
nullptr,
bool ownersFirst=
false,
954 double imbalanceTol = 1.1,
955 bool allowDistributedWells =
false)
957 auto ret = scatterGrid(method, ownersFirst, wells, possibleFutureConnections, serialPartitioning, transmissibilities,
958 addCornerCells, overlapLayers, partitionMethod, imbalanceTol, allowDistributedWells,
959 std::vector<int>{}, 0);
992 template<
class DataHandle>
993 std::pair<bool, std::vector<std::pair<std::string,bool> > >
995 const std::vector<cpgrid::OpmWellType> * wells,
996 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections = {},
997 bool ownersFirst=
false,
998 bool addCornerCells=
false,
int overlapLayers=1)
1002 possibleFutureConnections,
1023 template<
class DataHandle>
1029 bool ret =
loadBalance(overlapLayers, partitionMethod);
1048 bool loadBalance(
const std::vector<int>& parts,
bool ownersFirst=
false,
1049 bool addCornerCells=
false,
int overlapLayers=1)
1073 template<
class DataHandle>
1074 bool loadBalance(DataHandle& data,
const std::vector<int>& parts,
bool ownersFirst=
false,
1075 bool addCornerCells=
false,
int overlapLayers=1)
1077 bool ret =
loadBalance(parts, ownersFirst, addCornerCells, overlapLayers);
1100 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1101 const double* transmissibilities,
1103 const double imbalanceTol)
const;
1112 template<
class DataHandle>
1113 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir,
int )
const
1125 template<
class DataHandle>
1126 void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir)
const;
1142 typedef Dune::FieldVector<double, 3> Vector;
1145 const std::vector<double>& zcornData()
const;
1153 int numCells(
int level = -1)
const;
1159 int numFaces(
int level = -1)
const;
1181 int cellFace(
int cell,
int local_index,
int level = -1)
const;
1185 const cpgrid::OrientedEntityTable<0,1>::row_type
cellFaceRow(
int cell)
const;
1206 int faceCell(
int face,
int local_index,
int level = -1)
const;
1216 int numFaceVertices(
int face)
const;
1222 int faceVertex(
int face,
int local_index)
const;
1231 const Vector faceAreaNormalEcl(
int face)
const;
1265 :
public RandomAccessIteratorFacade<CentroidIterator<codim>,
1266 FieldVector<double, 3>,
1267 const FieldVector<double, 3>&, int>
1271 typedef typename std::vector<
cpgrid::Geometry<3-codim, 3> >::const_iterator
1279 const FieldVector<double,3>& dereference()
const
1281 return iter_->center();
1287 const FieldVector<double,3>& elementAt(
int n)
1289 return iter_[n]->center();
1291 void advance(
int n){
1317 int boundaryId(
int face)
const;
1325 template<
class Cell2FacesRowIterator>
1327 faceTag(
const Cell2FacesRowIterator& cell_face)
const;
1349 template<
class DataHandle>
1358 template<
class DataHandle>
1407 using ParallelIndexSet = cpgrid::CpGridDataTraits::ParallelIndexSet;
1409 using RemoteIndices = cpgrid::CpGridDataTraits::RemoteIndices;
1412 using CommunicationType = cpgrid::CpGridDataTraits::CommunicationType;
1417 const CommunicationType& cellCommunication()
const;
1419 ParallelIndexSet& getCellIndexSet();
1421 RemoteIndices& getCellRemoteIndices();
1423 const ParallelIndexSet& getCellIndexSet()
const;
1425 const RemoteIndices& getCellRemoteIndices()
const;
1463 std::pair<bool, std::vector<std::pair<std::string,bool> > >
1466 const std::vector<cpgrid::OpmWellType> * wells,
1467 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
1468 bool serialPartitioning,
1469 const double* transmissibilities,
1470 bool addCornerCells,
1473 double imbalanceTol = 1.1,
1474 bool allowDistributedWells =
true,
1475 const std::vector<int>& input_cell_part = {},
1482 std::vector<std::shared_ptr<cpgrid::CpGridData>> data_;
1484 std::vector<std::shared_ptr<cpgrid::CpGridData>> distributed_data_;
1486 std::vector<std::shared_ptr<cpgrid::CpGridData>>* current_data_;
1488 std::map<std::string,int> lgr_names_ = {{
"GLOBAL", 0}};
1494 std::shared_ptr<InterfaceMap> cell_scatter_gather_interfaces_;
1500 std::shared_ptr<InterfaceMap> point_scatter_gather_interfaces_;
1504 std::shared_ptr<cpgrid::GlobalIdSet> global_id_set_ptr_;
1510 std::map<std::string,std::string> partitioningParams;
1516#include <opm/grid/cpgrid/Entity.hpp>
1517#include <opm/grid/cpgrid/Iterators.hpp>
1518#include <opm/grid/cpgrid/CpGridData.hpp>
1524 namespace Capabilities
1530 static const bool v =
true;
1537 static const bool v =
true;
1543 static const bool v =
true;
1549 static const bool v =
true;
1556 static const bool v =
false;
1561 template<
class DataHandle>
1564 current_data_->back()->
communicate(data, iftype, dir);
1568 template<
class DataHandle>
1572 if (distributed_data_.empty()) {
1573 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balanced grid!");
1575 distributed_data_[0]->scatterData(handle, data_[0].get(),
1576 distributed_data_[0].get(),
1583 template<
class DataHandle>
1587 if (distributed_data_.empty()) {
1588 OPM_THROW(std::runtime_error,
"Moving Data only allowed with a load balance grid!");
1590 distributed_data_[0]->gatherData(handle, data_[0].get(), distributed_data_[0].get());
1596 template<
class Cell2FacesRowIterator>
1609 const int cell = cell_face.getCellIndex();
1610 const int face = *cell_face;
1611 assert (0 <= cell); assert (cell <
numCells());
1612 assert (0 <= face); assert (face <
numFaces());
1614 typedef cpgrid::OrientedEntityTable<1,0>::row_type F2C;
1617 const F2C& f2c = current_data_->back()->face_to_cell_[f];
1618 const face_tag tag = current_data_->back()->face_tag_[f];
1620 assert ((f2c.size() == 1) || (f2c.size() == 2));
1622 int inside_cell = 0;
1624 if ( f2c.size() == 2 )
1626 if ( f2c[1].index() == cell )
1631 const bool normal_is_in = ! f2c[inside_cell].orientation();
1636 return normal_is_in ? 0 : 1;
1639 return normal_is_in ? 2 : 3;
1643 return normal_is_in ? 4 : 5;
1648 OPM_THROW(std::logic_error,
"Unhandled face tag. This should never happen!");
1657#include <opm/grid/cpgrid/PersistentContainer.hpp>
1658#include <opm/grid/cpgrid/CartesianIndexMapper.hpp>
1659#include "cpgrid/Intersection.hpp"
1660#include "cpgrid/Geometry.hpp"
1661#include "cpgrid/Indexsets.hpp"
An iterator over the centroids of the geometry of the entities.
Definition CpGrid.hpp:1268
CentroidIterator(GeometryIterator iter)
Constructs a new iterator from an iterator over the geometries.
Definition CpGrid.hpp:1275
std::vector< cpgrid::Geometry< 3-codim, 3 > >::const_iterator GeometryIterator
The type of the iterator over the codim geometries.
Definition CpGrid.hpp:1272
[ 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:994
std::string name() const
Get the grid name.
Definition CpGrid.cpp:762
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:727
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:1112
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:1584
int numFaces(int level=-1) const
Get the number of faces.
Definition CpGrid.cpp:1154
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition CpGrid.cpp:1019
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:737
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData() const
Returns either data_ or distributed_data_(if non empty).
Definition CpGrid.cpp:662
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:1014
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:948
int faceTag(const Cell2FacesRowIterator &cell_face) const
Get the cartesian tag associated with a face tag.
Definition CpGrid.hpp:1598
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGrid.cpp:757
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:854
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:583
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:1090
cpgrid::CpGridDataTraits::InterfaceMap InterfaceMap
The type of the map describing communication interfaces.
Definition CpGrid.hpp:1362
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:897
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:1222
const Dune::cpgrid::CpGridData & currentLeafData() const
Returns current view data (the leaf grid).
Definition CpGrid.cpp:671
const std::vector< int > & globalCell() const
Retrieve mapping from internal ("compressed") active grid cells to external ("uncompressed") cells.
Definition CpGrid.cpp:681
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:753
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:720
int maxLevel() const
Return maximum level defined in this grid. Levels are 0 and 1, maxlevel = 1 (not counting leafview),...
Definition CpGrid.cpp:767
const CpGridTraits::Communication & comm() const
Get the collective communication object.
Definition CpGrid.cpp:1137
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:1024
void postAdapt()
Clean up refinement markers - set every element to the mark 0 which represents 'doing nothing'.
Definition CpGrid.cpp:2383
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:990
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the global grid.
Definition CpGrid.cpp:655
int faceVertex(int face, int local_index) const
Get the index identifying a vertex of a face.
Definition CpGrid.cpp:1232
bool mark(int refCount, const cpgrid::Entity< 0 > &element, bool throwOnFailure=false)
------------— Adaptivity (begin) ------------—
Definition CpGrid.cpp:1774
int numVertices() const
Get The number of vertices.
Definition CpGrid.cpp:1160
Dune::cpgrid::Intersection getParentIntersectionFromLgrBoundaryFace(const Dune::cpgrid::Intersection &intersection) const
Definition CpGrid.cpp:1237
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:747
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:1741
void syncDistributedGlobalCellIds()
Synchronizes cell global ids across processes after load balancing.
Definition CpGrid.cpp:2528
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGrid.cpp:752
int faceCell(int face, int local_index, int level=-1) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1184
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
void globalRefine(int refCount, bool throwOnFailure=false)
Refine the grid refCount times using the default refinement rule.
Definition CpGrid.cpp:1036
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim and PartitionIteratorType.
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:2591
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:1024
unsigned int ghostSize(int) const
Size of the ghost cell layer on the leaf level.
Definition CpGrid.cpp:1096
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:806
cpgrid::Entity< codim > entity(const cpgrid::Entity< codim > &seed) const
given an EntitySeed (or EntityPointer) return an entity object
Definition CpGrid.cpp:1080
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition CpGrid.cpp:1031
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:166
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:1810
void switchToDistributedView()
Switch to the distributed view.
Definition CpGrid.cpp:1665
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:1798
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
Definition CpGrid.cpp:1793
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:1074
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:1179
const InterfaceMap & cellScatterGatherInterface() const
Get an interface for gathering/scattering data attached to cells with communication.
Definition CpGrid.cpp:1650
bool refineAndUpdateGrid(bool throwOnFailure, 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:1845
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:1048
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:189
int numCells(int level=-1) const
Get the number of cells.
Definition CpGrid.cpp:1148
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 communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir, int) const
The new communication interface.
Definition CpGrid.hpp:1113
int cellFace(int cell, int local_index, int level=-1) const
Get a specific face of a cell.
Definition CpGrid.cpp:1172
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:2748
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition CpGrid.hpp:1569
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:98
Definition Entity.hpp:107
This class encapsulates geometry for vertices, intersections, and cells.
Definition Geometry.hpp:76
The global id set for Dune.
Definition Indexsets.hpp:483
Only needs to provide interface for doing nothing.
Definition Iterators.hpp:118
Definition Indexsets.hpp:199
Definition Indexsets.hpp:57
Definition Intersection.hpp:276
Definition Intersection.hpp:63
Iterator intended to be used as LeafIterator and LevelIterator (no difference due to no adaptivity) f...
Definition Iterators.hpp:60
Definition Indexsets.hpp:367
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:71
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