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
129template <typename ScalarType>
130ScalarType processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
131 const Dune::cpgrid::Entity<0>& element,
132 const Dune::CpGrid& grid);
133
145template <typename ScalarType>
147 int maxLevel,
148 const std::vector<ScalarType>& leafVector,
149 const std::vector<std::vector<int>>& toOutput_refinedLevels,
150 std::vector<std::vector<ScalarType>>& levelVectors);
151
173 const std::vector<std::vector<int>>& toOutput_refinedLevels,
174 const Opm::data::Solution& leafSolution,
175 std::vector<Opm::data::Solution>&);
176
188template <typename Grid>
189void extractRestartValueLevelGrids(const Grid& grid,
190 const Opm::RestartValue& leafRestartValue,
191 std::vector<Opm::RestartValue>& restartValue_levels);
192
193} // namespace Lgr
194} // namespace Opm
195
196template <typename Container>
197Container Opm::Lgr::reorderForOutput(const Container& simulatorContainer,
198 const std::vector<int>& toOutput)
199{
200 // Use toOutput to reorder simulatorContainer
201 Container outputContainer;
202 outputContainer.resize(toOutput.size());
203 for (std::size_t i = 0; i < toOutput.size(); ++i) {
204 outputContainer[i] = simulatorContainer[toOutput[i]];
205 }
206 return outputContainer;
207}
208
209template <typename ScalarType>
210ScalarType Opm::Lgr::processChildrenData(const std::vector<std::vector<ScalarType>>& levelVectors,
211 const Dune::cpgrid::Entity<0>& element,
212 const Dune::CpGrid& grid)
213{
214 const auto& [level, level_indices] = grid.currentData()[element.level()]->getChildrenLevelAndIndexList(element.index());
215
216 // For now, we compute field output properties as follows:
217 // - int fields: use the maximum value among the children data
218 // - double fields: use the average of the children data
219 //
220 // In the future, this could be extended to use the std::string field name
221 // to apply field-specific computation logic.
222 auto childrenDataFunc = [](){
223 if constexpr (std::is_same_v<ScalarType, double>)
225 else
227 }();
228
229 for (const auto& levelIdx : level_indices) {
230 childrenDataFunc(level, levelIdx, levelVectors);
231 }
232 return childrenDataFunc.getValue();
233}
234
235template <typename ScalarType>
237 int maxLevel,
238 const std::vector<ScalarType>& leafVector,
239 const std::vector<std::vector<int>>& toOutput_refinedLevels,
240 std::vector<std::vector<ScalarType>>& levelVectors)
241{
242 for (int level = 0; level <= maxLevel; ++level) {
243 levelVectors[level].resize(grid.levelGridView(level).size(0));
244 }
245 // For level cells that appear in the leaf, extract the data value from leafVector
246 // and assign it to the equivalent level cell.
247 for (const auto& element : Dune::elements(grid.leafGridView())) {
248 levelVectors[element.level()][element.getLevelElem().index()] = leafVector[element.index()];
249 }
250 // Note that all cells from maxLevel have assigned values at this point.
251 // Now, assign values for parent cells (for now, average of children values).
252 if (maxLevel) {
253 for (int level = maxLevel-1; level >= 0; --level) {
254 for (const auto& element : Dune::elements(grid.levelGridView(level))) {
255 if (!element.isLeaf()) {
256 levelVectors[level][element.index()] = Opm::Lgr::processChildrenData(levelVectors,
257 element,
258 grid);
259 }
260 }
261 }
262 }
263 // Use toOutput_levels to reorder in ascending level cartesian indices
264 for (int level = 1; level<=maxLevel; ++level) { // exclude level zero (does not need reordering)
265 levelVectors[level] = Opm::Lgr::reorderForOutput(levelVectors[level], toOutput_refinedLevels[level-1]);
266 }
267}
268
269template <typename Grid>
271 const Opm::RestartValue& leafRestartValue,
272 std::vector<Opm::RestartValue>& restartValue_levels)
273{
274 if constexpr (std::is_same_v<Grid, Dune::CpGrid>) {
275
276 int maxLevel = grid.maxLevel();
277 restartValue_levels.resize(maxLevel+1); // level 0, 1, ..., max level
278
279 // To store leafRestartValue.extra data in the order expected
280 // by outout files (increasing level Cartesian indices)
281 std::vector<std::vector<int>> toOutput_refinedLevels{};
282 toOutput_refinedLevels.resize(maxLevel); // exclude level zero (does not need reordering)
283
284 const Opm::LevelCartesianIndexMapper<Dune::CpGrid> levelCartMapp(grid);
285 for (int level = 1; level <= maxLevel; ++level) { // exclude level zero (does not need reordering)
286 toOutput_refinedLevels[level-1] = mapLevelIndicesToCartesianOutputOrder(grid, levelCartMapp, level);
287 }
288
289 std::vector<Opm::data::Solution> dataSolutionLevels{};
291 toOutput_refinedLevels,
292 leafRestartValue.solution,
293 dataSolutionLevels);
294
295
296
297 for (int level = 0; level <= maxLevel; ++level) {
298 restartValue_levels[level] = Opm::RestartValue(std::move(dataSolutionLevels[level]),
299 leafRestartValue.wells,
300 leafRestartValue.grp_nwrk,
301 leafRestartValue.aquifer,
302 level);
303 }
304
305 for (const auto& [rst_key, leafVector] : leafRestartValue.extra) {
306
307 std::vector<std::vector<double>> levelVectors{};
308 levelVectors.resize(maxLevel+1);
309
310 Opm::Lgr::populateDataVectorLevelGrids<double>(grid,
311 maxLevel,
312 leafVector,
313 toOutput_refinedLevels,
314 levelVectors);
315
316 for (int level = 0; level <= maxLevel; ++level) {
317 restartValue_levels[level].addExtra(rst_key.key, rst_key.dim, std::move(levelVectors[level]));
318 }
319 }
320 }
321}
322
323#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
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:236
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:270
ScalarType processChildrenData(const std::vector< std::vector< ScalarType > > &levelVectors, const Dune::cpgrid::Entity< 0 > &element, const Dune::CpGrid &grid)
Definition: LgrOutputHelpers.hpp:210
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:197
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