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;
76 const std::array<int, 3>& initial_split,
78 std::vector<int>& cell_part,
79 bool recursive =
false,
80 bool ensureConnectivity =
true);
92 void addOverlapLayer(
const CpGrid& grid,
93 const std::vector<int>& cell_part,
94 std::vector<std::set<int> >& cell_overlap,
113 int addOverlapLayer(
const CpGrid& grid,
114 const std::vector<int>& cell_part,
115 std::vector<std::tuple<int,int,char>>& exportList,
116 std::vector<std::tuple<int,int,char,int>>& importList,
117 const Communication<Dune::MPIHelper::MPICommunicator>& cc,
150 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
151 std::vector<std::tuple<int,int,char> >,
152 std::vector<std::tuple<int,int,char,int> >,
154 createListsFromParts(
const CpGrid& grid,
155 const std::vector<cpgrid::OpmWellType> * wells,
156 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
157 const double* transmissibilities,
158 const std::vector<int>& parts,
159 bool allowDistributedWells,
160 std::shared_ptr<cpgrid::CombinedGridWellGraph> gridAndWells =
nullptr,
183 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
184 std::vector<std::tuple<int,int,char> >,
185 std::vector<std::tuple<int,int,char,int> >,
187 vanillaPartitionGridOnRoot(
const CpGrid& grid,
188 const std::vector<cpgrid::OpmWellType> * wells,
189 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
190 const double* transmissibilities,
191 bool allowDistributedWells);
[ provides Dune::Grid ]
Definition CpGrid.hpp:203
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
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:198
Definition GridPartitioning.hpp:58