27#ifndef OPM_GENERIC_CPGRID_VANGUARD_HPP
28#define OPM_GENERIC_CPGRID_VANGUARD_HPP
30#include <opm/grid/CpGrid.hpp>
47 class ParallelEclipseState;
58 explicit MPIPartitionFromFile(
const std::filesystem::path& partitionFile)
59 : partitionFile_(partitionFile)
62 std::vector<int> operator()(
const Dune::CpGrid& grid)
const;
65 std::filesystem::path partitionFile_{};
77template<
class ElementMapper,
class Gr
idView,
class Scalar>
81 using Element =
typename GridView::template Codim<0>::Entity;
97 const Dune::CpGrid&
grid()
const
139 const std::vector<int>& cellPartition()
const
141 return this->cell_part_;
151 void doLoadBalance_(
const Dune::EdgeWeightMethod edgeWeightsMethod,
152 const bool ownersFirst,
155 const bool enableDistributedWells,
157 const GridView& gridView,
159 EclipseState& eclState,
160 FlowGenericVanguard::ParallelWellStruct& parallelWells,
161 const int numJacobiBlocks);
168 void distributeGrid(
const Dune::EdgeWeightMethod edgeWeightsMethod,
169 const bool ownersFirst,
172 const bool enableDistributedWells,
176 const std::vector<Well>& wells,
178 EclipseState& eclState,
179 FlowGenericVanguard::ParallelWellStruct& parallelWells);
181 void distributeGrid(
const Dune::EdgeWeightMethod edgeWeightsMethod,
182 const bool ownersFirst,
185 const bool enableDistributedWells,
189 const std::vector<Well>& wells,
191 ParallelEclipseState* eclState,
192 FlowGenericVanguard::ParallelWellStruct& parallelWells);
196 virtual const std::string&
metisParams()
const = 0;
202 void doCreateGrids_(EclipseState& eclState);
205 virtual void allocTrans() = 0;
206 virtual double getTransmissibility(
unsigned I,
unsigned J)
const = 0;
209 void doFilterConnections_(
Schedule& schedule);
211 Scalar computeCellThickness(
const Element& element)
const;
213 std::unique_ptr<Dune::CpGrid> grid_;
214 std::unique_ptr<Dune::CpGrid> equilGrid_;
215 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
216 std::unique_ptr<CartesianIndexMapper> equilCartesianIndexMapper_;
219 std::vector<int> cell_part_{};
Helper class for grid instantiation of ECL file-format using problems.
Definition CollectDataOnIORank.hpp:51
Definition GenericCpGridVanguard.hpp:78
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition GenericCpGridVanguard.cpp:138
const Dune::CpGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition GenericCpGridVanguard.cpp:560
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition GenericCpGridVanguard.cpp:568
const CartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition GenericCpGridVanguard.cpp:575
const Dune::CpGrid & grid() const
Return a reference to the simulation grid.
Definition GenericCpGridVanguard.hpp:97
static void setExternalLoadBalancer(const std::function< std::vector< int >(const Dune::CpGrid &)> &loadBalancer)
Sets a function that returns external load balancing information when passed the grid.
Definition GenericCpGridVanguard.hpp:123
Dune::CpGrid & grid()
Return a reference to the simulation grid.
Definition GenericCpGridVanguard.hpp:91
void allocCartMapper()
Distribute the simulation grid over multiple processes.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
std::optional< std::function< std::vector< int >(const Dune::CpGrid &)> > externalLoadBalancer
optional functor returning external load balancing information
Definition GenericCpGridVanguard.cpp:125
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242