27#ifndef OPM_ALUGRID_VANGUARD_HPP
28#define OPM_ALUGRID_VANGUARD_HPP
30#include <dune/alugrid/common/fromtogridfactory.hh>
31#include <dune/alugrid/dgf.hh>
32#include <dune/alugrid/grid.hh>
34#include <opm/common/OpmLog/OpmLog.hpp>
36#include <opm/grid/CpGrid.hpp>
37#include <opm/grid/cpgrid/LevelCartesianIndexMapper.hpp>
41#include <opm/simulators/flow/AluGridCartesianIndexMapper.hpp>
42#include <opm/simulators/flow/AluGridLevelCartesianIndexMapper.hpp>
45#include <opm/simulators/utils/ParallelEclipseState.hpp>
54template <
class TypeTag>
59namespace Opm::Properties {
63 using InheritsFrom = std::tuple<FlowBaseVanguard>;
68template<
class TypeTag>
72template<
class TypeTag>
75 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridMPIComm>;
77 using type = Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming, Dune::ALUGridNoComm>;
80template<
class TypeTag>
82 using type = Dune::CpGrid;
96template <
class TypeTag>
114 using Factory = Dune::FromToGridFactory<Grid>;
116 static constexpr int dimension = Grid::dimension;
117 static constexpr int dimensionworld = Grid::dimensionworld;
119 explicit AluGridVanguard(Simulator& simulator)
123 this->callImplementationInit();
148 {
return *equilGrid_; }
159 delete equilCartesianIndexMapper_;
160 equilCartesianIndexMapper_ =
nullptr;
163 equilGrid_ =
nullptr;
174 auto dataHandle = cartesianIndexMapper_->dataHandle(
gridView);
175 grid().loadBalance(*dataHandle);
183 grid().communicate(*dataHandle,
184 Dune::InteriorBorder_All_Interface,
185 Dune::ForwardCommunication );
190 globalTrans_ = std::make_unique<TransmissibilityType>(this->
eclState(),
196 == EnergyModules::FullyImplicitThermal,
202 globalTrans_->update(
false, TransmissibilityType::TransUpdateQuantities::Trans,
203 [&](
unsigned int i) {
return gridEquilIdxToGridIdx(i);});
213 template<
class DataHandle>
214 void scatterData(DataHandle& )
const
219 template<
class DataHandle>
220 void gatherData(DataHandle& )
const
225 template<
class DataHandle,
class InterfaceType,
class CommunicationDirection>
226 void communicate (DataHandle& , InterfaceType ,
227 CommunicationDirection )
const
239 globalTrans_.reset();
247 {
return *cartesianIndexMapper_; }
255 {
return LevelCartesianIndexMapper(*cartesianIndexMapper_); }
261 {
return *equilCartesianIndexMapper_; }
270 std::function<std::array<double,dimensionworld>(
int)>
276 const TransmissibilityType& globalTransmissibility()
const
278 assert( globalTrans_ !=
nullptr );
279 return *globalTrans_;
282 const std::vector<int>& globalCell()
284 return cartesianCellId_;
287 std::vector<int> cellPartition()
const
293 unsigned int gridEquilIdxToGridIdx(
unsigned int elemIndex)
const {
294 return equilGridToGrid_[elemIndex];
297 unsigned int gridIdxToEquilGridIdx(
unsigned int elemIndex)
const {
298 return ordering_[elemIndex];
312 const EclipseGrid* input_grid =
nullptr;
313 std::vector<double> global_porv;
319 input_grid = &this->
eclState().getInputGrid();
320 global_porv = this->
eclState().fieldProps().porv(
true);
321 OpmLog::info(
"\nProcessing grid");
327 this->equilGrid_ = std::make_unique<Dune::CpGrid>();
331 this->equilGrid_->processEclipseFormat(input_grid,
337 cartesianCellId_ = this->equilGrid_->globalCell();
339 for (
unsigned i = 0; i < dimension; ++i)
340 cartesianDimension_[i] = this->equilGrid_->logicalCartesianSize()[i];
342 equilCartesianIndexMapper_ = std::make_unique<EquilCartesianIndexMapper>(*equilGrid_);
348 factory_ = std::make_unique<Factory>();
349 grid_ = factory_->convert(*equilGrid_, cartesianCellId_, ordering_);
350 OpmLog::warning(
"Space Filling Curve (SFC) ordering is enabled: see flow_blackoil_alugrid for more informations on disabling/enabling SFC reordering");
351 equilGridToGrid_.resize(ordering_.size());
352 for (std::size_t index = 0; index < ordering_.size(); ++index) {
353 equilGridToGrid_[ordering_[index]] = index;
356 cartesianIndexMapper_ = std::make_unique<CartesianIndexMapper>(*grid_, cartesianDimension_, cartesianCellId_);
357 this->updateGridView_();
358 this->updateCartesianToCompressedMapping_();
359 this->updateCellDepths_();
360 this->updateCellThickness_();
363 void filterConnections_()
368 std::unique_ptr<Grid> grid_;
369 std::unique_ptr<EquilGrid> equilGrid_;
370 std::vector<int> cartesianCellId_;
371 std::vector<unsigned int> ordering_;
372 std::vector<unsigned int> equilGridToGrid_;
373 std::array<int,dimension> cartesianDimension_;
374 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
375 std::unique_ptr<EquilCartesianIndexMapper> equilCartesianIndexMapper_;
376 std::unique_ptr<Factory> factory_;
381 std::unique_ptr<TransmissibilityType> globalTrans_;
Helper class for grid instantiation of ECL file-format using problems.
Definition CollectDataOnIORank.hpp:49
Helper class for grid instantiation of ECL file-format using problems.
Definition AluGridVanguard.hpp:98
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition AluGridVanguard.hpp:171
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition AluGridVanguard.hpp:246
const LevelCartesianIndexMapper levelCartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition AluGridVanguard.hpp:254
std::function< std::array< double, dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition AluGridVanguard.hpp:271
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition AluGridVanguard.hpp:237
const Grid & grid() const
Return a reference to the simulation grid.
Definition AluGridVanguard.hpp:135
const EquilCartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition AluGridVanguard.hpp:260
Grid & grid()
Return a reference to the simulation grid.
Definition AluGridVanguard.hpp:129
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition AluGridVanguard.hpp:147
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition AluGridVanguard.hpp:157
void addLgrs()
Add LGRs to the grid, if any.
Definition basevanguard.hh:117
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition basevanguard.hh:70
Definition AluGridVanguard.hpp:55
FlowBaseVanguard(Simulator &simulator)
Create the grid for problem data files which use the ECL file format.
Definition FlowBaseVanguard.hpp:121
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:305
static Parallel::Communication & comm()
Obtain global communicator.
Definition FlowGenericVanguard.hpp:336
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition FlowGenericVanguard.hpp:168
Definition EclGenericWriter.hpp:50
Definition Transmissibility.hpp:54
Defines the common properties required by the porous medium multi-phase models.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:82
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Enable diffusive fluxes?
Definition multiphasebaseproperties.hh:91
Enable dispersive fluxes?
Definition multiphasebaseproperties.hh:95
Definition FlowBaseVanguard.hpp:70
The type of the DUNE grid.
Definition basicproperties.hh:104
Definition AluGridVanguard.hpp:62
Property which provides a Vanguard (manages grids).
Definition basicproperties.hh:100