27#ifndef OPM_FLOW_BASE_VANGUARD_HPP
28#define OPM_FLOW_BASE_VANGUARD_HPP
30#include <dune/grid/common/partitionset.hh>
32#include <opm/grid/common/GridEnums.hpp>
33#include <opm/grid/common/CartesianIndexMapper.hpp>
34#include <opm/grid/LookUpCellCentroid.hh>
36#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquiferCell.hpp>
37#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
45#include <opm/simulators/flow/BlackoilModelParameters.hpp>
51#include <unordered_map>
55template <
class TypeTag>
56class FlowBaseVanguard;
60namespace Opm::Properties {
69template<
class TypeTag,
class MyTypeTag>
81template <
class TypeTag>
98 static const int dimension = Grid::dimension;
99 static const int dimensionworld = Grid::dimensionworld;
100 using Element =
typename GridView::template Codim<0>::Entity;
109 FlowGenericVanguard::registerParameters_<Scalar>();
122 imbalanceTol_ = Parameters::Get<Parameters::ImbalanceTol<Scalar>>();
127 enableExperiments_ = enableExperiments;
132 const CartesianIndexMapper& cartesianMapper()
const
133 {
return asImp_().cartesianIndexMapper(); }
139 {
return asImp_().cartesianIndexMapper().cartesianDimensions(); }
145 {
return asImp_().cartesianIndexMapper().cartesianSize(); }
151 {
return asImp_().equilCartesianIndexMapper().cartesianSize(); }
166 for (
unsigned i = 1; i < dimension; ++i) {
212 throw std::runtime_error(
"compressedIndexForInteriorLGR not implemented");
222 {
return asImp_().cartesianIndexMapper().cartesianCoordinate(
cellIdx,
ijk); }
237 {
return asImp_().equilCartesianIndexMapper().cartesianCoordinate(
cellIdx,
ijk); }
251 const std::vector<Scalar>& cellCenterDepths()
const
276 const auto& grid = asImp_().grid();
277 if (grid.comm().size() == 1)
279 return grid.leafGridView().size(0);
281 const auto&
gridView = grid.leafGridView();
282 constexpr int codim = 0;
283 constexpr auto Part = Dune::Interior_Partition;
300 template<
class CartMapper>
301 std::function<std::array<double,dimensionworld>(
int)>
305 std::array<double,dimensionworld>
centroid;
306 const auto rank = this->
gridView().comm().rank();
311 centroid = this->
eclState().getInputGrid().getCellCenter(cartMapper.cartesianIndex(elemIdx));
322 void callImplementationInit()
324 asImp_().createGrids_();
325 asImp_().filterConnections_();
326 std::string outputDir = Parameters::Get<Parameters::OutputDir>();
329 const std::string&
dryRunString = Parameters::Get<Parameters::EnableDryRun>();
331 asImp_().finalizeInit_();
334 void updateCartesianToCompressedMapping_()
336 std::size_t
num_cells = asImp_().grid().leafGridView().size(0);
342 const auto elemIdx =
elemMapper.index(element);
345 if (element.partitionType() == Dune::InteriorEntity)
356 void updateCellDepths_()
358 int numCells = this->
gridView().size(0);
366 const unsigned int elemIdx =
elemMapper.index(element);
379 void updateCellThickness_()
381 if (!this->drsdtconEnabled())
399 typedef typename Element::Geometry Geometry;
400 static constexpr int zCoord = Element::dimension - 1;
403 const Geometry& geometry = element.geometry();
404 const int corners = geometry.corners();
405 for (
int i=0; i <
corners; ++i)
411 Scalar computeCellThickness(
const typename GridView::template Codim<0>::Entity& element)
const
413 typedef typename Element::Geometry Geometry;
414 static constexpr int zCoord = Element::dimension - 1;
418 const Geometry& geometry = element.geometry();
424 assert(geometry.corners() == 8);
425 for (
int i=0; i < 4; ++i){
434 Implementation& asImp_()
435 {
return *
static_cast<Implementation*
>(
this); }
437 const Implementation& asImp_()
const
438 {
return *
static_cast<const Implementation*
>(
this); }
448 mutable std::optional<std::vector<std::unordered_map<std::size_t, std::size_t>>>
lgrMappers_;
Helper class for grid instantiation of ECL file-format using problems.
Provides the base class for most (all?) simulator vanguards.
Definition CollectDataOnIORank.hpp:49
Provides the base class for most (all?) simulator vanguards.
Definition basevanguard.hh:50
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition basevanguard.hh:70
Helper class for grid instantiation of ECL file-format using problems.
Definition FlowBaseVanguard.hpp:84
unsigned cartesianIndex(const std::array< int, dimension > &coords) const
Return the index of the cells in the logical Cartesian grid.
Definition FlowBaseVanguard.hpp:162
int cartesianSize() const
Returns the overall number of cells of the logically Cartesian grid.
Definition FlowBaseVanguard.hpp:144
int compressedIndex(int cartesianCellIdx) const
Return compressed index from cartesian index.
Definition FlowBaseVanguard.hpp:180
std::vector< Scalar > cellCenterDepth_
Cell center depths.
Definition FlowBaseVanguard.hpp:451
const std::array< int, dimension > & cartesianDimensions() const
Returns the number of logically Cartesian cells in each direction.
Definition FlowBaseVanguard.hpp:138
unsigned equilCartesianIndex(unsigned compressedEquilCellIdx) const
Returns the Cartesian cell id given an element index for the grid used for equilibration.
Definition FlowBaseVanguard.hpp:227
FlowBaseVanguard(Simulator &simulator)
Create the grid for problem data files which use the ECL file format.
Definition FlowBaseVanguard.hpp:118
int compressedIndexForInterior(int cartesianCellIdx) const
Return compressed index from cartesian index only in interior.
Definition FlowBaseVanguard.hpp:195
unsigned cartesianIndex(unsigned compressedCellIdx) const
Returns the Cartesian cell id for identifaction with ECL data.
Definition FlowBaseVanguard.hpp:156
std::vector< int > is_interior_
Whether a cells is in the interior.
Definition FlowBaseVanguard.hpp:459
void cartesianCoordinate(unsigned cellIdx, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition FlowBaseVanguard.hpp:221
std::unordered_map< int, int > cartesianToCompressed_
Mapping between cartesian and compressed cells.
Definition FlowBaseVanguard.hpp:444
std::vector< Scalar > cellThickness_
Cell thickness.
Definition FlowBaseVanguard.hpp:455
std::optional< std::vector< std::unordered_map< std::size_t, std::size_t > > > lgrMappers_
Mapping between LGR cartesian and compressed cells.
Definition FlowBaseVanguard.hpp:448
Scalar cellCenterDepth(unsigned globalSpaceIdx) const
Returns the depth of a degree of freedom [m].
Definition FlowBaseVanguard.hpp:246
std::size_t globalNumCells() const
Get the number of cells in the global leaf grid view.
Definition FlowBaseVanguard.hpp:274
void equilCartesianCoordinate(unsigned cellIdx, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell of the grid used for EQUIL.
Definition FlowBaseVanguard.hpp:236
std::function< std::array< double, dimensionworld >(int)> cellCentroids_(const CartMapper &cartMapper, const bool &isCpGrid) const
Get function to query cell centroids for a distributed grid.
Definition FlowBaseVanguard.hpp:302
Scalar cellThickness(unsigned globalSpaceIdx) const
Returns the thickness of a degree of freedom [m].
Definition FlowBaseVanguard.hpp:263
static void registerParameters()
Register the common run-time parameters for all ECL simulator vanguards.
Definition FlowBaseVanguard.hpp:107
int equilCartesianSize() const
Returns the overall number of cells of the logically EquilCartesian grid.
Definition FlowBaseVanguard.hpp:150
Definition FlowGenericVanguard.hpp:108
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition FlowGenericVanguard.hpp:168
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition FlowBaseVanguard.hpp:57
Definition FlowBaseVanguard.hpp:70
Definition FlowBaseVanguard.hpp:64
a tag to mark properties as undefined
Definition propertysystem.hh:38