34#ifndef OPM_LOOKUPDATA_HH 
   35#define OPM_LOOKUPDATA_HH 
   37#include <dune/grid/common/mcmgmapper.hh> 
   42#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp> 
   63template <
typename Gr
id, 
typename Gr
idView>
 
   73    explicit LookUpData(
const GridView& gridView, 
bool isFieldPropInLgr = 
false) :
 
   91    template<
class ElementOrIndex, 
class FieldProperties>
 
   92    auto operator()(
const ElementOrIndex& elementOrIdx, 
const FieldProperties& fieldProp) 
const;
 
   96                                                     const std::string& propString) 
const;
 
   99    template<
typename IntType>
 
  101                                                   const std::string& propString,
 
  102                                                   const bool& needsTranslation,
 
  103                                                   std::function<
void(IntType, 
int)> valueCheck = [](IntType, 
int){}) 
const;
 
  106    template<
typename ElemOrIndex>
 
  108                           const std::string& propString,
 
  109                           const ElemOrIndex& elemOrIndex) 
const;
 
  112    template<
typename ElemOrIndex>
 
  113    int fieldPropInt(
const FieldPropsManager& fieldPropsManager,
 
  114                     const std::string& propString,
 
  115                     const ElemOrIndex& elemOrIndex) 
const;
 
  118    template<
typename ElementType>
 
  154    template<
typename Gr
idType, 
typename ElementType>
 
  169template<
typename Gr
id, 
typename Gr
idView>
 
  183                                 bool isFieldPropInLgr = 
false) :
 
  199    template<
class ElementOrIndex, 
class FieldProperties>
 
  200    auto operator()(
const ElementOrIndex& elemIdx, 
const FieldProperties& fieldProp) 
const;
 
  205                                                     const std::string& propString) 
const;
 
  208    template<
typename IntType>
 
  210                                                   const std::string& propString,
 
  211                                                   const bool& needsTranslation,
 
  212                                                   std::function<
void(IntType, 
int)> valueCheck = [](IntType, 
int){}) 
const;
 
  215    template<
typename ElemOrIndex>
 
  217                           const std::string& propString,
 
  218                           const ElemOrIndex& elemOrIndex) 
const;
 
  221    template<
typename ElemOrIndex>
 
  222    int fieldPropInt(
const FieldPropsManager& fieldPropsManager,
 
  223                     const std::string& propString,
 
  224                     const ElemOrIndex& elemOrIndex) 
const;
 
  227    template<
class ElementType>
 
  267    template<
class Gr
idType, 
class ElementType>
 
  282template<
typename Gr
id, 
typename Gr
idView>
 
  283template<
class ElementOrIndex, 
class FieldProperties>
 
  285                                                const FieldProperties& fieldProp)
 const 
  287    const auto& fieldPropIdx = this->getFieldPropIdx(elemIdx);
 
  288    assert(0 <= fieldPropIdx && 
static_cast<int>(fieldProp.size()) > fieldPropIdx);
 
  289    return fieldProp[fieldPropIdx];
 
  292template<
typename Gr
id, 
typename Gr
idView>
 
  294                                                                                 const std::string& propString)
 const 
  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)) {
 
  305        for (
const auto& element : elements(gridView_)) {
 
  306            const auto& elemIdx = this-> elemMapper_.index(element);
 
  307            const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); 
 
  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;
 
  314                fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
 
  319        for (
const auto& element : elements(gridView_)) {
 
  320            const auto& elemIdx = this-> elemMapper_.index(element);
 
  321            const auto& fieldPropIdx = this->getFieldPropIdx<Grid>(elemIdx); 
 
  322            fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx];
 
  325    return fieldPropOnLeaf;
 
  328template<
typename Gr
id, 
typename Gr
idView>
 
  329template<
typename IntType>
 
  331                                                                               const std::string& propString,
 
  332                                                                               const bool& needsTranslation,
 
  333                                                                               std::function<
void(IntType, 
int)> valueCheck)
 const 
  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); 
 
  342        fieldPropOnLeaf[elemIdx] = fieldProp[fieldPropIdx] - needsTranslation;
 
  343        valueCheck(fieldProp[fieldPropIdx], fieldPropIdx);
 
  345    return fieldPropOnLeaf;
 
  348template<
typename Gr
id, 
typename Gr
idView>
 
  349template<
typename ElemOrIndex>
 
  351                                                       const std::string& propString,
 
  352                                                       const ElemOrIndex& elemOrIndex)
 const 
  354    const auto& fieldPropVec = fieldPropsManager.get_double(propString);
 
  355    return this ->operator()(elemOrIndex,fieldPropVec);
 
  358template<
