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#if HAVE_OPM_COMMON
27#include <opm/input/eclipse/Units/UnitSystem.hpp>
28#include <opm/output/data/Cells.hpp>
29#include <opm/output/data/Solution.hpp>
30#include <opm/output/eclipse/RestartValue.hpp>
31#endif
32
33#include <algorithm> // for std::min/max
34#include <cstddef> // for std::size_t
35#include <utility> // for std::move
36#include <type_traits> // for std::is_same_v
37#include <vector>
38
39template<class ScalarType>
41 void operator()(int level,
42 int levelIdx,
43 const std::vector<std::vector<ScalarType>>& levelVectors)
44 {
45 min = std::min(min, levelVectors[ level ][ levelIdx ]);
46 }
47
48 ScalarType getValue()
49 {
50 return min;
51 }
52
53 ScalarType min{};
54};
55
56template<class ScalarType>
58 void operator()(int level,
59 int levelIdx,
60 const std::vector<std::vector<ScalarType>>& levelVectors)
61 {
62 max = std::max(max, levelVectors[ level ][ levelIdx ]);
63 }
64
65 ScalarType getValue()
66 {
67 return max;
68 }
69
70 ScalarType max{};
71};
72
73template<class ScalarType>
75 void operator()(int level,
76 int levelIdx,
77 const std::vector<std::vector<ScalarType>>& levelVectors)
78 {
79 partialSum+= levelVectors[ level ][ levelIdx ];
80 ++count;
81 }
82
83 ScalarType getValue()
84 {
85 return partialSum/count;
86 }
87
88 std::size_t count{};
89 ScalarType partialSum{};
90};
91
92namespace Opm
93{
94namespace Lgr
95{
96
109 int level);
110
118template <typename Container>
119Container reorderForOutput(const Container& simulatorContainer,
120 const std::vector<int>& toOutput);
121
128std::vector<std::unordered_map<int,int>>
131
141template <typename ScalarType>
142ScalarType processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
143 const Dune::cpgrid::Entity<0>& element,
144 const Dune::CpGrid& grid);
145
157template <typename ScalarType>
159 int maxLevel,
160 const std::vector<ScalarType>& leafVector,
161 const std::vector<std::vector<int>>& toOutput_refinedLevels,
162 std::vector<std::vector<ScalarType>>& levelVectors);
163
184#if HAVE_OPM_COMMON
185void extractSolutionLevelGrids(const Dune::CpGrid& grid,
186 const std::vector<std::vector<int>>& toOutput_refinedLevels,
187 const Opm::data::Solution& leafSolution,
188 std::vector<Opm::data::Solution>&);
189
201template <typename Grid>
202void extractRestartValueLevelGrids(const Grid& grid,
203 const Opm::RestartValue& leafRestartValue,
204 std::vector<Opm::RestartValue>& restartValue_levels);
205#endif
206
207} // namespace Lgr
208} // namespace Opm
209
210template <typename Container>
211Container Opm::Lgr::reorderForOutput(const Container& simulatorContainer,
212 const std::vector<int>& toOutput)
213{
214 // Use toOutput to reorder simulatorContainer
215 Container outputContainer;
216 outputContainer.resize(toOutput.size());
217 for (std::size_t i = 0; i < toOutput.size(); ++i) {
218 outputContainer[i] = simulatorContainer[toOutput[i]];
219 }
220 return outputContainer;
221}
222
223template <typename ScalarType>
224ScalarType Opm::Lgr::processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
225 const Dune::cpgrid::Entity<0>& element,
226 const Dune::CpGrid& grid)
227{
228 const auto& [level, level_indices] = grid.currentData()[element.level()]->getChildrenLevelAndIndexList(element.index());
229
230 // For now, we compute field output properties as follows:
231 // - int fields: use the maximum value among the children data
232 // - double fields: use the average of the children data
233 //
234 // In the future, this could be extended to use the std::string field name
235 // to apply field-specific computation logic.
236 auto childrenDataFunc = [](){
237 if constexpr (std::is_same_v<ScalarType, double>)
239 else
241 }();
242
243 for (const auto& levelIdx : level_indices) {
244 childrenDataFunc(level, levelIdx, levelVectors);
245 }
246 return childrenDataFunc.getValue();
247}
248
249template <typename ScalarType>
251 int maxLevel,
252 const std::vector<ScalarType>& leafVector,
253 const std::vector<std::vector<int>>& toOutput_refinedLevels,
254 std::vector<std::vector<ScalarType>>& levelVectors)
255{
256 for (int level = 0; level <= maxLevel; ++level) {
257 levelVectors[level].resize(grid.levelGridView(level).size(0));
258 }
259 // For level cells that appear in the leaf, extract the data value from leafVector
260 // and assign it to the equivalent level cell.
261 for (const auto& element : Dune::elements(grid.leafGridView())) {
262 levelVectors[element.level()][element.getLevelElem().index()] = leafVector[element.index()];
263 }
264 // Note that all cells from maxLevel have assigned values at this point.
265 // Now, assign values for parent cells (for now, average of children values).
266 if (maxLevel) {
267 for (int level = maxLevel-1; level >= 0; --level) {
268 for (const auto& element : Dune::elements(grid.levelGridView(level))) {
269 if (!element.isLeaf()) {
270 levelVectors[level][element.index()] = Opm::Lgr::processChildrenData(levelVectors,
271 element,
272 grid);
273 }
274 }
275 }
276 }
277 // Use toOutput_levels to reorder in ascending level cartesian indices
278 for (int level = 1; level<=maxLevel; ++level) { // exclude level zero (does not need reordering)
279 levelVectors[level] = Opm::Lgr::reorderForOutput(levelVectors[level], toOutput_refinedLevels[level-1]);
280 }
281}
282
283#if HAVE_OPM_COMMON
284template <typename Grid>
285void Opm::Lgr::extractRestartValueLevelGrids(const Grid& grid,
286 const Opm::RestartValue& leafRestartValue,
287 std::vector<Opm::RestartValue>& restartValue_levels)
288{
289 if constexpr (std::is_same_v<Grid, Dune::CpGrid>) {
290
291 int maxLevel = grid.maxLevel();
292 restartValue_levels.resize(maxLevel+1); // level 0, 1, ..., max level
293
294 // To store leafRestartValue.extra data in the order expected
295 // by outout files (increasing level Cartesian indices)
296 std::vector<std::vector<int>> toOutput_refinedLevels{};
297 toOutput_refinedLevels.resize(maxLevel); // exclude level zero (does not need reordering)
298
299 const Opm::LevelCartesianIndexMapper<Dune::CpGrid> levelCartMapp(grid);
300 for (int level = 1; level <= maxLevel; ++level) { // exclude level zero (does not need reordering)
301 toOutput_refinedLevels[level-1] = mapLevelIndicesToCartesianOutputOrder(grid, levelCartMapp, level);
302 }
303
304 std::vector<Opm::data::Solution> dataSolutionLevels{};
305 extractSolutionLevelGrids(grid,
306 toOutput_refinedLevels,
307 leafRestartValue.solution,
308 dataSolutionLevels);
309
310
311
312 for (int level = 0; level <= maxLevel; ++level) {
313 restartValue_levels[level] = Opm::RestartValue(std::move(dataSolutionLevels[level]),
314 leafRestartValue.wells,
315 leafRestartValue.grp_nwrk,
316 leafRestartValue.aquifer,
317 level);
318 }
319
320 for (const auto& [rst_key, leafVector] : leafRestartValue.extra) {
321
322 std::vector<std::vector<double>> levelVectors{};
323 levelVectors.resize(maxLevel+1);
324
325 if (rst_key.key == "OPMEXTRA") {
326 // For OPMEXTRA, leafVector has size 1 instead of
327 // grid.leafGridView().size(0)
328 continue; // skip it
329 }
330 Opm::Lgr::populateDataVectorLevelGrids<double>(grid,
331 maxLevel,
332 leafVector,
333 toOutput_refinedLevels,
334 levelVectors);
335 for (int level = 0; level <= maxLevel; ++level) {
336 restartValue_levels[level].addExtra(rst_key.key, rst_key.dim, std::move(levelVectors[level]));
337 }
338 }
339 }
340}
341#endif
342
343#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:473
int index() const
The (positive) index of an entity. Not a Dune interface method.
Definition: EntityRep.hpp:125
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 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:250
ScalarType processChildrenData(const std::vector< std::vector< ScalarType > > &levelVectors, const Dune::cpgrid::Entity< 0 > &element, const Dune::CpGrid &grid)
Definition: LgrOutputHelpers.hpp:224
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:211
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:74
ScalarType partialSum
Definition: LgrOutputHelpers.hpp:89
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:83
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:75
std::size_t count
Definition: LgrOutputHelpers.hpp:88
Definition: LgrOutputHelpers.hpp:57
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:65
ScalarType max
Definition: LgrOutputHelpers.hpp:70
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:58
Definition: LgrOutputHelpers.hpp:40
void operator()(int level, int levelIdx, const std::vector< std::vector< ScalarType > > &levelVectors)
Definition: LgrOutputHelpers.hpp:41
ScalarType getValue()
Definition: LgrOutputHelpers.hpp:48
ScalarType min
Definition: LgrOutputHelpers.hpp:53