Make test on one process pass by taking the possible future connections into account when creating the sparsity pattern of the jacobian and when creating the partitions for the solver

This commit is contained in:
Lisa Julia Nebel 2024-08-05 14:27:38 +02:00
parent f17ad653ec
commit 9a97a8d658
4 changed files with 70 additions and 36 deletions

View File

@ -953,9 +953,10 @@ private:
? param.num_local_domains_
: num_cells / default_cells_per_domain;
auto possibleFutureConnectionSet = this->model_.simulator().vanguard().schedule().getPossibleFutureConnections();
return ::Opm::partitionCells(param.local_domain_partition_method_,
num_domains,
grid.leafGridView(), wells, zoltan_ctrl);
grid.leafGridView(), wells, &possibleFutureConnectionSet, zoltan_ctrl);
}
std::vector<int> reconstitutePartitionVector() const

View File

@ -106,9 +106,10 @@ public:
/// \param[in] zoltan_ctrl Control parameters for on-rank subdomain
/// partitioning.
template <class GridView, class Element>
void buildLocalGraph(const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl);
void buildLocalGraph(const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl);
/// Partition rank's interior cells into non-overlapping domains using
/// the Zoltan graph partitioning software package.
@ -178,9 +179,10 @@ private:
/// \param[in] g2l Mapping from globally unique cell IDs to local,
/// on-rank active cell IDs. Return value from \c connectElements().
template <typename Comm>
void connectWells(const Comm comm,
const std::vector<Opm::Well>& wells,
const std::unordered_map<int, int>& g2l);
void connectWells(const Comm comm,
const std::vector<Opm::Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const std::unordered_map<int, int>& g2l);
};
// Note: "grid_view.size(0)" is intentional here. It is not an error. The
@ -194,11 +196,12 @@ ZoltanPartitioner::ZoltanPartitioner(const GridView&
{}
template <class GridView, class Element>
void ZoltanPartitioner::buildLocalGraph(const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl)
void ZoltanPartitioner::buildLocalGraph(const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl)
{
this->connectWells(grid_view.comm(), wells, this->connectElements(grid_view, zoltan_ctrl));
this->connectWells(grid_view.comm(), wells, possibleFutureConnections, this->connectElements(grid_view, zoltan_ctrl));
}
template <class GridView, class Element>
@ -282,9 +285,10 @@ ZoltanPartitioner::connectElements(const GridView&
}
template <typename Comm>
void ZoltanPartitioner::connectWells(const Comm comm,
const std::vector<Opm::Well>& wells,
const std::unordered_map<int, int>& g2l)
void ZoltanPartitioner::connectWells(const Comm comm,
const std::vector<Opm::Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const std::unordered_map<int, int>& g2l)
{
auto distributedWells = 0;
@ -301,6 +305,19 @@ void ZoltanPartitioner::connectWells(const Comm comm,
cellIx.push_back(locPos->second);
}
if (possibleFutureConnections) {
const auto possibleFutureConnectionSetIt = possibleFutureConnections->find(well.name());
if (possibleFutureConnectionSetIt != possibleFutureConnections->end()) {
for (auto& global_index : possibleFutureConnectionSetIt->second) {
auto locPos = g2l.find(global_index);
if (locPos == g2l.end()) {
++otherProc;
continue;
}
cellIx.push_back(locPos->second);
}
}
}
if ((otherProc > 0) && !cellIx.empty()) {
++distributedWells;
@ -338,10 +355,11 @@ This is not supported in the current NLDD implementation.)",
template <class GridView, class Element>
std::pair<std::vector<int>, int>
partitionCellsZoltan(const int num_domains,
const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl)
partitionCellsZoltan(const int num_domains,
const GridView& grid_view,
const std::vector<Opm::Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const Opm::ZoltanPartitioningControl<Element>& zoltan_ctrl)
{
if (num_domains <= 1) { // No partitioning => every cell in domain zero.
const auto num_interior_cells =
@ -355,7 +373,7 @@ partitionCellsZoltan(const int num_domains,
}
auto partitioner = ZoltanPartitioner { grid_view, zoltan_ctrl.local_to_global };
partitioner.buildLocalGraph(grid_view, wells, zoltan_ctrl);
partitioner.buildLocalGraph(grid_view, wells, possibleFutureConnections, zoltan_ctrl);
return partitioner.partition(num_domains, grid_view, zoltan_ctrl);
}
@ -574,13 +592,14 @@ std::pair<std::vector<int>, int>
Opm::partitionCells(const std::string& method,
const int num_local_domains,
const GridView& grid_view,
[[maybe_unused]] const std::vector<Well>& wells,
[[maybe_unused]] const ZoltanPartitioningControl<Element>& zoltan_ctrl)
[[maybe_unused]] const std::vector<Well>& wells,
[[maybe_unused]] const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
[[maybe_unused]] const ZoltanPartitioningControl<Element>& zoltan_ctrl)
{
if (method == "zoltan") {
#if HAVE_MPI && HAVE_ZOLTAN
return partitionCellsZoltan(num_local_domains, grid_view, wells, zoltan_ctrl);
return partitionCellsZoltan(num_local_domains, grid_view, wells, possibleFutureConnections, zoltan_ctrl);
#else // !HAVE_MPI || !HAVE_ZOLTAN
@ -665,15 +684,16 @@ Opm::partitionCellsSimple(const int num_cells, const int num_domains)
// Deliberately placed at end of file. No other code beyond this separator.
// ---------------------------------------------------------------------------
#define InstantiatePartitionCells(Grid) \
template std::pair<std::vector<int>, int> \
Opm::partitionCells(const std::string&, \
const int, \
const std::remove_cv_t<std::remove_reference_t< \
decltype(std::declval<Grid>().leafGridView())>>&, \
const std::vector<Opm::Well>&, \
const Opm::ZoltanPartitioningControl< \
typename std::remove_cv_t<std::remove_reference_t< \
#define InstantiatePartitionCells(Grid) \
template std::pair<std::vector<int>, int> \
Opm::partitionCells(const std::string&, \
const int, \
const std::remove_cv_t<std::remove_reference_t< \
decltype(std::declval<Grid>().leafGridView())>>&, \
const std::vector<Opm::Well>&, \
const std::unordered_map<std::string, std::set<int>>*, \
const Opm::ZoltanPartitioningControl< \
typename std::remove_cv_t<std::remove_reference_t< \
decltype(std::declval<Grid>().leafGridView())>>::template Codim<0>::Entity>&)
// ---------------------------------------------------------------------------

View File

@ -21,6 +21,7 @@
#define OPM_ASPINPARTITION_HEADER_INCLUDED
#include <functional>
#include <set>
#include <string>
#include <utility>
#include <vector>
@ -88,11 +89,12 @@ struct ZoltanPartitioningControl
/// on current rank.
template <class GridView, class Element>
std::pair<std::vector<int>, int>
partitionCells(const std::string& method,
const int num_local_domains,
const GridView& grid_view,
const std::vector<Well>& wells,
const ZoltanPartitioningControl<Element>& zoltan_ctrl);
partitionCells(const std::string& method,
const int num_local_domains,
const GridView& grid_view,
const std::vector<Well>& wells,
const std::unordered_map<std::string, std::set<int>>* possibleFutureConnections,
const ZoltanPartitioningControl<Element>& zoltan_ctrl);
/// Read a partitioning from file, assumed to contain one number per cell, its partition number.
/// \return pair containing a partition vector (partition number for each cell), and the number of partitions.

View File

@ -194,11 +194,22 @@ namespace Opm {
// Create cartesian to compressed mapping
const auto& schedule_wells = this->schedule().getWellsatEnd();
const auto& possibleFutureConnections = this->schedule().getPossibleFutureConnections();
// initialize the additional cell connections introduced by wells.
for (const auto& well : schedule_wells)
{
std::vector<int> wellCells = this->getCellsForConnections(well);
// Now add the cells of the possible future connections
const auto possibleFutureConnectionSetIt = possibleFutureConnections.find(well.name());
if (possibleFutureConnectionSetIt != possibleFutureConnections.end()) {
for (auto& global_index : possibleFutureConnectionSetIt->second) {
int compressed_idx = compressedIndexForInterior(global_index);
if (compressed_idx >= 0) { // Ignore connections in inactive/remote cells.
wellCells.push_back(compressed_idx);
}
}
}
for (int cellIdx : wellCells) {
neighbors[cellIdx].insert(wellCells.begin(),
wellCells.end());