20#ifndef OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
21#define OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
23#include <opm/grid/CpGrid.hpp>
24#include <opm/grid/cpgrid/LevelCartesianIndexMapper.hpp>
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>
39template<
class ScalarType>
41 void operator()(
int level,
43 const std::vector<std::vector<ScalarType>>& levelVectors)
45 min = std::min(min, levelVectors[ level ][ levelIdx ]);
56template<
class ScalarType>
58 void operator()(
int level,
60 const std::vector<std::vector<ScalarType>>& levelVectors)
62 max = std::max(max, levelVectors[ level ][ levelIdx ]);
73template<
class ScalarType>
75 void operator()(
int level,
77 const std::vector<std::vector<ScalarType>>& levelVectors)
79 partialSum+= levelVectors[ level ][ levelIdx ];
85 return partialSum/count;
89 ScalarType partialSum{};
107std::vector<int> mapLevelIndicesToCartesianOutputOrder(
const Dune::CpGrid& grid,
118template <
typename Container>
119Container reorderForOutput(
const Container& simulatorContainer,
120 const std::vector<int>& toOutput);
128std::vector<std::unordered_map<int,int>>
129levelCartesianToLevelCompressedMaps(
const Dune::CpGrid& grid,
141template <
typename ScalarType>
142ScalarType processChildrenData(
const std::vector<std::vector<ScalarType>>& levelVectors,
157template <
typename ScalarType>
158void populateDataVectorLevelGrids(
const Dune::CpGrid& grid,
160 const std::vector<ScalarType>& leafVector,
161 const std::vector<std::vector<int>>& toOutput_refinedLevels,
162 std::vector<std::vector<ScalarType>>& levelVectors);
186 const std::vector<std::vector<int>>& toOutput_refinedLevels,
187 const Opm::data::Solution& leafSolution,
188 std::vector<Opm::data::Solution>&);
201template <
typename Gr
id>
202void extractRestartValueLevelGrids(
const Grid& grid,
203 const Opm::RestartValue& leafRestartValue,
204 std::vector<Opm::RestartValue>& restartValue_levels);
210template <
typename Container>
211Container Opm::Lgr::reorderForOutput(
const Container& simulatorContainer,
212 const std::vector<int>& toOutput)
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]];
220 return outputContainer;
223template <
typename ScalarType>
224ScalarType Opm::Lgr::processChildrenData(
const std::vector<std::vector<ScalarType>>& levelVectors,
228 const auto& [level, level_indices] = grid.
currentData()[element.
level()]->getChildrenLevelAndIndexList(element.
index());
236 auto childrenDataFunc = [](){
237 if constexpr (std::is_same_v<ScalarType, double>)
243 for (
const auto& levelIdx : level_indices) {
244 childrenDataFunc(level, levelIdx, levelVectors);
246 return childrenDataFunc.getValue();
249template <
typename ScalarType>
250void Opm::Lgr::populateDataVectorLevelGrids(
const Dune::CpGrid& grid,
252 const std::vector<ScalarType>& leafVector,
253 const std::vector<std::vector<int>>& toOutput_refinedLevels,
254 std::vector<std::vector<ScalarType>>& levelVectors)
256 for (
int level = 0; level <= maxLevel; ++level) {
257 levelVectors[level].resize(grid.levelGridView(level).
size(0));
261 for (
const auto& element : Dune::elements(grid.leafGridView())) {
267 for (
int level = maxLevel-1; level >= 0; --level) {
268 for (
const auto& element : Dune::elements(grid.levelGridView(level))) {
270 levelVectors[level][element.
index()] = Opm::Lgr::processChildrenData(levelVectors,
278 for (
int level = 1; level<=maxLevel; ++level) {
279 levelVectors[level] = Opm::Lgr::reorderForOutput(levelVectors[level], toOutput_refinedLevels[level-1]);
284template <
typename Gr
id>
285void Opm::Lgr::extractRestartValueLevelGrids(
const Grid& grid,
286 const Opm::RestartValue& leafRestartValue,
287 std::vector<Opm::RestartValue>& restartValue_levels)
289 if constexpr (std::is_same_v<Grid, Dune::CpGrid>) {
291 int maxLevel = grid.maxLevel();
292 restartValue_levels.resize(maxLevel+1);
296 std::vector<std::vector<int>> toOutput_refinedLevels{};
297 toOutput_refinedLevels.resize(maxLevel);
300 for (
int level = 1; level <= maxLevel; ++level) {
301 toOutput_refinedLevels[level-1] = mapLevelIndicesToCartesianOutputOrder(grid, levelCartMapp, level);
304 std::vector<Opm::data::Solution> dataSolutionLevels{};
305 extractSolutionLevelGrids(grid,
306 toOutput_refinedLevels,
307 leafRestartValue.solution,
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,
320 for (
const auto& [rst_key, leafVector] : leafRestartValue.extra) {
322 std::vector<std::vector<double>> levelVectors{};
323 levelVectors.resize(maxLevel+1);
325 if (rst_key.key ==
"OPMEXTRA") {
330 Opm::Lgr::populateDataVectorLevelGrids<double>(grid,
333 toOutput_refinedLevels,
335 for (
int level = 0; level <= maxLevel; ++level) {
336 restartValue_levels[level].addExtra(rst_key.key, rst_key.dim, std::move(levelVectors[level]));
[ 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).
Definition CpGrid.cpp:662
int size(int level, int codim) const
Number of grid entities per level and codim.
Definition CpGrid.cpp:990
int index() const
The (positive) index of an entity.
Definition EntityRep.hpp:125
Definition Entity.hpp:107
Entity< 0 > getLevelElem() const
Get equivalent element on the level grid where the entity was born, if grid = leaf-grid-view....
Definition Entity.hpp:629
int level() const
Return the level of the entity in the grid hierarchy. Level = 0 represents the coarsest grid.
Definition Entity.hpp:473
bool isLeaf() const
Check if the entity is in the leafview.
Definition Entity.hpp:497
Definition LevelCartesianIndexMapper.hpp:45
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:71
Definition LgrOutputHelpers.hpp:74
Definition LgrOutputHelpers.hpp:57
Definition LgrOutputHelpers.hpp:40