typename Gr
id, 
typename Gr
idView>
 
  359template<
typename ElemOrIndex>
 
  361                                                 const std::string& propString,
 
  362                                                 const ElemOrIndex& elemOrIndex)
 const 
  364    const auto& fieldPropVec = fieldPropsManager.get_int(propString);
 
  365    return this ->operator()(elemOrIndex,fieldPropVec);
 
  368template<
typename Gr
id, 
typename Gr
idView>
 
  369template<
typename ElementType>
 
  371    return this->
template getFieldPropIdx<Grid, ElementType>(elementOrIndex);
 
  374template<
typename Gr
id, 
typename Gr
idView>
 
  375template<
typename Gr
idType, 
typename IndexType>
 
  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) {
 
  383            if (isFieldPropInLgr_ && elem.level()) { 
 
  385                return elem.getLevelElem().index();
 
  388                return elem.getOrigin().index();
 
  391            static_assert(std::is_same_v<IndexType, Dune::cpgrid::Entity<0>>);
 
  392            if (isFieldPropInLgr_ && elementOrIndex.level()) { 
 
  394                return elementOrIndex.getLevelElem().index();
 
  397                return elementOrIndex.getOrigin().index();
 
  401        static_assert(std::is_same_v<Grid,GridType>);
 
  402        if constexpr (isIntegral) {
 
  404            assert(gridView_.grid().maxLevel() == 0);
 
  405            return elementOrIndex;
 
  407            assert(elementOrIndex.level() == 0); 
 
  408            return this-> elemMapper_.index(elementOrIndex);
 
  416template<
typename Gr
id, 
typename Gr
idView>
 
  417template<
typename IndexType, 
typename FieldPropType>
 
  419                                                         const FieldPropType& fieldProp)
 const 
  422    const auto fieldPropCartIdx = this->getFieldPropCartesianIdx<Grid>(elementOrIndex);
 
  423    assert(0 <=  fieldPropCartIdx && (
static_cast<int>(fieldProp.size()) > fieldPropCartIdx));
 
  424    return fieldProp[fieldPropCartIdx];
 
  427template<
typename Gr
id, 
typename Gr
idView>
 
  429                                                                                          const std::string& propString)
 const 
  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];
 
  439    return fieldPropOnLeaf;
 
  442template<
typename Gr
id, 
typename Gr
idView>
 
  443template<
typename IntType>
 
  445                                                                                        const std::string& propString,
 
  446                                                                                        const bool& needsTranslation,
 
  447                                                                                        std::function<
void(IntType, 
int)> valueCheck)
 const 
  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);
 
  458    return fieldPropOnLeaf;
 
  461template<
typename Gr
id, 
typename Gr
idView>
 
  462template<
typename ElemOrIndex>
 
  464                                                                const std::string& propString,
 
  465                                                                const ElemOrIndex& elemOrIndex)
 const 
  467    const auto& fieldPropVec = fieldPropsManager.get_double(propString);
 
  468    return this ->operator()(elemOrIndex,fieldPropVec);
 
  471template<
typename Gr
id, 
typename Gr
idView>
 
  472template<
typename ElemOrIndex>
 
  474                                                          const std::string& propString,
 
  475                                                          const ElemOrIndex& elemOrIndex)
 const 
  477    const auto& fieldPropVec = fieldPropsManager.get_int(propString);
 
  478    return this ->operator()(elemOrIndex,fieldPropVec);
 
  481template<
typename Gr
id, 
typename Gr
idView>
 
  482template<
class ElementType>
 
  484    return this->getFieldPropCartesianIdx<Grid, ElementType>(elemIdx);
 
  487template<
typename Gr
id, 
typename Gr
idView>
 
  488template<
typename Gr
idType, 
typename ElementType>
 
  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) {
 
  496            return this -> getFieldPropCartesianIdx<Dune::CpGrid>(elem);
 
  498            if (isFieldPropInLgr_ && elementOrIndex.level()) { 
 
  499                return elementOrIndex.getLevelCartesianIdx();
 
  502                return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
 
  506        if constexpr (isIntegral) {
 
  507            static_assert(std::is_same_v<Grid,GridType>);
 
  508            return  cartMapper_-> cartesianIndex(elementOrIndex);
 
  510            static_assert(std::is_same_v<Grid,GridType>);
 
  512            assert(gridView_.grid().maxLevel() == 0);
 
  513            return cartMapper_-> cartesianIndex(this->elemMapper_.index(elementOrIndex));
 
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