Skip to content

Commit 1affdfa

Browse files
author
Lisa Julia Nebel
authored
Merge pull request #874 from michal-toth/feature/graph-of-grid-second-try
Allow distributed wells with partition method `zoltanwell`
2 parents 7970508 + 643e81c commit 1affdfa

3 files changed

Lines changed: 80 additions & 29 deletions

File tree

opm/grid/GraphOfGridWrappers.cpp

Lines changed: 76 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -27,6 +27,7 @@
2727
#include "GraphOfGridWrappers.hpp"
2828
#include <opm/grid/common/CommunicationUtils.hpp>
2929
#include <opm/grid/common/ZoltanPartition.hpp> // function scatterExportInformation
30+
#include <opm/grid/common/ZoltanGraphFunctions.hpp> // makeImportAndExportLists when allowDistributedWells==true
3031

3132
namespace Opm {
3233

@@ -524,6 +525,7 @@ zoltanPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
524525
Dune::EdgeWeightMethod edgeWeightMethod,
525526
int root,
526527
const double zoltanImbalanceTol,
528+
bool allowDistributedWells,
527529
const std::map<std::string, std::string>& params,
528530
int level)
529531
{
@@ -562,8 +564,11 @@ zoltanPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
562564
assert(gog.size()==0 || !partitionIsEmpty);
563565
auto wellConnections = partitionIsEmpty || !wells ? Dune::cpgrid::WellConnections()
564566
: Dune::cpgrid::WellConnections(*wells, possibleFutureConnections, grid);
565-
addWellConnections(gog, wellConnections);
566-
gog.addNeighboringCellsToWells(layers);
567+
if (!allowDistributedWells){
568+
// skip cell contraction if wells can be distributed over multiple processes
569+
addWellConnections(gog, wellConnections);
570+
gog.addNeighboringCellsToWells(layers);
571+
}
567572

568573
// call partitioner
569574
setGraphOfGridZoltanGraphFunctions(zz, gog, partitionIsEmpty);
@@ -590,30 +595,63 @@ zoltanPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
590595
}
591596

592597
// arrange output into tuples and add well cells
593-
auto importExportLists = makeImportAndExportLists(gog,
594-
cc,
595-
wells,
596-
wellConnections,
597-
root,
598-
numExport,
599-
numImport,
600-
exportLocalGids,
601-
exportGlobalGids,
602-
exportProcs,
603-
importGlobalGids,
604-
level);
598+
auto prepareIELists = [&]() {
599+
if (allowDistributedWells) {
600+
// wells can be split among several processes
601+
using CombinedGridWellGraph = Dune::cpgrid::CombinedGridWellGraph;
602+
std::shared_ptr<CombinedGridWellGraph> gridAndWells;
603+
if (wells) {
604+
gridAndWells.reset(new CombinedGridWellGraph(grid,
605+
wells,
606+
possibleFutureConnections,
607+
transmissibilities,
608+
partitionIsEmpty,
609+
edgeWeightMethod));
610+
}
611+
auto result = makeImportAndExportLists(grid,
612+
cc,
613+
wells,
614+
possibleFutureConnections,
615+
gridAndWells.get(),
616+
root,
617+
numExport,
618+
numImport,
619+
exportGlobalGids, // function uses Local GIDs that are identical to global GIDs
620+
exportGlobalGids,
621+
exportProcs,
622+
importGlobalGids,
623+
allowDistributedWells);
624+
std::sort(std::get<2>(result).begin(), std::get<2>(result).end());
625+
std::sort(std::get<3>(result).begin(), std::get<3>(result).end());
626+
return result;
627+
} else {
628+
// each well is guaranteed to be on a single process
629+
auto partResult = makeImportAndExportLists(gog,
630+
cc,
631+
wells,
632+
wellConnections,
633+
root,
634+
numExport,
635+
numImport,
636+
exportLocalGids,
637+
exportGlobalGids,
638+
exportProcs,
639+
importGlobalGids,
640+
level);
641+
return std::tuple(std::move(std::get<0>(partResult)),
642+
std::move(std::get<1>(partResult)),
643+
std::move(std::get<2>(partResult)),
644+
std::move(std::get<3>(partResult)),
645+
std::move(wellConnections));
646+
}
647+
};
648+
auto importExportLists = prepareIELists();
605649

606650
Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, &exportProcs, &exportToPart);
607651
Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, &importProcs, &importToPart);
608652
Zoltan_Destroy(&zz);
609653

610-
// add wellConnections to the importExportLists and return it
611-
auto result = std::tuple(std::move(std::get<0>(importExportLists)),
612-
std::move(std::get<1>(importExportLists)),
613-
std::move(std::get<2>(importExportLists)),
614-
std::move(std::get<3>(importExportLists)),
615-
std::move(wellConnections));
616-
return result;
654+
return importExportLists;
617655
}
618656

