20#ifndef OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
21#define OPM_GRID_CPGRID_LGROUTPUTHELPERS_HEADER_INCLUDED
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>
37template<
class ScalarType>
41 const std::vector<std::vector<ScalarType>>& levelVectors)
43 min = std::min(
min, levelVectors[ level ][ levelIdx ]);
54template<
class ScalarType>
58 const std::vector<std::vector<ScalarType>>& levelVectors)
60 max = std::max(
max, levelVectors[ level ][ levelIdx ]);
71template<
class ScalarType>
75 const std::vector<std::vector<ScalarType>>& levelVectors)
77 partialSum+= levelVectors[ level ][ levelIdx ];
116template <
typename Container>
118 const std::vector<int>& toOutput);
126std::vector<std::unordered_map<int,int>>
139template <
typename ScalarType>
155template <
typename ScalarType>
158 const std::vector<ScalarType>& leafVector,
159 const std::vector<std::vector<int>>& toOutput_refinedLevels,
160 std::vector<std::vector<ScalarType>>& levelVectors);
183 const std::vector<std::vector<int>>& toOutput_refinedLevels,
184 const Opm::data::Solution& leafSolution,
185 std::vector<Opm::data::Solution>&);
198template <
typename Gr
id>
200 const Opm::RestartValue& leafRestartValue,
201 std::vector<Opm::RestartValue>& restartValue_levels);
206template <
typename Container>
208 const std::vector<int>& toOutput)
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]];
216 return outputContainer;
219template <
typename ScalarType>
224 const auto& [level, level_indices] = grid.
currentData()[element.
level()]->getChildrenLevelAndIndexList(element.
index());
232 auto childrenDataFunc = [](){
233 if constexpr (std::is_same_v<ScalarType, double>)
239 for (
const auto& levelIdx : level_indices) {
240 childrenDataFunc(level, levelIdx, levelVectors);
242 return childrenDataFunc.getValue();
245template <
typename ScalarType>
248 const std::vector<ScalarType>& leafVector,
249 const std::vector<std::vector<int>>& toOutput_refinedLevels,
250 std::vector<std::vector<ScalarType>>& levelVectors)
252 for (
int level = 0; level <= maxLevel; ++level) {
253 levelVectors[level].resize(grid.levelGridView(level).
size(0));
257 for (
const auto& element : Dune::elements(grid.leafGridView())) {
258 levelVectors[element.level()][element.getLevelElem().index()] = leafVector[element.index()];
263 for (
int level = maxLevel-1; level >= 0; --level) {
264 for (
const auto& element : Dune::elements(grid.levelGridView(level))) {
265 if (!element.isLeaf()) {
274 for (
int level = 1; level<=maxLevel; ++level) {
279template <
typename Gr
id>
281 const Opm::RestartValue& leafRestartValue,
282 std::vector<Opm::RestartValue>& restartValue_levels)
284 if constexpr (std::is_same_v<Grid, Dune::CpGrid>) {
286 int maxLevel = grid.maxLevel();
287 restartValue_levels.resize(maxLevel+1);
291 std::vector<std::vector<int>> toOutput_refinedLevels{};
292 toOutput_refinedLevels.resize(maxLevel);
295 for (
int level = 1; level <= maxLevel; ++level) {
299 std::vector<Opm::data::Solution> dataSolutionLevels{};
301 toOutput_refinedLevels,
302 leafRestartValue.solution,
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,
315 for (
const auto& [rst_key, leafVector] : leafRestartValue.extra) {
317 std::vector<std::vector<double>> levelVectors{};
318 levelVectors.resize(maxLevel+1);
320 if (rst_key.key ==
"OPMEXTRA") {
325 Opm::Lgr::populateDataVectorLevelGrids<double>(grid,
328 toOutput_refinedLevels,
330 for (
int level = 0; level <= maxLevel; ++level) {
331 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).
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