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
30
33#include <opm/grid/common/ZoltanGraphFunctions.hpp> // defines Zoltan and null-callback-functions
34
35namespace Opm {
36/*
37 This file contains wrappers for GraphOfGrid that satisfy interface
38 requirements of graph partitioners like Zoltan and (TODO!) Metis.
39
40 Additionally, parsing wells is done here.
41*/
42
43#if HAVE_MPI
47int getGraphOfGridNumVertices(void* pGraph, int *err);
48
54 [[maybe_unused]] int dimGlobalID,
55 [[maybe_unused]] int dimLocalID,
56 ZOLTAN_ID_PTR gIDs,
57 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
58 int weightDim,
59 float *objWeights,
60 int *err);
61
66void getGraphOfGridNumEdges(void *pGraph,
67 [[maybe_unused]] int dimGlobalID,
68 [[maybe_unused]] int dimLocalID,
69 int numCells,
70 ZOLTAN_ID_PTR gIDs,
71 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
72 int *numEdges,
73 int *err);
74
82void getGraphOfGridEdgeList(void *pGraph,
83 [[maybe_unused]] int dimGlobalID,
84 [[maybe_unused]] int dimLocalID,
85 int numCells,
86 ZOLTAN_ID_PTR gIDs,
87 [[maybe_unused]] ZOLTAN_ID_PTR lIDs,
88 int *numEdges,
89 ZOLTAN_ID_PTR nborGIDs,
90 int *nborProc,
91 int weightDim,
92 float *edgeWeights,
93 int *err);
94
96template<typename Zoltan_Struct>
99 bool pretendNull);
100#endif
101
112 const std::unordered_map<std::string, std::set<int>>& wells,
113 bool checkWellIntersections=true);
114
123 bool checkWellIntersections=true);
124
131 std::vector<int>& gIDtoRank,
132 const int& root = -1);
133
134#if HAVE_MPI
135namespace Impl{
139void extendAndSortImportList(std::vector<std::tuple<int,int,char,int>>& importList,
140 const std::vector<int>& extraCells);
141
149std::vector<std::vector<int>>
151 std::vector<std::tuple<int,int,char>>& exportList,
152 int root,
153 const std::vector<int>& gIDtoRank);
172std::vector<int> communicateExportedCells(const std::vector<std::vector<int>>& exportedCells,
174 int root);
175} // end namespace Impl
176
187 int root,
188 std::vector<std::tuple<int, int, char>>& exportList,
189 std::vector<std::tuple<int, int, char, int>>& importList,
190 const std::vector<int>& gIDtoRank = {});
191#endif // HAVE_MPI
192
198std::vector<int> getWellRanks(const std::vector<int>& gIDtoRank,
199 const Dune::cpgrid::WellConnections& wellConnections);
200
201#if HAVE_MPI
214std::vector<std::pair<std::string, bool>>
215wellsOnThisRank(const std::vector<Dune::cpgrid::OpmWellType>& wells,
216 const std::vector<int>& wellRanks,
218 int root);
219
239template<class Id>
240std::tuple<std::vector<int>,
241 std::vector<std::pair<std::string, bool>>,
242 std::vector<std::tuple<int,int,char> >,
243 std::vector<std::tuple<int,int,char,int> > >
245 const Dune::Communication<MPI_Comm>& cc,
246 const std::vector<Dune::cpgrid::OpmWellType> * wells,
247 const Dune::cpgrid::WellConnections& wellConnections,
248 int root,
249 int numExport,
250 int numImport,
251 [[maybe_unused]] const Id* exportLocalGids,
252 const Id* exportGlobalGids,
253 const int* exportToPart,
254 const Id* importGlobalGids,
255 int level);
256
261std::tuple<std::vector<int>, std::vector<std::pair<std::string, bool>>,
262 std::vector<std::tuple<int,int,char> >,
263 std::vector<std::tuple<int,int,char,int> >,
266 const std::vector<Dune::cpgrid::OpmWellType> * wells,
267 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
268 const double* transmissibilities,
270 Dune::EdgeWeightMethod edgeWeightMethod,
271 int root,
272 const double zoltanImbalanceTol,
273 bool allowDistributedWells,
274 const std::map<std::string,std::string>& params,
275 int level);
276
284std::vector<std::vector<int> >
285makeExportListsFromGIDtoRank(const std::vector<int>& gIDtoRank, int ccsize);
286
291std::tuple<std::vector<int>, std::vector<std::pair<std::string, bool>>,
292 std::vector<std::tuple<int,int,char> >,
293 std::vector<std::tuple<int,int,char,int> >,
296 const std::vector<Dune::cpgrid::OpmWellType> * wells,
297 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
298 const double* transmissibilities,
300 Dune::EdgeWeightMethod edgeWeightMethod,
301 int root,
302 const double zoltanImbalanceTol,
303 bool allowDistributedWells,
304 const std::map<std::string,std::string>& params);
305#endif // HAVE_MPI
306
307} // end namespace Opm
308
309#endif // GRAPH_OF_GRID_WRAPPERS_HEADER
[ provides Dune::Grid ]
Definition: CpGrid.hpp:203
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