619657
std::vector<std::vector<int> >
@@ -640,6 +678,7 @@ applySerialZoltan (const Dune::CpGrid& grid,
640678
Dune::EdgeWeightMethod edgeWeightMethod,
641679
int root,
642680
const double zoltanImbalanceTol,
681+
bool allowDistributedWells,
643682
const std::map<std::string, std::string>& params)
644683
{
645684
int rc = ZOLTAN_OK;
@@ -673,8 +712,11 @@ applySerialZoltan (const Dune::CpGrid& grid,
673712

674713
// prepare graph and contract well cells
675714
GraphOfGrid gog(grid, transmissibilities, edgeWeightMethod);
676-
addWellConnections(gog, wellConnections);
677-
gog.addNeighboringCellsToWells(layers);
715+
if (!allowDistributedWells){
716+
// skip cell contraction if wells can be distributed over multiple processes
717+
addWellConnections(gog, wellConnections);
718+
gog.addNeighboringCellsToWells(layers);
719+
}
678720

679721
// call partitioner
680722
setGraphOfGridZoltanGraphFunctions(zz, gog, false);
@@ -723,6 +765,7 @@ zoltanSerialPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
723765
Dune::EdgeWeightMethod edgeWeightMethod,
724766
int root,
725767
const double zoltanImbalanceTol,
768+
bool allowDistributedWells,
726769
const std::map<std::string, std::string>& params)
727770
{
728771
// root process has the whole grid, other ranks nothing
@@ -744,6 +787,7 @@ zoltanSerialPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
744787
edgeWeightMethod,
745788
root,
746789
zoltanImbalanceTol,
790+
allowDistributedWells,
747791
params);
748792
}
749793

@@ -788,8 +832,15 @@ zoltanSerialPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
788832
// get the distribution of wells
789833
std::vector<std::pair<std::string, bool>> parallel_wells;
790834
if (wells) {
791-
auto wellRanks = getWellRanks(gIDtoRank, wellConnections);
792-
parallel_wells = wellsOnThisRank(*wells, wellRanks, cc, root);
835+
if (allowDistributedWells) {
836+
// wells can be split among several processes
837+
auto wellsOnProc = Dune::cpgrid::perforatingWellIndicesOnProc(gIDtoRank, *wells, possibleFutureConnections, grid);
838+
parallel_wells = Dune::cpgrid::computeParallelWells(wellsOnProc, *wells, cc, root);
839+
} else {
840+
// each well is guaranteed to be on a single process
841+
auto wellRanks = getWellRanks(gIDtoRank, wellConnections);
842+
parallel_wells = wellsOnThisRank(*wells, wellRanks, cc, root);
843+
}
793844
}
794845

795846
return std::make_tuple(std::move(gIDtoRank),

opm/grid/GraphOfGridWrappers.hpp

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -269,6 +269,7 @@ zoltanPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
269269
Dune::EdgeWeightMethod edgeWeightMethod,
270270
int root,
271271
const double zoltanImbalanceTol,
272+
bool allowDistributedWells,
272273
const std::map<std::string,std::string>& params,
273274
int level);
274275

@@ -298,6 +299,7 @@ zoltanSerialPartitioningWithGraphOfGrid(const Dune::CpGrid& grid,
298299
Dune::EdgeWeightMethod edgeWeightMethod,
299300
int root,
300301
const double zoltanImbalanceTol,
302+
bool allowDistributedWells,
301303
const std::map<std::string,std::string>& params);
302304
#endif // HAVE_MPI
303305

opm/grid/cpgrid/CpGrid.cpp

Lines changed: 2 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -174,7 +174,6 @@ CpGrid::CpGrid()
174174
current_view_data_ = data_[0].get();
175175
current_data_ = &data_;
176176
global_id_set_ptr_ = std::make_shared<cpgrid::GlobalIdSet>(*current_view_data_);
177-
178177
}
179178

180179
CpGrid::CpGrid(MPIHelper::MPICommunicator comm)
@@ -188,7 +187,6 @@ CpGrid::CpGrid(MPIHelper::MPICommunicator comm)
188187
current_view_data_ = data_[0].get();
189188
current_data_ = &data_;
190189
global_id_set_ptr_ = std::make_shared<cpgrid::GlobalIdSet>(*current_view_data_);
191-
192190
}
193191

194192
std::vector<int>
@@ -381,8 +379,8 @@ CpGrid::scatterGrid(EdgeWeightMethod method,
381379
#ifdef HAVE_ZOLTAN
382380
std::tie(computedCellPart, wells_on_proc, exportList, importList, wellConnections)
383381
= serialPartitioning
384-
? Opm::zoltanSerialPartitioningWithGraphOfGrid(*this, wells, possibleFutureConnections, transmissibilities, cc, method, 0, imbalanceTol, partitioningParams)
385-
: Opm::zoltanPartitioningWithGraphOfGrid(*this, wells, possibleFutureConnections, transmissibilities, cc, method, 0, imbalanceTol, partitioningParams, level);
382+
? Opm::zoltanSerialPartitioningWithGraphOfGrid(*this, wells, possibleFutureConnections, transmissibilities, cc, method, 0, imbalanceTol, allowDistributedWells, partitioningParams)
383+
: Opm::zoltanPartitioningWithGraphOfGrid(*this, wells, possibleFutureConnections, transmissibilities, cc, method, 0, imbalanceTol, allowDistributedWells, partitioningParams, level);
386384
#else
387385
OPM_THROW(std::runtime_error, "Parallel runs depend on ZOLTAN if useZoltan is true. Please install!");
388386
#endif // HAVE_ZOLTAN

0 commit comments

Comments
 (0)