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
39#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
41
42#include <functional>
43#include <string>
44#include <type_traits>
45#include <vector>
46
47namespace Dune
48{
49class CpGrid;
50}
51
52namespace Opm
53{
54
61template <typename Grid, typename GridView>
63{
64public:
71 explicit LookUpData(const GridView& gridView, bool isFieldPropInLgr = false) :
72 gridView_(gridView),
73 elemMapper_(gridView, Dune::mcmgElementLayout()),
74 isFieldPropInLgr_(isFieldPropInLgr)
75 {
76 }
77
89 template<class ElementOrIndex, class FieldProperties>
90 auto operator()(const ElementOrIndex& elementOrIdx, const FieldProperties& fieldProp) const;
91
93 std::vector<double> assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
94 const std::string& propString) const;
95
97 template<typename IntType>
98 std::vector<IntType> assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
99 const std::string& propString,
100 const bool& needsTranslation,
101 std::function<void(IntType, int)> valueCheck = [](IntType, int){}) const;
102
104 template<typename ElemOrIndex>
105 double fieldPropDouble(const FieldPropsManager& fieldPropsManager,
106 const std::string& propString,
107 const ElemOrIndex& elemOrIndex) const;
108
110 template<typename ElemOrIndex>
111 int fieldPropInt(const FieldPropsManager& fieldPropsManager,
112 const std::string& propString,
113 const ElemOrIndex& elemOrIndex) const;
114
116 template<typename ElementType>
117 auto getFieldPropIdx(const ElementType& elem) const;
118
152 template<typename GridType, typename ElementType>
153 auto getFieldPropIdx(const ElementType& elem) const;
154
155protected:
156 const GridView& gridView_;
157 Dune::MultipleCodimMultipleGeomTypeMapper<GridView> elemMapper_;
159}; // end LookUpData class
160
167template<typename Grid, typename GridView>
169{
170public:
179 explicit LookUpCartesianData(const GridView& gridView,
181 bool isFieldPropInLgr = false) :
182 gridView_(gridView),
183 elemMapper_(gridView, Dune::mcmgElementLayout()),
184 cartMapper_(&mapper),
185 isFieldPropInLgr_(isFieldPropInLgr)
186 {
187 }
188
197 template<class ElementOrIndex, class FieldProperties>
198 auto operator()(const ElementOrIndex& elemIdx, const FieldProperties& fieldProp) const;
199
200
202 std::vector<double> assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
203 const std::string& propString) const;
204
206 template<typename IntType>
207 std::vector<IntType> assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
208 const std::string& propString,
209 const bool& needsTranslation,
210 std::function<void(IntType, int)> valueCheck = [](IntType, int){}) const;
211
213 template<typename ElemOrIndex>
214 double fieldPropDouble(const FieldPropsManager& fieldPropsManager,
215 const std::string& propString,
216 const ElemOrIndex& elemOrIndex) const;
217
219 template<typename ElemOrIndex>
220 int fieldPropInt(const FieldPropsManager& fieldPropsManager,
221 const std::string& propString,
222 const ElemOrIndex& elemOrIndex) const;
223
225 template<class ElementType>
226 auto getFieldPropCartesianIdx(const ElementType& elemIdx) const;
227
265 template<class GridType, class ElementType>
266 auto getFieldPropCartesianIdx(const ElementType& elemIdx) const;
267
268protected:
269 const GridView& gridView_;
270 Dune::MultipleCodimMultipleGeomTypeMapper<GridView> elemMapper_;
273}; // end LookUpCartesianData class
274}
275// end namespace Opm
276
277
278
280template<typename Grid, typename GridView>
281template<class ElementOrIndex, class FieldProperties>
282auto Opm::LookUpData<Grid,GridView>::operator()(const ElementOrIndex& elemIdx,
283 const FieldProperties& fieldProp) const
284{
285 const auto& fieldPropIdx = this->getFieldPropIdx(elemIdx);
286 assert(0 <= fieldPropIdx && static_cast<int>(fieldProp.size()) > fieldPropIdx);
287 return fieldProp[fieldPropIdx];
288}
289
290template<typename Grid, typename GridView>
291std::vector<double> Opm::LookUpData<Grid,GridView>::assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
292 const std::string& propString) const
293{
294 std::vector<double> fieldPropOnLeaf;
295 unsigned int numElements = gridView_.size(0);
296 fieldPropOnLeaf.resize(numElements);
297 const auto& fieldProp = fieldPropsManager.get_double(propString);
298 if ( (propString == "PORV") && (gridView_.grid().maxLevel() > 0)) {
299 // PORV poreVolume. LGRs supported (so far) only for CpGrid.
300 // For CpGrid with LGRs, poreVolume of a cell on the leaf grid view which has a parent cell on level 0,
301 // is computed as porv[parent] * leafCellVolume / parentCellVolume. In this way, the sum of the pore
302 // volume of a parent cell coincides with the sum of the pore volume of its children.
303 for (const auto& element : elements(gridView_)) {
304 const auto& elemIdx = this-> elemMapper_.index(element);
305 const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
306 if (element.hasFather()) {
307 const auto fatherVolume = element.father().geometry().volume();
308 const auto& elemVolume = element.geometry().volume();
309 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx] * elemVolume / fatherVolume;
310 }
311 else {
312 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
313 }
314 }
315 }
316 else {
317 for (const auto& element : elements(gridView_)) {
318 const auto& elemIdx = this-> elemMapper_.index(element);
319 const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
320 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
321 }
322 }
323 return fieldPropOnLeaf;
324}
325
326template<typename Grid, typename GridView>
327template<typename IntType>
328std::vector<IntType> Opm::LookUpData<Grid,GridView>::assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
329 const std::string& propString,
330 const bool& needsTranslation,
331 std::function<void(IntType, int)> valueCheck) const
332{
333 std::vector<IntType> fieldPropOnLeaf;
334 unsigned int numElements = gridView_.size(0);
335 fieldPropOnLeaf.resize(numElements);
336 const auto& fieldProp = fieldPropsManager.get_int(propString);
337 for (const auto& element : elements(gridView_)) {
338 const auto& elemIdx = this-> elemMapper_.index(element);
339 const auto& fieldPropIdx = this->getFieldPropIdx(elemIdx); // gets parentIdx (or (lgr)levelIdx) for CpGrid with LGRs
340 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx] - needsTranslation;
341 valueCheck(fieldProp[fieldPropIdx], fieldPropIdx);
342 }
343 return fieldPropOnLeaf;
344}
345
346template<typename Grid, typename GridView>
347template<typename ElemOrIndex>
348double Opm::LookUpData<Grid,GridView>::fieldPropDouble(const FieldPropsManager& fieldPropsManager,
349 const std::string& propString,
350 const ElemOrIndex& elemOrIndex) const
351{
352 const auto& fieldPropVec = fieldPropsManager.get_double(propString);
353 return this ->operator()(elemOrIndex,fieldPropVec);
354}
355
356template<typename Grid, typename GridView>
357template<typename ElemOrIndex>
358int Opm::LookUpData<Grid,GridView>::fieldPropInt(const FieldPropsManager& fieldPropsManager,
359 const std::string& propString,
360 const ElemOrIndex& elemOrIndex) const
361{
362 const auto& fieldPropVec = fieldPropsManager.get_int(propString);
363 return this ->operator()(elemOrIndex,fieldPropVec);
364}
365
366template<typename Grid, typename GridView>
367template<typename ElementType>
368auto Opm::LookUpData<Grid, GridView>::getFieldPropIdx(const ElementType& elementOrIndex) const {
369 return this->template getFieldPropIdx<Grid, ElementType>(elementOrIndex);
370}
371
372template<typename Grid, typename GridView>
373template<typename GridType, typename IndexType>
374auto Opm::LookUpData<Grid,GridView>::getFieldPropIdx(const IndexType& elementOrIndex) const
375{
376 constexpr static bool isIntegral = std::is_integral_v<IndexType>;
377 if constexpr (std::is_same_v<GridType, Dune::CpGrid>) {
378 static_assert(std::is_same_v<Grid,GridType>);
379 if constexpr (isIntegral) {
380 const auto& elem = Dune::cpgrid::Entity<0>(*(gridView_.grid().currentData().back()), elementOrIndex, true);
381 if (isFieldPropInLgr_ && elem.level()) { // level > 0 == true ; level == 0 == false
382 // In case some LGRs do not have refined field properties, the next line need to be modified.
383 return elem.getLevelElem().index();
384 }
385 else {
386 return elem.getOrigin().index();
387 }
388 } else {
389 static_assert(std::is_same_v<IndexType, Dune::cpgrid::Entity<0>>);
390 if (isFieldPropInLgr_ && elementOrIndex.level()) { // level > 0 == true ; level == 0 == false
391 // In case some LGRs do not have refined field properties, the next line need to be modified.
392 return elementOrIndex.getLevelElem().index();
393 }
394 else {
395 return elementOrIndex.getOrigin().index();
396 }
397 }
398 } else {
399 static_assert(std::is_same_v<Grid,GridType>);
400 if constexpr (isIntegral) {
401 // Check there are no LGRs. LGRs (level>0) only supported for CpGrid.
402 assert(gridView_.grid().maxLevel() == 0);
403 return elementOrIndex;
404 } else {
405 assert(elementOrIndex.level() == 0); // LGRs (level>0) only supported for CpGrid.
406 return this-> elemMapper_.index(elementOrIndex);
407 }
408 }
409
410}
411
413
414template<typename Grid, typename GridView>
415template<typename IndexType, typename FieldPropType>
416auto Opm::LookUpCartesianData<Grid,GridView>::operator()(const IndexType& elementOrIndex,
417 const FieldPropType& fieldProp) const
418{
419 assert(cartMapper_);
420 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elementOrIndex);
421 assert(0 <= fieldPropCartIdx && (static_cast<int>(fieldProp.size()) > fieldPropCartIdx));
422 return fieldProp[fieldPropCartIdx];
423}
424
425template<typename Grid, typename GridView>
426std::vector<double> Opm::LookUpCartesianData<Grid,GridView>::assignFieldPropsDoubleOnLeaf(const FieldPropsManager& fieldPropsManager,
427 const std::string& propString) const
428{
429 std::vector<double> fieldPropOnLeaf;
430 unsigned int numElements = gridView_.size(0);
431 fieldPropOnLeaf.resize(numElements);
432 const auto& fieldProp = fieldPropsManager.get_double(propString);
433 for (unsigned int elemIdx = 0; elemIdx < numElements; ++elemIdx) {
434 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elemIdx);
435 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropCartIdx];
436 }
437 return fieldPropOnLeaf;
438}
439
440template<typename Grid, typename GridView>
441template<typename IntType>
442std::vector<IntType> Opm::LookUpCartesianData<Grid,GridView>::assignFieldPropsIntOnLeaf(const FieldPropsManager& fieldPropsManager,
443 const std::string& propString,
444 const bool& needsTranslation,
445 std::function<void(IntType, int)> valueCheck) const
446{
447 std::vector<IntType> fieldPropOnLeaf;
448 unsigned int numElements = gridView_.size(0);
449 fieldPropOnLeaf.resize(numElements);
450 const auto& fieldProp = fieldPropsManager.get_int(propString);
451 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
452 const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elemIdx);
453 fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropCartIdx] - needsTranslation;
454 valueCheck(fieldProp[fieldPropCartIdx], fieldPropCartIdx);
455 }
456 return fieldPropOnLeaf;
457}
458
459template<typename Grid, typename GridView>
460template<typename ElemOrIndex>
461double Opm::LookUpCartesianData<Grid,GridView>::fieldPropDouble(const FieldPropsManager& fieldPropsManager,
462 const std::string& propString,
463 const ElemOrIndex& elemOrIndex) const
464{
465 const auto& fieldPropVec = fieldPropsManager.get_double(propString);
466 return this ->operator()(elemOrIndex,fieldPropVec);
467}
468
469template<typename Grid, typename GridView>
470template<typename ElemOrIndex>
471int Opm::LookUpCartesianData<Grid,GridView>::fieldPropInt(const FieldPropsManager& fieldPropsManager,
472 const std::string& propString,
473 const ElemOrIndex& elemOrIndex) const
474{
475 const auto& fieldPropVec = fieldPropsManager.get_int(propString);
476 return this ->operator()(elemOrIndex,fieldPropVec);
477}
478
479template<typename Grid, typename GridView>
480template<class ElementType>
482 return this->getFieldPropCartesianIdx<Grid, ElementType>(elemIdx);
483}
484
485template<typename Grid, typename GridView>
486template<typename GridType, typename ElementType>
487auto Opm::LookUpCartesianData<Grid,GridView>::getFieldPropCartesianIdx(const ElementType& elementOrIndex) const
488{
489 constexpr static bool isIntegral = std::is_integral_v<ElementType>;
490 if constexpr (std::is_same_v<GridType, Dune::CpGrid>) {
491 static_assert(std::is_same_v<Grid,GridType>);
492 if constexpr (isIntegral) {
493 const auto& elem = Dune::cpgrid::Entity<0>(*(gridView_.grid().currentData().back()), elementOrIndex, true);
494 return this -> getFieldPropCartesianIdx<Dune::CpGrid>(elem);
495 } else {
496 if (isFieldPropInLgr_ && elementOrIndex.level()) { // level == 0 false; level > 0 true
497 return elementOrIndex.getLevelCartesianIdx();
498 }
499 else {
500 return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
501 }
502 }
503 } else {
504 if constexpr (isIntegral) {
505 static_assert(std::is_same_v<Grid,GridType>);
506 return cartMapper_-> cartesianIndex(elementOrIndex);
507 } else {
508 static_assert(std::is_same_v<Grid,GridType>);
509 // Check there are no LGRs. LGRs (level>0) only supported for CpGrid.
510 assert(gridView_.grid().maxLevel() == 0);
511 return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
512 }
513 }
514
515}
516#endif
Interface class to access the logical Cartesian grid as used in industry standard simulator decks.
Definition: common/CartesianIndexMapper.hpp:16
Definition: LookUpData.hh:169
const GridView & gridView_
Definition: LookUpData.hh:269
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:426
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:471
LookUpCartesianData(const GridView &gridView, const Dune::CartesianIndexMapper< Grid > &mapper, bool isFieldPropInLgr=false)
: Constructor taking a GridView, a CartesianIndexMapper, and a bool
Definition: LookUpData.hh:179
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:461
auto getFieldPropCartesianIdx(const ElementType &elemIdx) const
Calls getFieldPropCartesianIdx<Grid, ElementType>(elemIdx)
Definition: LookUpData.hh:481
const Dune::CartesianIndexMapper< Grid > * cartMapper_
Definition: LookUpData.hh:271
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:442
bool isFieldPropInLgr_
Definition: LookUpData.hh:272
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > elemMapper_
Definition: LookUpData.hh:270
Definition: LookUpData.hh:63
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:282
bool isFieldPropInLgr_
Definition: LookUpData.hh:158
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:291
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > elemMapper_
Definition: LookUpData.hh:157
const GridView & gridView_
Definition: LookUpData.hh:156
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:71
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:328
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:358
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:348
auto getFieldPropIdx(const ElementType &elem) const
calls getFieldPropIdx<Grid, ElementType>(elem)
Definition: LookUpData.hh:368
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