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>
58 bool operator()(
const std::pair<int,int>& o,
const std::pair<int,int>& v)
60 return o.first < v.first;
73 const std::array<int, 3>& initial_split,
75 std::vector<int>& cell_part,
76 bool recursive =
false,
77 bool ensureConnectivity =
true);
87 void addOverlapLayer(
const CpGrid& grid,
88 const std::vector<int>& cell_part,
89 std::vector<std::set<int> >& cell_overlap,
90 int mypart,
int overlapLayers,
bool all=
false);
104 int addOverlapLayer(
const CpGrid& grid,
const std::vector<int>& cell_part,
105 std::vector<std::tuple<int,int,char>>& exportList,
106 std::vector<std::tuple<int,int,char,int>>& importList,
107 const Communication<Dune::MPIHelper::MPICommunicator>& cc,
108 bool addCornerCells,
const double* trans,
int layers = 1);
131 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
132 std::vector<std::tuple<int,int,char> >,
133 std::vector<std::tuple<int,int,char,int> >,
135 createZoltanListsFromParts(
const CpGrid& grid,
const std::vector<cpgrid::OpmWellType> * wells,
136 const double* transmissibilities,
const std::vector<int>& parts,
137 bool allowDistributedWells);
156 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
157 std::vector<std::tuple<int,int,char> >,
158 std::vector<std::tuple<int,int,char,int> >,
160 vanillaPartitionGridOnRoot(
const CpGrid& grid,
161 const std::vector<cpgrid::OpmWellType> * wells,
162 const double* transmissibilities,
163 bool allowDistributedWells);
[ provides Dune::Grid ]
Definition CpGrid.hpp:238
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:57