35#include <dune/grid/common/mcmgmapper.hh>
36#include <opm/input/eclipse/EclipseState/Grid/EclipseGrid.hpp>
58template <
typename Gr
id,
typename Gr
idView>
67 const Opm::EclipseGrid* eclgrid) :
69 cartMapper_(&cartMapper),
84 template<
typename Gr
idType = Gr
id>
85 typename std::enable_if_t<!std::is_same_v<GridType,Dune::CpGrid>, std::array<double,3>>
98 template<
typename Gr
idType = Gr
id>
99 typename std::enable_if_t<std::is_same_v<GridType,Dune::CpGrid>, std::array<double,3>>
103 const GridView& gridView_;
105 const Opm::EclipseGrid* eclGrid_;
112template<
typename Gr
id,
typename Gr
idView>
113template<
typename Gr
idType>
114typename std::enable_if_t<!std::is_same_v<GridType,Dune::CpGrid>, std::array<double,3>>
117 static_assert(std::is_same_v<Grid,GridType>);
118 return this -> eclGrid_ -> getCellCenter(
this -> cartMapper_->cartesianIndex(elemIdx));
121template<
typename Gr
id,
typename Gr
idView>
122template<
typename Gr
idType>
123typename std::enable_if_t<std::is_same_v<GridType,Dune::CpGrid>,std::array<double,3>>
126 static_assert(std::is_same_v<Grid,GridType>);
127 return this -> gridView_.grid().getEclCentroid(elemIdx);
Interface class to access the logical Cartesian grid as used in industry standard simulator decks.
Definition CartesianIndexMapper.hpp:16
Copyright 2019 Equinor AS.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68
LookUpCellCentroid struct - To search cell centroids via element index.
Definition LookUpCellCentroid.hh:60
std::enable_if_t<!std::is_same_v< GridType, Dune::CpGrid >, std::array< double, 3 > > operator()(std::size_t elemIdx) const
: Call operator
Definition LookUpCellCentroid.hh:115
LookUpCellCentroid(const GridView &gridView, const Dune::CartesianIndexMapper< Grid > &cartMapper, const Opm::EclipseGrid *eclgrid)
: Constructor taking a GridView, CartesianMapper
Definition LookUpCellCentroid.hh:65