38#ifndef OPM_GRIDPARTITIONING_HEADER
39#define OPM_GRIDPARTITIONING_HEADER
46#include <dune/common/parallel/mpihelper.hh>
48#include <opm/grid/utility/OpmWellType.hpp>
49#include <opm/grid/common/WellConnections.hpp>
50#include <opm/grid/common/ZoltanGraphFunctions.hpp>
59 bool operator()(
const std::pair<int,int>& o,
const std::pair<int,int>& v)
61 return o.first < v.first;
74 const std::array<int, 3>& initial_split,
76 std::vector<int>& cell_part,
77 bool recursive =
false,
78 bool ensureConnectivity =
true);
88 void addOverlapLayer(
const CpGrid& grid,
89 const std::vector<int>& cell_part,
90 std::vector<std::set<int> >& cell_overlap,
91 int mypart,
int overlapLayers,
bool all=
false);
105 int addOverlapLayer(
const CpGrid& grid,
const std::vector<int>& cell_part,
106 std::vector<std::tuple<int,int,char>>& exportList,
107 std::vector<std::tuple<int,int,char,int>>& importList,
108 const Communication<Dune::MPIHelper::MPICommunicator>& cc,
109 bool addCornerCells,
const double* trans,
int layers = 1);
135 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
136 std::vector<std::tuple<int,int,char> >,
137 std::vector<std::tuple<int,int,char,int> >,
139 createListsFromParts(
const CpGrid& grid,
const std::vector<cpgrid::OpmWellType> * wells,
140 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
141 const double* transmissibilities,
const std::vector<int>& parts,
142 bool allowDistributedWells, std::shared_ptr<cpgrid::CombinedGridWellGraph> gridAndWells =
nullptr);
164 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
165 std::vector<std::tuple<int,int,char> >,
166 std::vector<std::tuple<int,int,char,int> >,
168 vanillaPartitionGridOnRoot(
const CpGrid& grid,
169 const std::vector<cpgrid::OpmWellType> * wells,
170 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
171 const double* transmissibilities,
172 bool allowDistributedWells);
[ provides Dune::Grid ]
Definition CpGrid.hpp:201
Copyright 2019 Equinor AS.
Definition CartesianIndexMapper.hpp:10
void partition(const CpGrid &grid, const coord_t &initial_split, int &num_part, std::vector< int > &cell_part, bool recursive, bool ensureConnectivity)
Partition a CpGrid based on (ijk) coordinates, with splitting to ensure that each partition is connec...
Definition GridPartitioning.cpp:195
Definition GridPartitioning.hpp:58