GraphOfGridWrappers.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 Copyright 2024 Equinor ASA.
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
25
26#ifndef GRAPH_OF_GRID_WRAPPERS_HEADER
27#define GRAPH_OF_GRID_WRAPPERS_HEADER
28
29#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/grid/common/ZoltanGraphFunctions.hpp> // defines Zoltan and null-callback-functions
33
34namespace Opm {
35/*
36 This file contains wrappers for GraphOfGrid that satisfy interface
37 requirements of graph partitioners like Zoltan and (TODO!) Metis.
38
39 Additionally, parsing wells is done here.
40*/
41
42#if HAVE_MPI
46int getGraphOfGridNumVertices(void* pGraph, int *err);
47
53 [[maybe_unused]] int dimGlobalID,
54 [[maybe_unused]] int dimLocalID,
55 ZOLTAN_ID_PTR gIDs,
56 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
57 int weightDim,
58 float *objWeights,
59 int *err);
60
65void getGraphOfGridNumEdges(void *pGraph,
66 [[maybe_unused]] int dimGlobalID,
67 [[maybe_unused]] int dimLocalID,
68 int numCells,
69 ZOLTAN_ID_PTR gIDs,
70 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
71 int *numEdges,
72 int *err);
73
81void getGraphOfGridEdgeList(void *pGraph,
82 [[maybe_unused]] int dimGlobalID,
83 [[maybe_unused]] int dimLocalID,
84 int numCells,
85 ZOLTAN_ID_PTR gIDs,
86 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
87 int *numEdges,
88 ZOLTAN_ID_PTR nborGIDs,
89 int *nborProc,
90 int weightDim,
91 float *edgeWeights,
92 int *err);
93
95template<typename Zoltan_Struct>
98 bool pretendNull);
99#endif
100
111 const std::unordered_map<std::string, std::set<int>>& wells,
112 bool checkWellIntersections=true);
113
122 bool checkWellIntersections=true);
123
130 std::vector<int>& gIDtoRank,
131 const int& root = -1);
132
133#if HAVE_MPI
134namespace Impl{
138void extendAndSortImportList(std::vector<std::tuple<int,int,char,int>>& importList,
139 const std::vector<int>& extraCells);
140
148std::vector<std::vector<int>>
150 std::vector<std::tuple<int,int,char>>& exportList,
151 int root,
152 const std::vector<int>& gIDtoRank);
171std::vector<int> communicateExportedCells(const std::vector<std::vector<int>>& exportedCells,
173 int root);
174} // end namespace Impl
175
186 int root,
187 std::vector<std::tuple<int, int, char>>& exportList,
188 std::vector<std::tuple<int, int, char, int>>& importList,
189 const std::vector<int>& gIDtoRank = {});
190#endif // HAVE_MPI
191
197std::vector<int> getWellRanks(const std::vector<int>& gIDtoRank,
198 const Dune::cpgrid::WellConnections& wellConnections);
199
200#if HAVE_MPI
213std::vector<std::pair<std::string, bool>>
214wellsOnThisRank(const std::vector<Dune::cpgrid::OpmWellType>& wells,
215 const std::vector<int>& wellRanks,
217 int root);
218
238template<class Id>
239std::tuple<std::vector<int>,
240 std::vector<std::pair<std::string, bool>>,
241 std::vector<std::tuple<int,int,char> >,
242 std::vector<std::tuple<int,int,char,int> > >
244 const Dune::Communication<MPI_Comm>& cc,
245 const std::vector<Dune::cpgrid::OpmWellType> * wells,
246 const Dune::cpgrid::WellConnections& wellConnections,
247 int root,
248 int numExport,
249 int numImport,
250 [[maybe_unused]] const Id* exportLocalGids,
251 const Id* exportGlobalGids,
252 const int* exportToPart,
253 const Id* importGlobalGids,
254 int level);
255
260std::tuple<std::vector<int>, std::vector<std::pair<std::string, bool>>,
261 std::vector<std::tuple<int,int,char> >,
262 std::vector<std::tuple<int,int,char,int> >,
265 const std::vector<Dune::cpgrid::OpmWellType> * wells,
266 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
267 const double* transmissibilities,
269 Dune::EdgeWeightMethod edgeWeightMethod,
270 int root,
271 const double zoltanImbalanceTol,
272 bool allowDistributedWells,
273 const std::map<std::string,std::string>& params,
274 int level);
275
283std::vector<std::vector<int> >
284makeExportListsFromGIDtoRank(const std::vector<int>& gIDtoRank, int ccsize);
285
290std::tuple<std::vector<int>, std::vector<std::pair<std::string, bool>>,
291 std::vector<std::tuple<int,int,char> >,
292 std::vector<std::tuple<int,int,char,int> >,
295 const std::vector<Dune::cpgrid::OpmWellType> * wells,
296 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
297 const double* transmissibilities,
299 Dune::EdgeWeightMethod edgeWeightMethod,
300 int root,
301 const double zoltanImbalanceTol,
302 bool allowDistributedWells,
303 const std::map<std::string,std::string>& params);
304#endif // HAVE_MPI
305
306} // end namespace Opm
307
308#endif // GRAPH_OF_GRID_WRAPPERS_HEADER
[ provides Dune::Grid ]
Definition: CpGrid.hpp:198
A class calculating and representing all connections of wells.
Definition: WellConnections.hpp:51
A class storing a graph representation of the grid.
Definition: GraphOfGrid.hpp:43
EdgeWeightMethod
enum for choosing Methods for weighting graph-edges correspoding to cell interfaces in Zoltan's or Me...
Definition: GridEnums.hpp:34
void extendAndSortImportList(std::vector< std::tuple< int, int, char, int > > &importList, const std::vector< int > &extraCells)
Add cells to the import list and sort it.
std::vector< std::vector< int > > extendRootExportList(const GraphOfGrid< Dune::CpGrid > &gog, std::vector< std::tuple< int, int, char > > &exportList, int root, const std::vector< int > &gIDtoRank)
Add well cells' global IDs to the root's export list and output cells missing in other rank's import ...
std::vector< int > communicateExportedCells(const std::vector< std::vector< int > > &exportedCells, const Dune::cpgrid::CpGridDataTraits::Communication &cc, int root)
Communicate cells exported from root, needed for extending other rank's import lists.
int numCells(const Dune::CpGrid &grid)
Get the number of cells of a grid.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
std::vector< int > getWellRanks(const std::vector< int > &gIDtoRank, const Dune::cpgrid::WellConnections &wellConnections)
Find to which ranks wells were assigned.
std::vector< std::vector< int > > makeExportListsFromGIDtoRank(const std::vector< int > &gIDtoRank, int ccsize)
Make complete export lists from a vector holding destination rank for each global ID.
std::tuple< std::vector< int >, std::vector< std::pair< std::string, bool > >, std::vector< std::tuple< int, int, char > >, std::vector< std::tuple< int, int, char, int > > > makeImportAndExportLists(const GraphOfGrid< Dune::CpGrid > &gog, const Dune::Communication< MPI_Comm > &cc, const std::vector< Dune::cpgrid::OpmWellType > *wells, const Dune::cpgrid::WellConnections &wellConnections, int root, int numExport, int numImport, const Id *exportLocalGids, const Id *exportGlobalGids, const int *exportToPart, const Id *importGlobalGids, int level)
Transform Zoltan output into tuples.
void addWellConnections(GraphOfGrid< Dune::CpGrid > &gog, const Dune::cpgrid::WellConnections &wells, bool checkWellIntersections=true)
Add WellConnections to the GraphOfGrid.
std::tuple< std::vector< int >, std::vector< std::pair< std::string, bool > >, std::vector< std::tuple< int, int, char > >, std::vector< std::tuple< int, int, char, int > >, Dune::cpgrid::WellConnections > zoltanPartitioningWithGraphOfGrid(const Dune::CpGrid &grid, const std::vector< Dune::cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const double *transmissibilities, const Dune::cpgrid::CpGridDataTraits::Communication &cc, Dune::EdgeWeightMethod edgeWeightMethod, int root, const double zoltanImbalanceTol, bool allowDistributedWells, const std::map< std::string, std::string > &params, int level)
Call Zoltan partitioner on GraphOfGrid.
void extendAndSortExportAndImportLists(const GraphOfGrid< Dune::CpGrid > &gog, const Dune::cpgrid::CpGridDataTraits::Communication &cc, int root, std::vector< std::tuple< int, int, char > > &exportList, std::vector< std::tuple< int, int, char, int > > &importList, const std::vector< int > &gIDtoRank={})
Add well cells' global IDs to the root's export and others' import list.
void addFutureConnectionWells(GraphOfGrid< Dune::CpGrid > &gog, const std::unordered_map< std::string, std::set< int > > &wells, bool checkWellIntersections=true)
Adds well to the GraphOfGrid.
std::tuple< std::vector< int >, std::vector< std::pair< std::string, bool > >, std::vector< std::tuple< int, int, char > >, std::vector< std::tuple< int, int, char, int > >, Dune::cpgrid::WellConnections > zoltanSerialPartitioningWithGraphOfGrid(const Dune::CpGrid &grid, const std::vector< Dune::cpgrid::OpmWellType > *wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const double *transmissibilities, const Dune::cpgrid::CpGridDataTraits::Communication &cc, Dune::EdgeWeightMethod edgeWeightMethod, int root, const double zoltanImbalanceTol, bool allowDistributedWells, const std::map< std::string, std::string > &params)
Call serial Zoltan partitioner on GraphOfGrid.
void getGraphOfGridVerticesList(void *pGraph, int dimGlobalID, int dimLocalID, ZOLTAN_ID_PTR gIDs, ZOLTAN_ID_PTR lIDs, int weightDim, float *objWeights, int *err)
callback function for ZOLTAN_OBJ_LIST_FN
void extendGIDtoRank(const GraphOfGrid< Dune::CpGrid > &gog, std::vector< int > &gIDtoRank, const int &root=-1)
Correct gIDtoRank's data about well cells.
void getGraphOfGridEdgeList(void *pGraph, int dimGlobalID, int dimLocalID, int numCells, ZOLTAN_ID_PTR gIDs, ZOLTAN_ID_PTR lIDs, int *numEdges, ZOLTAN_ID_PTR nborGIDs, int *nborProc, int weightDim, float *edgeWeights, int *err)
callback function for ZOLTAN_EDGE_LIST_MULTI_FN
void setGraphOfGridZoltanGraphFunctions(Zoltan_Struct *zz, GraphOfGrid< Dune::CpGrid > &gog, bool pretendNull)
Register callback functions to Zoltan.
void getGraphOfGridNumEdges(void *pGraph, int dimGlobalID, int dimLocalID, int numCells, ZOLTAN_ID_PTR gIDs, ZOLTAN_ID_PTR lIDs, int *numEdges, int *err)
callback function for ZOLTAN_NUM_EDGES_MULTI_FN
int getGraphOfGridNumVertices(void *pGraph, int *err)
callback function for ZOLTAN_NUM_OBJ_FN
std::vector< std::pair< std::string, bool > > wellsOnThisRank(const std::vector< Dune::cpgrid::OpmWellType > &wells, const std::vector< int > &wellRanks, const Dune::cpgrid::CpGridDataTraits::Communication &cc, int root)
Get rank-specific information about which wells are present.
Dune::Communication< MPICommunicator > Communication
Definition: CpGridDataTraits.hpp:58