NlddReporting.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025, SINTEF Digital
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_NLDD_REPORTING_HEADER_INCLUDED
21#define OPM_NLDD_REPORTING_HEADER_INCLUDED
22
23#include <dune/grid/common/gridenums.hh>
24#include <dune/grid/common/partitionset.hh>
25
28
29#include <algorithm>
30#include <cmath>
31#include <filesystem>
32#include <fstream>
33#include <limits>
34#include <numeric>
35#include <string_view>
36#include <vector>
37
38#include <fmt/format.h>
39
40namespace Opm {
41
42namespace details {
43
48{
49 std::vector<int> all_owned;
50 std::vector<int> all_overlap;
51 std::vector<int> all_wells;
52 std::vector<int> all_domains;
53
54 explicit DomainInfo(std::size_t size)
55 : all_owned(size)
56 , all_overlap(size)
57 , all_wells(size)
58 , all_domains(size)
59 {}
60};
61
67
74void writeNlddFile(const std::string& file,
75 std::string_view header,
76 const std::vector<int>& data);
77
78} // namespace details
79
88void reportNlddStatistics(const std::vector<SimulatorReport>& domain_reports,
89 const SimulatorReport& local_report,
90 const bool output_cout,
91 const Parallel::Communication& comm);
92
103template <class Grid, class Domain, class ElementMapper, class CartMapper>
105 const std::filesystem::path& odir,
106 const std::vector<Domain>& domains,
107 const std::vector<SimulatorReport>& domain_reports,
108 const Grid& grid,
109 const ElementMapper& elementMapper,
110 const CartMapper& cartMapper)
111{
112 const auto& dims = cartMapper.cartesianDimensions();
113 const auto total_size = dims[0] * dims[1] * dims[2];
114 const auto& comm = grid.comm();
115 const int rank = comm.rank();
116
117 // Create a cell-to-iterations mapping for this process
118 std::vector<int> cell_iterations(grid.size(0), 0);
119
120 // Populate the mapping with iteration counts for each domain
121 const auto ds = domains.size();
122 for (auto domain_idx = 0*ds; domain_idx < ds; ++domain_idx) {
123 const auto& domain = domains[domain_idx];
124 const auto& report = domain_reports[domain_idx];
125 const int iterations = report.success.total_newton_iterations +
126 report.failure.total_newton_iterations;
127
128 for (const int cell_idx : domain.cells) {
129 cell_iterations[cell_idx] = iterations;
130 }
131 }
132
133 // Create a full-sized vector initialized with zeros (indicating inactive cells)
134 auto full_iterations = std::vector<int>(total_size, 0);
135
136 // Convert local cell indices to cartesian indices
137 const auto& gridView = grid.leafGridView();
138 for (const auto& cell : elements(gridView, Dune::Partitions::interior)) {
139 const int cell_idx = elementMapper.index(cell);
140 const int cart_idx = cartMapper.cartesianIndex(cell_idx);
141 full_iterations[cart_idx] = cell_iterations[cell_idx];
142 }
143
144 // Gather all iteration data using max operation
145 comm.max(full_iterations.data(), full_iterations.size());
146
147 // Only rank 0 writes the file
148 if (rank == 0) {
149 details::writeNlddFile(odir / "ResInsight_nonlinear_iterations.txt",
150 "NLDD_ITER", full_iterations);
151 }
152}
153
163template <class Grid, class Domain, class ElementMapper, class CartMapper>
165 const std::filesystem::path& odir,
166 const std::vector<Domain>& domains,
167 const Grid& grid,
168 const ElementMapper& elementMapper,
169 const CartMapper& cartMapper)
170{
171 const auto& dims = cartMapper.cartesianDimensions();
172 const auto total_size = dims[0] * dims[1] * dims[2];
173 const auto& comm = grid.comm();
174 const int rank = comm.rank();
175
176 const auto& partition_vector = reconstitutePartitionVector(domains, grid);
177
178 // Create a full-sized vector initialized with -1 (indicating inactive cells)
179 auto full_partition = std::vector<int>(total_size, -1);
180
181 // Fill active cell values for this rank
182 auto i = 0;
183 for (const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
184 full_partition[cartMapper.cartesianIndex(elementMapper.index(cell))] = partition_vector[i++];
185 }
186
187 // Gather all partitions using max operation
188 comm.max(full_partition.data(), full_partition.size());
189
190 // Only rank 0 writes the file
191 if (rank == 0) {
192 details::writeNlddFile(odir / "ResInsight_compatible_partition.txt",
193 "NLDD_DOM", full_partition);
194 }
195
196 const auto nDigit = 1 + static_cast<int>(std::floor(std::log10(comm.size())));
197 auto partition_fname = odir / fmt::format("{1:0>{0}}", nDigit, rank);
198 std::ofstream pfile { partition_fname };
199
200 auto cell_index = 0;
201 for (const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
202 pfile << rank << ' '
203 << cartMapper.cartesianIndex(elementMapper.index(cell)) << ' '
204 << partition_vector[cell_index++] << '\n';
205 }
206}
207
218template <class Grid, class Domain>
220 const std::vector<int>& partition_vector,
221 const std::vector<Domain>& domains,
222 SimulatorReport& local_reports_accumulated,
223 std::vector<SimulatorReport>& domain_reports_accumulated,
224 const Grid& grid,
225 int num_wells)
226{
227 const auto& gridView = grid.leafGridView();
228 const auto& comm = grid.comm();
229 const int rank = comm.rank();
230
231 const int num_domains = domains.size();
232 const int owned_cells = partition_vector.size();
233
234 // Count overlap cells using grid view iteration
235 int overlap_cells = std::count_if(elements(gridView).begin(), elements(gridView).end(),
236 [](const auto& cell) { return cell.partitionType() == Dune::OverlapEntity; });
237
238 // Store data for summary output
239 local_reports_accumulated.success.num_wells = num_wells;
240 local_reports_accumulated.success.num_domains = num_domains;
241 local_reports_accumulated.success.num_overlap_cells = overlap_cells;
242 local_reports_accumulated.success.num_owned_cells = owned_cells;
243
244 // Set statistics for each domain report
245 const auto dr_size = domain_reports_accumulated.size();
246 for (auto i = 0*dr_size; i < dr_size; ++i) {
247 auto& domain_report = domain_reports_accumulated[i];
248 domain_report.success.num_domains = 1;
249 domain_report.success.num_overlap_cells = 0;
250 domain_report.success.num_owned_cells = domains[i].cells.size();
251 }
252
253 // Gather data from all ranks
254 details::DomainInfo info(comm.size());
255
256 comm.gather(&owned_cells, info.all_owned.data(), 1, 0);
257 comm.gather(&overlap_cells, info.all_overlap.data(), 1, 0);
258 comm.gather(&num_wells, info.all_wells.data(), 1, 0);
259 comm.gather(&num_domains, info.all_domains.data(), 1, 0);
260
261 if (rank == 0) {
263 }
264}
265
274template <class Grid, class Domain>
276 const std::vector<Domain>& domains,
277 const Grid& grid)
278{
279 const auto& comm = grid.comm();
280 const int rank = comm.rank();
281
282 auto numD = std::vector<int>(comm.size() + 1, 0);
283 numD[rank + 1] = static_cast<int>(domains.size());
284 comm.sum(numD.data(), numD.size());
285 std::partial_sum(numD.begin(), numD.end(), numD.begin());
286
287 auto p = std::vector<int>(grid.size(0));
288 auto maxCellIdx = std::numeric_limits<int>::min();
289
290 auto d = numD[rank];
291 for (const auto& domain : domains) {
292 for (const auto& cell : domain.cells) {
293 p[cell] = d;
294 if (cell > maxCellIdx) {
295 maxCellIdx = cell;
296 }
297 }
298
299 ++d;
300 }
301
302 p.erase(p.begin() + maxCellIdx + 1, p.end());
303 return p;
304}
305
306} // namespace Opm
307
308#endif // OPM_NLDD_REPORTING_HEADER_INCLUDED
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
void printDistributionSummary(const DomainInfo &info)
void writeNlddFile(const std::string &file, std::string_view header, const std::vector< int > &data)
Definition: blackoilboundaryratevector.hh:39
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)
Definition: NlddReporting.hpp:219
void reportNlddStatistics(const std::vector< SimulatorReport > &domain_reports, const SimulatorReport &local_report, const bool output_cout, const Parallel::Communication &comm)
std::vector< int > reconstitutePartitionVector(const std::vector< Domain > &domains, const Grid &grid)
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)
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)
Definition: NlddReporting.hpp:104
Definition: SimulatorReport.hpp:122
SimulatorReportSingle success
Definition: SimulatorReport.hpp:123
int num_wells
Definition: SimulatorReport.hpp:65
int num_owned_cells
Definition: SimulatorReport.hpp:67
int num_overlap_cells
Definition: SimulatorReport.hpp:66
int num_domains
Definition: SimulatorReport.hpp:64
Definition: NlddReporting.hpp:48
std::vector< int > all_owned
All owned cells.
Definition: NlddReporting.hpp:49
DomainInfo(std::size_t size)
Definition: NlddReporting.hpp:54
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