27#ifndef OPM_POLYHEDRAL_GRID_VANGUARD_HPP
28#define OPM_POLYHEDRAL_GRID_VANGUARD_HPP
30#include <opm/grid/polyhedralgrid.hh>
31#include <opm/grid/polyhedralgrid/levelcartesianindexmapper.hh>
42#include <unordered_set>
45template <
class TypeTag>
46class PolyhedralGridVanguard;
49namespace Opm::Properties {
53 using InheritsFrom = std::tuple<FlowBaseVanguard>;
58template<
class TypeTag>
62template<
class TypeTag>
64 using type = Dune::PolyhedralGrid<3, 3>;
66template<
class TypeTag>
82template <
class TypeTag>
99 static constexpr int dimension = Grid::dimension;
100 static constexpr int dimensionworld = Grid::dimensionworld;
103 using GridPointer = Grid*;
104 using EquilGridPointer = EquilGrid*;
112 , simulator_(simulator)
114 this->callImplementationInit();
176 {
return *cartesianIndexMapper_; }
193 {
return *cartesianIndexMapper_; }
195 const std::vector<int>& globalCell()
200 unsigned int gridEquilIdxToGridIdx(
unsigned int elemIndex)
const {
204 unsigned int gridIdxToEquilGridIdx(
unsigned int elemIndex)
const {
217 std::unordered_set<std::string> defunctWellNames()
const
218 {
return defunctWellNames_; }
220 const TransmissibilityType& globalTransmissibility()
const
222 return simulator_.problem().eclTransmissibilities();
232 std::function<std::array<double,FlowBaseVanguard<TypeTag>::dimensionworld>(
int)>
238 std::vector<int> cellPartition()
const
247 this->grid_ = std::make_unique<Grid>
249 this->
eclState().fieldProps().porv(
true));
251 this->cartesianIndexMapper_ =
252 std::make_unique<CartesianIndexMapper>(*this->grid_);
254 this->updateGridView_();
255 this->updateCartesianToCompressedMapping_();
256 this->updateCellDepths_();
259 void filterConnections_()
264 Simulator& simulator_;
266 std::unique_ptr<Grid> grid_;
267 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
270 std::unordered_set<std::string> defunctWellNames_;
271 std::vector<int> globalcell_;
Helper class for grid instantiation of ECL file-format using problems.
Definition CollectDataOnIORank.hpp:49
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition basevanguard.hh:69
Helper class for grid instantiation of ECL file-format using problems.
Definition FlowBaseVanguard.hpp:83
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:297
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition FlowGenericVanguard.hpp:167
Definition RelpermDiagnostics.hpp:31
Helper class for grid instantiation of ECL file-format using problems.
Definition PolyhedralGridVanguard.hpp:84
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition PolyhedralGridVanguard.hpp:155
const CartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition PolyhedralGridVanguard.hpp:192
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition PolyhedralGridVanguard.hpp:213
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition PolyhedralGridVanguard.hpp:163
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition PolyhedralGridVanguard.hpp:145
Grid & grid()
Return a reference to the simulation grid.
Definition PolyhedralGridVanguard.hpp:127
const Grid & grid() const
Return a reference to the simulation grid.
Definition PolyhedralGridVanguard.hpp:133
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition PolyhedralGridVanguard.hpp:175
const LevelCartesianIndexMapper levelCartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition PolyhedralGridVanguard.hpp:183
std::function< std::array< double, FlowBaseVanguard< TypeTag >::dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition PolyhedralGridVanguard.hpp:233
Definition Transmissibility.hpp:54
Defines the common properties required by the porous medium multi-phase models.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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:235
Definition FlowBaseVanguard.hpp:69
The type of the DUNE grid.
Definition basicproperties.hh:100
Definition PolyhedralGridVanguard.hpp:52
Property which provides a Vanguard (manages grids)
Definition basicproperties.hh:96