27#ifndef OPM_FLOW_GENERIC_VANGUARD_HPP
28#define OPM_FLOW_GENERIC_VANGUARD_HPP
30#include <dune/common/parallel/communication.hh>
32#include <opm/grid/common/GridEnums.hpp>
34#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
36#include <opm/simulators/utils/ParallelCommunication.hpp>
42#include <unordered_map>
48namespace Action {
class State; }
51struct NumericalAquiferCell;
62 using ParallelWellStruct = std::vector<std::pair<std::string,bool>>;
66 std::unique_ptr<UDQState> udqState_;
67 std::unique_ptr<Action::State> actionState_;
68 std::unique_ptr<SummaryState> summaryState_;
69 std::unique_ptr<WellTestState> wtestState_;
70 std::shared_ptr<EclipseState> eclState_;
71 std::shared_ptr<Schedule> eclSchedule_;
72 std::shared_ptr<SummaryConfig> eclSummaryConfig_;
104 {
return setupTime_; }
110 static void readDeck(
const std::string& filename);
121 {
return *eclState_; }
124 {
return *eclState_; }
130 {
return *eclSchedule_; }
133 {
return *eclSchedule_; }
140 {
return *eclSummaryConfig_; }
150 {
return *summaryState_; }
153 {
return *summaryState_; }
161 {
return *actionState_; }
164 {
return *actionState_; }
172 {
return *udqState_; }
175 {
return *udqState_; }
177 WellTestState transferWTestState() {
178 return *this->wtestState_.release();
189 {
return caseName_; }
195 {
return edgeWeightsMethod_; }
202#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
203 return numJacobiBlocks_;
213 {
return ownersFirst_; }
220 bool serialPartitioning()
const
221 {
return serialPartitioning_; }
226 double zoltanImbalanceTol()
const
227 {
return zoltanImbalanceTol_; }
229 const std::string& externalPartitionFile()
const
231 return this->externalPartitionFile_;
239 {
return enableDistributedWells_; }
252 { comm_ = std::move(
comm); }
255 static Parallel::Communication&
comm()
263 template<
class Serializer>
264 void serializeOp(Serializer& serializer);
270 void updateOutputDir_(std::string outputDir,
271 bool enableEclCompatFile);
273 bool drsdtconEnabled()
const;
275 std::unordered_map<std::size_t, const NumericalAquiferCell*> allAquiferCells()
const;
282 static std::unique_ptr<Parallel::Communication> comm_;
284 std::string caseName_;
285 std::string fileName_;
286 Dune::EdgeWeightMethod edgeWeightsMethod_;
288#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
289 int numJacobiBlocks_{0};
294 bool serialPartitioning_;
295 double zoltanImbalanceTol_;
296 std::string zoltanParams_;
297 std::string externalPartitionFile_{};
299 bool enableDistributedWells_;
300 std::string ignoredKeywords_;
301 std::optional<int> outputInterval_;
302 bool useMultisegmentWell_;
303 bool enableExperiments_;
305 std::unique_ptr<SummaryState> summaryState_;
306 std::unique_ptr<UDQState> udqState_;
307 std::unique_ptr<Action::State> actionState_;
313 std::unique_ptr<WellTestState> wtestState_;
317 std::shared_ptr<Python> python;
319 std::shared_ptr<EclipseState> eclState_;
320 std::shared_ptr<Schedule> eclSchedule_;
321 std::shared_ptr<SummaryConfig> eclSummaryConfig_;
Definition FlowGenericVanguard.hpp:60
static Parallel::Communication & comm()
Obtain global communicator.
Definition FlowGenericVanguard.hpp:255
int numJacobiBlocks() const
Number of blocks in the Block-Jacobi preconditioner.
Definition FlowGenericVanguard.hpp:200
const Schedule & schedule() const
Return a reference to the object that managages the ECL schedule.
Definition FlowGenericVanguard.hpp:129
double setupTime()
Returns the wall time required to set up the simulator before it was born.
Definition FlowGenericVanguard.hpp:103
static void readDeck(const std::string &filename)
Read a deck.
Definition FlowGenericVanguard.cpp:126
const SummaryConfig & summaryConfig() const
Return a reference to the object that determines which quantities ought to be put into the ECL summar...
Definition FlowGenericVanguard.hpp:139
ParallelWellStruct parallelWells_
Information about wells in parallel.
Definition FlowGenericVanguard.hpp:328
static void setCommunication(std::unique_ptr< Opm::Parallel::Communication > comm)
Set global communication.
Definition FlowGenericVanguard.hpp:251
bool enableDistributedWells() const
Whether perforations of a well might be distributed.
Definition FlowGenericVanguard.hpp:238
static std::string canonicalDeckPath(const std::string &caseName)
Returns the canonical path to a deck file.
Definition FlowGenericVanguard.cpp:143
~FlowGenericVanguard()
Destructor.
Action::State & actionState()
Returns the action state.
Definition FlowGenericVanguard.hpp:160
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition FlowGenericVanguard.hpp:120
FlowGenericVanguard()
Constructor.
Definition FlowGenericVanguard.cpp:90
void defineSimulationModel(SimulationModelParams &¶ms)
Set the simulation configuration objects.
Definition FlowGenericVanguard.cpp:114
const std::string & caseName() const
Returns the name of the case.
Definition FlowGenericVanguard.hpp:188
bool ownersFirst() const
Parameter that decide if cells owned by rank are ordered before ghost cells.
Definition FlowGenericVanguard.hpp:212
Dune::EdgeWeightMethod edgeWeightsMethod() const
Parameter deciding the edge-weight strategy of the load balancer.
Definition FlowGenericVanguard.hpp:194
UDQState & udqState()
Returns the udq state.
Definition FlowGenericVanguard.hpp:171
const ParallelWellStruct & parallelWells() const
Returns vector with name and whether the has local perforated cells for all wells.
Definition FlowGenericVanguard.hpp:247
SummaryState & summaryState()
Returns the summary state.
Definition FlowGenericVanguard.hpp:149
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Definition FlowGenericVanguard.hpp:64