LgrOutputHelpers.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA.
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_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
21#define OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
22
23#include <opm/grid/CpGrid.hpp>
25
26#include <opm/input/eclipse/Units/UnitSystem.hpp>
27#include <opm/output/data/Cells.hpp>
28#include <opm/output/data/Solution.hpp>
29#include <opm/output/eclipse/RestartValue.hpp>
30
31#include <algorithm> // for std::min/max
32#include <cstddef> // for std::size_t
33#include <utility> // for std::move
34#include <type_traits> // for std::is_same_v
35#include <vector>
36
37template<class ScalarType>
39 void operator()(int level,
40 int levelIdx,
41 const std::vector<std::vector<ScalarType>>& levelVectors)
42 {
43 min = std::min(min, levelVectors[ level ][ levelIdx ]);
44 }
45
46 ScalarType getValue()
47 {
48 return min;
49 }
50
51 ScalarType min{};
52};
53
54template<class ScalarType>
56 void operator()(int level,
57 int levelIdx,
58 const std::vector<std::vector<ScalarType>>& levelVectors)
59 {
60 max = std::max(max, levelVectors[ level ][ levelIdx ]);
61 }
62
63 ScalarType getValue()
64 {
65 return max;
66 }
67
68 ScalarType max{};
69};
70
71template<class ScalarType>
73 void operator()(int level,
74 int levelIdx,
75 const std::vector<std::vector<ScalarType>>& levelVectors)
76 {
77 partialSum+= levelVectors[ level ][ levelIdx ];
78 ++count;
79 }
80
81 ScalarType getValue()
82 {
83 return partialSum/count;
84 }
85
86 std::size_t count{};
87 ScalarType partialSum{};
88};
89
90namespace Opm
91{
92namespace Lgr
93{
94
107 int level);
108
116template <typename Container>
117Container reorderForOutput(const Container& simulatorContainer,
118 const std::vector<int>& toOutput);
119
126std::vector<std::unordered_map<int,int>>
129
139template <typename ScalarType>
140ScalarType processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
141 const Dune::cpgrid::Entity<0>& element,
142 const Dune::CpGrid& grid);
143
155template <typename ScalarType>
157 int maxLevel,
158 const std::vector<ScalarType>& leafVector,
159 const std::vector<std::vector<int>>& toOutput_refinedLevels,
160 std::vector<std::vector<ScalarType>>& levelVectors);
161
183 const std::vector<std::vector<int>>& toOutput_refinedLevels,
184 const Opm::data::Solution& leafSolution,
185 std::vector<Opm::data::Solution>&);
186
198template <typename Grid>
199void extractRestartValueLevelGrids(const Grid& grid,
200 const Opm::RestartValue& leafRestartValue,
201 std::vector<Opm::RestartValue>& restartValue_levels);
202
203} // namespace Lgr
204} // namespace Opm
205
206template <typename Container>
207Container Opm::Lgr::reorderForOutput(const Container& simulatorContainer,
208 const std::vector<int>& toOutput)
209{
210 // Use toOutput to reorder simulatorContainer
211 Container outputContainer;
212 outputContainer.resize(toOutput.size());
213 for (std::size_t i = 0; i < toOutput.size(); ++i) {
214 outputContainer[i] = simulatorContainer[toOutput[i]];
215 }
216 return outputContainer;
217}
218
219template <typename ScalarType>
220ScalarType Opm::Lgr::processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
221 const Dune::cpgrid::Entity<0>& element,
222 const Dune::CpGrid& grid)
223{
224 const auto& [level, level_indices] = grid.currentData()[element.level()]->getChildrenLevelAndIndexList(element.index());
225
226 // For now, we compute field output properties as follows:
227 // - int fields: use the maximum value among the children data
228 // - double fields: use the average of the children data
229 //
230 // In the future, this could be extended to use the std::string field name
231 // to apply field-specific computation logic.
232 auto childrenDataFunc = [](){
233 if constexpr (std::is_same_v<ScalarType, double>)
235 else
237 }();
238
239 for (const auto& levelIdx : level_indices) {
240 childrenDataFunc(level, levelIdx, levelVectors);
241 }
242 return childrenDataFunc.getValue();
243}
244
245template <typename ScalarType>
247 int maxLevel,
248 const std::vector<ScalarType>& leafVector,
249 const std::vector<std::vector<int>>& toOutput_refinedLevels,
250 std::vector<std::vector<ScalarType>>& levelVectors)
251{
252 for (int level = 0; level <= maxLevel; ++level) {
253 levelVectors[level].resize(grid.levelGridView(level).size(0));
254 }
255 // For level cells that appear in the leaf, extract the data value from leafVector
256 // and assign it to the equivalent level cell.
257 for (const auto& element : Dune::elements(grid.leafGridView())) {
258 levelVectors[element.level()][element.getLevelElem().index()] = leafVector[element.index()];
259 }
260 // Note that all cells from maxLevel have assigned values at this point.
261 // Now, assign values for parent cells (for now, average of children values).
262 if (maxLevel) {
263 for (int level = maxLevel-1; level >= 0; --level) {
264 for (const auto& element : Dune::elements(grid.levelGridView(level))) {
265 if (!element.isLeaf()) {
266 levelVectors[level][element.index()] = Opm::Lgr::processChildrenData(levelVectors,
267 element,
268 grid);
269 }
270 }
271 }
272 }
273 // Use toOutput_levels to reorder in ascending level cartesian indices
274 for (int level = 1; level<=maxLevel; ++level) { // exclude level zero (does not need reordering)
275 levelVectors[level] = Opm::Lgr::reorderForOutput(levelVectors[level], toOutput_refinedLevels[level-1]);
276 }
277}
278
279template <typename Grid>
281 const Opm::RestartValue& leafRestartValue,
282 std::vector<Opm::RestartValue>& restartValue_levels)
283{
284 if constexpr (std::is_same_v<Grid, Dune::CpGrid>) {
285
286 int maxLevel = grid.maxLevel();
287 restartValue_levels.resize(maxLevel+1); // level 0, 1, ..., max level
288
289 // To store leafRestartValue.extra data in the order expected
290 // by outout files (increasing level Cartesian indices)
291 std::vector<std::vector<int>> toOutput_refinedLevels{};
292 toOutput_refinedLevels.resize(maxLevel); // exclude level zero (does not need reordering)
293
294 const Opm::LevelCartesianIndexMapper<Dune::CpGrid> levelCartMapp(grid);
295 for (int level = 1; level <= maxLevel; ++level) { // exclude level zero (does not need reordering)
296 toOutput_refinedLevels[level-1] = mapLevelIndicesToCartesianOutputOrder(grid, levelCartMapp, level);
297 }
298
299 std::vector<Opm::data::Solution> dataSolutionLevels{};
301 toOutput_refinedLevels,
302 leafRestartValue.solution,
303 dataSolutionLevels);
304
305
306
307 for (int level = 0; level <= maxLevel; ++level) {
308 restartValue_levels[level] = Opm::RestartValue(std::move(dataSolutionLevels[level]),
309 leafRestartValue.wells,
310 leafRestartValue.grp_nwrk,
311 leafRestartValue.aquifer,
312 level);
313 }
314
315 for (const auto& [rst_key, leafVector] : leafRestartValue.extra) {
316
317 std::vector<std::vector<double>> levelVectors{};
318 levelVectors.resize(maxLevel+1);
319
320 if (rst_key.key == "OPMEXTRA") {
321 // For OPMEXTRA, leafVector has size 1 instead of
322 // grid.leafGridView().size(0)
323 continue; // skip it
324 }
325 Opm::Lgr::populateDataVectorLevelGrids<double>(grid,
326 maxLevel,
327 leafVector,
328 toOutput_refinedLevels,
329 levelVectors);
330 for (int level = 0; level <= maxLevel; ++level) {
331 restartValue_levels[level].addExtra(rst_key.key, rst_key.dim, std::move(levelVectors[level]));
332 }
333 }
334 }
335}
336
337#endif // OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
[ provides Dune::Grid ]
Definition: CpGrid.hpp:203
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & currentData() const
Returns either data_ or distributed_data_(if non empty).
int size(int level, int codim) const
Number of grid entities per level and codim.
int level() const
Return the level of the entity in the grid hierarchy. Level = 0 represents the coarsest grid.
Definition: Entity.hpp:442
int index() const
The (positive) index of an entity. Not a Dune interface method.
Definition: EntityRep.hpp:126
Definition: cpgrid/LevelCartesianIndexMapper.hpp:53
std::vector< std::unordered_map< int, int > > levelCartesianToLevelCompressedMaps(const Dune::CpGrid &grid, const Opm::LevelCartesianIndexMapper< Dune::CpGrid > &levelCartMapp)
Map level Cartesian index to level compressed index (active cell)
void extractSolutionLevelGrids(const Dune::CpGrid &grid, const std::vector< std::vector< int > > &toOutput_refinedLevels, const Opm::data::Solution &leafSolution, std::vector< Opm::data::Solution > &)
Extracts and organizes solution data for all grid refinement levels.
void populateDataVectorLevelGrids(const Dune::CpGrid &grid, int maxLevel, const std::vector< ScalarType > &leafVector, const std::vector< std::vector< int > > &toOutput_refinedLevels, std::vector< std::vector< ScalarType > > &levelVectors)
Populate level data vectors based on leaf vector, for a specific named data field.
Definition: LgrOutputHelpers.hpp:246
void extractRestartValueLevelGrids(const Grid &grid, const Opm::RestartValue &leafRestartValue, std::vector< Opm::RestartValue > &restartValue_levels)
Constructs restart-value containers for all grid refinement levels.
Definition: LgrOutputHelpers.hpp:280
ScalarType processChildrenData(const std::vector< std::vector< ScalarType > > &levelVectors, const Dune::cpgrid::Entity< 0 > &element, const Dune::CpGrid &grid)
Definition: LgrOutputHelpers.hpp:220
Container reorderForOutput(const Container &simulatorContainer, const std::vector< int > &toOutput)
Reorder data from a simulation container into the order assumed by output for refined level grids.
Definition: LgrOutputHelpers.hpp:207
std::vector< int > mapLevelIndicesToCartesianOutputOrder(const Dune::CpGrid &grid, const Opm::LevelCartesianIndexMapper< Dune::CpGrid > &levelCartMapp, int level)
Builds a mapping from level element indices to the Cartesian ordering required by output files.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
Definition: LgrOutputHelpers.hpp:72
ScalarType partialSum
Definition: LgrOutputHelpers.hpp:87
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:81
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:73
std::size_t count
Definition: LgrOutputHelpers.hpp:86
Definition: LgrOutputHelpers.hpp:55
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:63
ScalarType max
Definition: LgrOutputHelpers.hpp:68
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:56
Definition: LgrOutputHelpers.hpp:38
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:39
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:46
ScalarType min
Definition: LgrOutputHelpers.hpp:51