20#ifndef OPM_NLDD_REPORTING_HEADER_INCLUDED
21#define OPM_NLDD_REPORTING_HEADER_INCLUDED
23#include <dune/grid/common/gridenums.hh>
24#include <dune/grid/common/partitionset.hh>
26#include <opm/simulators/timestepping/SimulatorReport.hpp>
27#include <opm/simulators/utils/ParallelCommunication.hpp>
38#include <fmt/format.h>
66void printDistributionSummary(
const DomainInfo& info);
74void writeNlddFile(
const std::string&
file,
76 const std::vector<int>& data);
91 const Parallel::Communication& comm);
103template <
class Gr
id,
class Domain,
class ElementMapper,
class CartMapper>
105 const std::filesystem::path&
odir,
106 const std::vector<Domain>&
domains,
109 const ElementMapper& elementMapper,
114 const auto& comm = grid.comm();
115 const int rank = comm.rank();
125 const int iterations = report.success.total_newton_iterations +
126 report.failure.total_newton_iterations;
137 const auto& gridView = grid.leafGridView();
138 for (
const auto& cell :
elements(gridView, Dune::Partitions::interior)) {
139 const int cell_idx = elementMapper.index(cell);
149 details::writeNlddFile(
odir /
"ResInsight_nonlinear_iterations.txt",
163template <
class Gr
id,
class Domain,
class ElementMapper,
class CartMapper>
165 const std::filesystem::path&
odir,
166 const std::vector<Domain>&
domains,
168 const ElementMapper& elementMapper,
173 const auto& comm = grid.comm();
174 const int rank = comm.rank();
183 for (
const auto& cell :
elements(grid.leafGridView(), Dune::Partitions::interior)) {
192 details::writeNlddFile(
odir /
"ResInsight_compatible_partition.txt",
196 const auto nDigit = 1 +
static_cast<int>(std::floor(std::log10(comm.size())));
201 for (
const auto& cell :
elements(grid.leafGridView(), Dune::Partitions::interior)) {
203 <<
cartMapper.cartesianIndex(elementMapper.index(cell)) <<
' '
218template <
class Gr
id,
class Domain>
221 const std::vector<Domain>&
domains,
227 const auto& gridView = grid.leafGridView();
228 const auto& comm = grid.comm();
229 const int rank = comm.rank();
231 const int num_domains =
domains.size();
236 [](
const auto& cell) {
return cell.partitionType() == Dune::OverlapEntity; });
258 comm.gather(&num_wells, info.
all_wells.data(), 1, 0);
259 comm.gather(&num_domains, info.
all_domains.data(), 1, 0);
262 details::printDistributionSummary(info);
274template <
class Gr
id,
class Domain>
276 const std::vector<Domain>&
domains,
279 const auto& comm = grid.comm();
280 const int rank = comm.rank();
282 auto numD = std::vector<int>(comm.size() + 1, 0);
285 std::partial_sum(
numD.begin(),
numD.end(),
numD.begin());
287 auto p = std::vector<int>(grid.size(0));
288 auto maxCellIdx = std::numeric_limits<int>::min();
292 for (
const auto& cell :
domain.cells) {
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
void printDomainDistributionSummary(const std::vector< int > &partition_vector, const std::vector< Domain > &domains, SimulatorReport &local_reports_accumulated, std::vector< SimulatorReport > &domain_reports_accumulated, const Grid &grid, int num_wells)
Prints a summary of domain distribution across ranks.
Definition NlddReporting.hpp:219
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void reportNlddStatistics(const std::vector< SimulatorReport > &domain_reports, const SimulatorReport &local_report, const bool output_cout, const Parallel::Communication &comm)
Reports NLDD statistics after simulation.
Definition NlddReporting.cpp:97
std::vector< int > reconstitutePartitionVector(const std::vector< Domain > &domains, const Grid &grid)
Reconstructs the partition vector that maps each grid cell to its corresponding domain ID,...
Definition NlddReporting.hpp:275
void writePartitions(const std::filesystem::path &odir, const std::vector< Domain > &domains, const Grid &grid, const ElementMapper &elementMapper, const CartMapper &cartMapper)
Writes the partition vector to a file in ResInsight compatible format and a partition file for each r...
Definition NlddReporting.hpp:164
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir, const std::vector< Domain > &domains, const std::vector< SimulatorReport > &domain_reports, const Grid &grid, const ElementMapper &elementMapper, const CartMapper &cartMapper)
Writes the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition NlddReporting.hpp:104
Definition SimulatorReport.hpp:122
Struct holding information about domains used for printing a summary.
Definition NlddReporting.hpp:48
std::vector< int > all_owned
All owned cells.
Definition NlddReporting.hpp:49
std::vector< int > all_domains
All domains.
Definition NlddReporting.hpp:52
std::vector< int > all_wells
All wells.
Definition NlddReporting.hpp:51
std::vector< int > all_overlap
All overlap cells.
Definition NlddReporting.hpp:50