LookUpData.hh
Go to the documentation of this file.
1//===========================================================================
2//
3// File: LookUpData.hh
4//
5// Created: Tue May 23 14:44:00 2023
6//
7// Author(s): Antonella Ritorto <antonella.ritorto@opm-op.com>
8//
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2023 Equinor ASA.
18
19 This file is part of The Open Porous Media project (OPM).
20
21 OPM is free software: you can redistribute it and/or modify
22 it under the terms of the GNU General Public License as published by
23 the Free Software Foundation, either version 3 of the License, or
24 (at your option) any later version.
25
26 OPM is distributed in the hope that it will be useful,
27 but WITHOUT ANY WARRANTY; without even the implied warranty of
28 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
29 GNU General Public License for more details.
30
31 You should have received a copy of the GNU General Public License
32 along with OPM. If not, see <http://www.gnu.org/licenses/>.
33*/
34#ifndef OPM_LOOKUPDATA_HH
35#define OPM_LOOKUPDATA_HH
36
37#include <dune/grid/common/mcmgmapper.hh>
38
41
42#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
43
44#include <functional>
45#include <string>
46#include <type_traits>
47#include <vector>
48
49namespace Dune
50{
51class CpGrid;
52}
53
54namespace Opm
55{
56
63template <typename Grid, typename GridView>
65{
66public:
73 explicit LookUpData(const GridView& gridView, bool isFieldPropInLgr = false) :
74 gridView_(gridView),
75 elemMapper_(gridView, Dune::mcmgElementLayout()),
76 isFieldPropInLgr_(isFieldPropInLgr)
77 {
78 }
79
91 template<class ElementOrIndex, class FieldProperties>
92 auto operator()(const ElementOrIndex& elementOrIdx, const FieldProperties& fieldProp) const;
93
95 std::vector<double> assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
96 const std::string& propString) const;
97
99 template<typename IntType>
100 std::vector<IntType> assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
101 const std::string& propString,
102 const bool& needsTranslation,
103 std::function<void(IntType, int)> valueCheck = [](IntType, int){}) const;
104
106 template<typename ElemOrIndex>
107 double fieldPropDouble(const FieldPropsManager& fieldPropsManager,
108 const std::string& propString,
109 const ElemOrIndex& elemOrIndex) const;
110
112 template<typename ElemOrIndex>
113 int fieldPropInt(const FieldPropsManager& fieldPropsManager,
114 const std::string& propString,
115 const ElemOrIndex& elemOrIndex) const;
116
118 template<typename ElementType>
119 auto getFieldPropIdx(const ElementType& elem) const;
120
154 template<typename GridType, typename ElementType>
155 auto getFieldPropIdx(const ElementType& elem) const;
156
157protected:
158 const GridView& gridView_;
159 Dune::MultipleCodimMultipleGeomTypeMapper<GridView> elemMapper_;
161}; // end LookUpData class
162
169template<typename Grid, typename GridView>
171{
172public:
181 explicit LookUpCartesianData(const GridView& gridView,
183 bool isFieldPropInLgr = false) :
184 gridView_(gridView),
185 elemMapper_(gridView, Dune::mcmgElementLayout()),
186 cartMapper_(&mapper),
187 isFieldPropInLgr_(isFieldPropInLgr)
188 {
189 }
190
199 template<class ElementOrIndex, class FieldProperties>
200 auto operator()(const ElementOrIndex& elemIdx, const FieldProperties& fieldProp) const;
201
202
204 std::vector<double> assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
205 const std::string& propString) const;
206
208 template<typename IntType>
209 std::vector<IntType> assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
210 const std::string& propString,
211 const bool& needsTranslation,
212 std::function<void(IntType, int)> valueCheck = [](IntType, int){}) const;
213
215 template<typename ElemOrIndex>
216 double fieldPropDouble(const FieldPropsManager& fieldPropsManager,
217 const std::string& propString,
218 const ElemOrIndex& elemOrIndex) const;
219
221 template<typename ElemOrIndex>
222 int fieldPropInt(const FieldPropsManager& fieldPropsManager,
223 const std::string& propString,
224 const ElemOrIndex& elemOrIndex) const;
225
227 template<class ElementType>
228 auto getFieldPropCartesianIdx(const ElementType& elemIdx) const;
229
267 template<class GridType, class ElementType>
268 auto getFieldPropCartesianIdx(const ElementType& elemIdx) const;
269
270protected:
271 const GridView& gridView_;
272 Dune::MultipleCodimMultipleGeomTypeMapper<GridView> elemMapper_;
275}; // end LookUpCartesianData class
276}
277// end namespace Opm
278
279
280
282template<typename Grid, typename GridView>
283template<class ElementOrIndex, class FieldProperties>
284auto Opm::LookUpData<Grid,GridView>::operator()(const ElementOrIndex& elemIdx,
285 const FieldProperties& fieldProp) const
286{
287 const auto& fieldPropIdx = this->getFieldPropIdx(elemIdx);
288 assert(0 <= fieldPropIdx && static_cast<int>(fieldProp.size()) > fieldPropIdx);
289 return fieldProp[fieldPropIdx];
290}
291
292template<typename Grid, typename GridView>
293std::vector<double> Opm::LookUpData<Grid,GridView>::assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
294 const std::string& propString) const
295{
296 std::vector<double> fieldPropOnLeaf;
297 unsigned int numElements = gridView_.size(0);
298 fieldPropOnLeaf.resize(numElements);
299 const auto& fieldProp = fieldPropsManager.get_double(propString);
300 if ( (propString == "PORV") && (gridView_.grid().maxLevel() > 0)) {
301 // PORV poreVolume. LGRs supported (so far) only for CpGrid.
302 // For CpGrid with LGRs, poreVolume of a cell on the leaf grid view which has a parent cell on level 0,
303 // is computed as porv[parent] * leafCellVolume / parentCellVolume. In this way, the sum of the pore
304 // volume of a parent cell coincides with the sum of the pore volume of its children.
305 for (const auto& element : elements(gridView_)) {
306 const auto& elemIdx = this-> elemMapper_.index(element);
307 const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
308 if (element.hasFather()) {
309 const auto fatherVolume = element.father().geometry().volume();
310 const auto& elemVolume = element.geometry().volume();
311 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx] * elemVolume / fatherVolume;
312 }
313 else {
314 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
315 }
316 }
317 }
318 else {
319 for (const auto& element : elements(gridView_)) {
320 const auto& elemIdx = this-> elemMapper_.index(element);
321 const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
322 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
323 }
324 }
325 return fieldPropOnLeaf;
326}
327
328template<typename Grid, typename GridView>
329template<typename IntType>
330std::vector<IntType> Opm::LookUpData<Grid,GridView>::assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
331 const std::string& propString,
332 const bool& needsTranslation,
333 std::function<void(IntType, int)> valueCheck) const
334{
335 std::vector<IntType> fieldPropOnLeaf;
336 unsigned int numElements = gridView_.size(0);
337 fieldPropOnLeaf.resize(numElements);
338 const auto& fieldProp = fieldPropsManager.get_int(propString);
339 for (const auto& element : elements(gridView_)) {
340 const auto& elemIdx = this-> elemMapper_.index(element);
341 const auto& fieldPropIdx = this->getFieldPropIdx(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
342 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx] - needsTranslation;
343 valueCheck(fieldProp[fieldPropIdx], fieldPropIdx);
344 }
345 return fieldPropOnLeaf;
346}
347
348template<typename Grid, typename GridView>
349template<typename ElemOrIndex>
350double Opm::LookUpData<Grid,GridView>::fieldPropDouble(const FieldPropsManager& fieldPropsManager,
351 const std::string& propString,
352 const ElemOrIndex& elemOrIndex) const
353{
354 const auto& fieldPropVec = fieldPropsManager.get_double(propString);
355 return this ->operator()(elemOrIndex,fieldPropVec);
356}
357
358template<typename Grid, typename GridView>
359template<typename ElemOrIndex>
360int Opm::LookUpData<Grid,GridView>::fieldPropInt(const FieldPropsManager& fieldPropsManager,
361 const std::string& propString,
362 const ElemOrIndex& elemOrIndex) const
363{
364 const auto& fieldPropVec = fieldPropsManager.get_int(propString);
365 return this ->operator()(elemOrIndex,fieldPropVec);
366}
367
368template<typename Grid, typename GridView>
369template<typename ElementType>
370auto Opm::LookUpData<Grid, GridView>::getFieldPropIdx(const ElementType& elementOrIndex) const {
371 return this->template getFieldPropIdx<Grid, ElementType>(elementOrIndex);
372}
373
374template<typename Grid, typename GridView>
375template<typename GridType, typename IndexType>
376auto Opm::LookUpData<Grid,GridView>::getFieldPropIdx(const IndexType& elementOrIndex) const
377{
378 constexpr static bool isIntegral = std::is_integral_v<IndexType>;
379 if constexpr (std::is_same_v<GridType, Dune::CpGrid>) {
380 static_assert(std::is_same_v<Grid,GridType>);
381 if constexpr (isIntegral) {
382 const auto& elem = Dune::cpgrid::Entity<0>(*(gridView_.grid().currentData().back()), elementOrIndex, true);
383 if (isFieldPropInLgr_ && elem.level()) { // level > 0 == true ; level == 0 == false
384 // In case some LGRs do not have refined field properties, the next line need to be modified.
385 return elem.getLevelElem().index();
386 }
387 else {
388 return elem.getOrigin().index();
389 }
390 } else {
391 static_assert(std::is_same_v<IndexType, Dune::cpgrid::Entity<0>>);
392 if (isFieldPropInLgr_ && elementOrIndex.level()) { // level > 0 == true ; level == 0 == false
393 // In case some LGRs do not have refined field properties, the next line need to be modified.
394 return elementOrIndex.getLevelElem().index();
395 }
396 else {
397 return elementOrIndex.getOrigin().index();
398 }
399 }
400 } else {
401 static_assert(std::is_same_v<Grid,GridType>);
402 if constexpr (isIntegral) {
403 // Check there are no LGRs. LGRs (level>0) only supported for CpGrid.
404 assert(gridView_.grid().maxLevel() == 0);
405 return elementOrIndex;
406 } else {
407 assert(elementOrIndex.level() == 0); // LGRs (level>0) only supported for CpGrid.
408 return this-> elemMapper_.index(elementOrIndex);
409 }
410 }
411
412}
413
415
416template<typename Grid, typename GridView>
417template<typename IndexType, typename FieldPropType>
418auto Opm::LookUpCartesianData<Grid,GridView>::operator()(const IndexType& elementOrIndex,
419 const FieldPropType& fieldProp) const
420{
421 assert(cartMapper_);
422 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elementOrIndex);
423 assert(0 <= fieldPropCartIdx && (static_cast<int>(fieldProp.size()) > fieldPropCartIdx));
424 return fieldProp[fieldPropCartIdx];
425}
426
427template<typename Grid, typename GridView>
428std::vector<double> Opm::LookUpCartesianData<Grid,GridView>::assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
429 const std::string& propString) const
430{
431 std::vector<double> fieldPropOnLeaf;
432 unsigned int numElements = gridView_.size(0);
433 fieldPropOnLeaf.resize(numElements);
434 const auto& fieldProp = fieldPropsManager.get_double(propString);
435 for (unsigned int elemIdx = 0; elemIdx < numElements; ++elemIdx) {
436 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elemIdx);
437 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropCartIdx];
438 }
439 return fieldPropOnLeaf;
440}
441
442template<typename Grid, typename GridView>
443template<typename IntType>
444std::vector<IntType> Opm::LookUpCartesianData<Grid,GridView>::assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
445 const std::string& propString,
446 const bool& needsTranslation,
447 std::function<void(IntType, int)> valueCheck) const
448{
449 std::vector<IntType> fieldPropOnLeaf;
450 unsigned int numElements = gridView_.size(0);
451 fieldPropOnLeaf.resize(numElements);
452 const auto& fieldProp = fieldPropsManager.get_int(propString);
453 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
454 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elemIdx);
455 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropCartIdx] - needsTranslation;
456 valueCheck(fieldProp[fieldPropCartIdx], fieldPropCartIdx);
457 }
458 return fieldPropOnLeaf;
459}
460
461template<typename Grid, typename GridView>
462template<typename ElemOrIndex>
463double Opm::LookUpCartesianData<Grid,GridView>::fieldPropDouble(const FieldPropsManager& fieldPropsManager,
464 const std::string& propString,
465 const ElemOrIndex& elemOrIndex) const
466{
467 const auto& fieldPropVec = fieldPropsManager.get_double(propString);
468 return this ->operator()(elemOrIndex,fieldPropVec);
469}
470
471template<typename Grid, typename GridView>
472template<typename ElemOrIndex>
473int Opm::LookUpCartesianData<Grid,GridView>::fieldPropInt(const FieldPropsManager& fieldPropsManager,
474 const std::string& propString,
475 const ElemOrIndex& elemOrIndex) const
476{
477 const auto& fieldPropVec = fieldPropsManager.get_int(propString);
478 return this ->operator()(elemOrIndex,fieldPropVec);
479}
480
481template<typename Grid, typename GridView>
482template<class ElementType>
484 return this->getFieldPropCartesianIdx<Grid, ElementType>(elemIdx);
485}
486
487template<typename Grid, typename GridView>
488template<typename GridType, typename ElementType>
489auto Opm::LookUpCartesianData<Grid,GridView>::getFieldPropCartesianIdx(const ElementType& elementOrIndex) const
490{
491 constexpr static bool isIntegral = std::is_integral_v<ElementType>;
492 if constexpr (std::is_same_v<GridType, Dune::CpGrid>) {
493 static_assert(std::is_same_v<Grid,GridType>);
494 if constexpr (isIntegral) {
495 const auto& elem = Dune::cpgrid::Entity<0>(*(gridView_.grid().currentData().back()), elementOrIndex, true);
496 return this -> getFieldPropCartesianIdx<Dune::CpGrid>(elem);
497 } else {
498 if (isFieldPropInLgr_ && elementOrIndex.level()) { // level == 0 false; level > 0 true
499 return elementOrIndex.getLevelCartesianIdx();
500 }
501 else {
502 return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
503 }
504 }
505 } else {
506 if constexpr (isIntegral) {
507 static_assert(std::is_same_v<Grid,GridType>);
508 return cartMapper_-> cartesianIndex(elementOrIndex);
509 } else {
510 static_assert(std::is_same_v<Grid,GridType>);
511 // Check there are no LGRs. LGRs (level>0) only supported for CpGrid.
512 assert(gridView_.grid().maxLevel() == 0);
513 return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
514 }
515 }
516
517}
518#endif
Interface class to access the logical Cartesian grid as used in industry standard simulator decks.
Definition: common/CartesianIndexMapper.hpp:16
Definition: LookUpData.hh:171
const GridView & gridView_
Definition: LookUpData.hh:271
auto operator()(const ElementOrIndex &elemIdx, const FieldProperties &fieldProp) const
: Get field property for an element in the leaf grid view, from a vector, via Cartesian Index.
std::vector< double > assignFieldPropsDoubleOnLeaf(const FieldPropsManager &fieldPropsManager, const std::string &propString) const
: Get field property of type double from field properties manager by name.
Definition: LookUpData.hh:428
int fieldPropInt(const FieldPropsManager &fieldPropsManager, const std::string &propString, const ElemOrIndex &elemOrIndex) const
: Get property of type int from field properties manager by name, via element or its index.
Definition: LookUpData.hh:473
LookUpCartesianData(const GridView &gridView, const Dune::CartesianIndexMapper< Grid > &mapper, bool isFieldPropInLgr=false)
: Constructor taking a GridView, a CartesianIndexMapper, and a bool
Definition: LookUpData.hh:181
double fieldPropDouble(const FieldPropsManager &fieldPropsManager, const std::string &propString, const ElemOrIndex &elemOrIndex) const
: Get property of type double from field properties manager by name, via element or its index.
Definition: LookUpData.hh:463
auto getFieldPropCartesianIdx(const ElementType &elemIdx) const
Calls getFieldPropCartesianIdx<Grid, ElementType>(elemIdx)
Definition: LookUpData.hh:483
const Dune::CartesianIndexMapper< Grid > * cartMapper_
Definition: LookUpData.hh:273
std::vector< IntType > assignFieldPropsIntOnLeaf(const FieldPropsManager &fieldPropsManager, const std::string &propString, const bool &needsTranslation, std::function< void(IntType, int)> valueCheck=[](IntType, int){}) const
: Get field property of type int from field properties manager by name.
Definition: LookUpData.hh:444
bool isFieldPropInLgr_
Definition: LookUpData.hh:274
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > elemMapper_
Definition: LookUpData.hh:272
Definition: LookUpData.hh:65
auto operator()(const ElementOrIndex &elementOrIdx, const FieldProperties &fieldProp) const
: Get field propertry for an element or index in the leaf grid view, from a vector and element index.
Definition: LookUpData.hh:284
bool isFieldPropInLgr_
Definition: LookUpData.hh:160
std::vector< double > assignFieldPropsDoubleOnLeaf(const FieldPropsManager &fieldPropsManager, const std::string &propString) const
: Get field property of type double from field properties manager by name.
Definition: LookUpData.hh:293
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > elemMapper_
Definition: LookUpData.hh:159
const GridView & gridView_
Definition: LookUpData.hh:158
auto getFieldPropIdx(const ElementType &elem) const
Return the index used for retrieving field properties, depending on whether the grid is a CpGrid or a...
LookUpData(const GridView &gridView, bool isFieldPropInLgr=false)
: Constructor taking a GridView and a bool
Definition: LookUpData.hh:73
std::vector< IntType > assignFieldPropsIntOnLeaf(const FieldPropsManager &fieldPropsManager, const std::string &propString, const bool &needsTranslation, std::function< void(IntType, int)> valueCheck=[](IntType, int){}) const
: Get field property of type int from field properties manager by name.
Definition: LookUpData.hh:330
int fieldPropInt(const FieldPropsManager &fieldPropsManager, const std::string &propString, const ElemOrIndex &elemOrIndex) const
: Get property of type int from field properties manager by name, via element.
Definition: LookUpData.hh:360
double fieldPropDouble(const FieldPropsManager &fieldPropsManager, const std::string &propString, const ElemOrIndex &elemOrIndex) const
: Get property of type double from field properties manager by name, via element or its index.
Definition: LookUpData.hh:350
auto getFieldPropIdx(const ElementType &elem) const
calls getFieldPropIdx<Grid, ElementType>(elem)
Definition: LookUpData.hh:370
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26