27#ifndef OPM_FLOW_GENERIC_VANGUARD_HPP
28#define OPM_FLOW_GENERIC_VANGUARD_HPP
30#include <dune/common/parallel/communication.hh>
32#include <opm/common/OpmLog/OpmLog.hpp>
34#include <opm/grid/common/GridEnums.hpp>
36#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
38#include <opm/simulators/utils/ParallelCommunication.hpp>
44#include <unordered_map>
48namespace Opm::Parameters {
65struct MetisParams {
static constexpr auto value =
"default"; };
67#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
68struct NumJacobiBlocks {
static constexpr int value = 0; };
78struct AddCorners {
static constexpr bool value =
false; };
96namespace Action {
class State; }
99struct NumericalAquiferCell;
110 using ParallelWellStruct = std::vector<std::pair<std::string,bool>>;
114 std::unique_ptr<UDQState> udqState_;
115 std::unique_ptr<Action::State> actionState_;
116 std::unique_ptr<SummaryState> summaryState_;
117 std::unique_ptr<WellTestState> wtestState_;
118 std::shared_ptr<EclipseState> eclState_;
119 std::shared_ptr<Schedule> eclSchedule_;
120 std::shared_ptr<SummaryConfig> eclSummaryConfig_;
152 {
return setupTime_; }
158 static void readDeck(
const std::string& filename);
169 {
return *eclState_; }
172 {
return *eclState_; }
178 {
return *eclSchedule_; }
181 {
return *eclSchedule_; }
188 {
return *eclSummaryConfig_; }
198 {
return *summaryState_; }
201 {
return *summaryState_; }
209 {
return *actionState_; }
212 {
return *actionState_; }
220 {
return *udqState_; }
223 {
return *udqState_; }
225 std::unique_ptr<WellTestState> transferWTestState() {
226 return std::move(this->wtestState_);
237 {
return caseName_; }
243 {
return edgeWeightsMethod_; }
250#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
251 return numJacobiBlocks_;
261 {
return ownersFirst_; }
263 bool edgeConformal()
const
264 {
return edgeConformal_; }
267 bool addCorners()
const
268 {
return addCorners_; }
270 int numOverlap()
const
271 {
return numOverlap_; }
276 Dune::PartitionMethod partitionMethod()
const
277 {
return partitionMethod_; }
283 bool serialPartitioning()
const
284 {
return serialPartitioning_; }
290 double imbalanceTol()
const
292 if (zoltanImbalanceTolSet_) {
293 OpmLog::info(
"The parameter --zoltan-imbalance-tol is deprecated "
294 "and has been renamed to --imbalance-tol, please "
295 "adjust your calls and scripts!");
296 return zoltanImbalanceTol_;
298 return imbalanceTol_;
302 const std::string& externalPartitionFile()
const
304 return this->externalPartitionFile_;
312 {
return enableDistributedWells_; }
319 {
return enableEclOutput_; }
333 { comm_ = std::move(
comm); }
336 static Parallel::Communication&
comm()
338 assert(comm_ !=
nullptr);
344 template<
class Serializer>
345 void serializeOp(Serializer& serializer);
351 void updateOutputDir_(std::string outputDir,
352 bool enableEclCompatFile);
354 void updateNOSIM_(std::string_view enableDryRun);
356 bool drsdtconEnabled()
const;
358 std::unordered_map<std::size_t, const NumericalAquiferCell*> allAquiferCells()
const;
362 template<
class Scalar>
363 static void registerParameters_();
368 static std::unique_ptr<Parallel::Communication> comm_;
370 std::string caseName_;
371 std::string fileName_;
372 Dune::EdgeWeightMethod edgeWeightsMethod_;
374#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
375 int numJacobiBlocks_{0};
385 Dune::PartitionMethod partitionMethod_;
386 bool serialPartitioning_;
387 double imbalanceTol_;
389 bool zoltanImbalanceTolSet_;
390 double zoltanImbalanceTol_;
391 double zoltanPhgEdgeSizeThreshold_;
392 std::string zoltanParams_;
394 std::string metisParams_;
396 std::string externalPartitionFile_{};
399 bool enableDistributedWells_;
400 bool enableEclOutput_;
401 bool allow_splitting_inactive_wells_;
403 std::string ignoredKeywords_;
404 std::optional<int> outputInterval_;
405 bool useMultisegmentWell_;
406 bool enableExperiments_;
408 std::unique_ptr<SummaryState> summaryState_;
409 std::unique_ptr<UDQState> udqState_;
410 std::unique_ptr<Action::State> actionState_;
416 std::unique_ptr<WellTestState> wtestState_;
420 std::shared_ptr<Python> python;
422 std::shared_ptr<EclipseState> eclState_;
423 std::shared_ptr<Schedule> eclSchedule_;
424 std::shared_ptr<SummaryConfig> eclSummaryConfig_;
Definition FlowGenericVanguard.hpp:108
static Parallel::Communication & comm()
Obtain global communicator.
Definition FlowGenericVanguard.hpp:336
int numJacobiBlocks() const
Number of blocks in the Block-Jacobi preconditioner.
Definition FlowGenericVanguard.hpp:248
const Schedule & schedule() const
Return a reference to the object that managages the ECL schedule.
Definition FlowGenericVanguard.hpp:177
double setupTime()
Returns the wall time required to set up the simulator before it was born.
Definition FlowGenericVanguard.hpp:151
static void readDeck(const std::string &filename)
Read a deck.
Definition FlowGenericVanguard.cpp:210
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:187
ParallelWellStruct parallelWells_
Information about wells in parallel.
Definition FlowGenericVanguard.hpp:431
static void setCommunication(std::unique_ptr< Opm::Parallel::Communication > comm)
Set global communication.
Definition FlowGenericVanguard.hpp:332
bool enableDistributedWells() const
Whether perforations of a well might be distributed.
Definition FlowGenericVanguard.hpp:311
static std::string canonicalDeckPath(const std::string &caseName)
Returns the canonical path to a deck file.
Definition FlowGenericVanguard.cpp:227
~FlowGenericVanguard()
Destructor.
Action::State & actionState()
Returns the action state.
Definition FlowGenericVanguard.hpp:208
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition FlowGenericVanguard.hpp:168
FlowGenericVanguard()
Constructor.
Definition FlowGenericVanguard.cpp:108
void defineSimulationModel(SimulationModelParams &¶ms)
Set the simulation configuration objects.
Definition FlowGenericVanguard.cpp:198
bool enableEclOutput() const
Whether or not to emit result files that are compatible with a commercial reservoir simulator.
Definition FlowGenericVanguard.hpp:318
const std::string & caseName() const
Returns the name of the case.
Definition FlowGenericVanguard.hpp:236
bool ownersFirst() const
Parameter that decide if cells owned by rank are ordered before ghost cells.
Definition FlowGenericVanguard.hpp:260
Dune::EdgeWeightMethod edgeWeightsMethod() const
Parameter deciding the edge-weight strategy of the load balancer.
Definition FlowGenericVanguard.hpp:242
UDQState & udqState()
Returns the udq state.
Definition FlowGenericVanguard.hpp:219
const ParallelWellStruct & parallelWells() const
Retrieve collection (a vector of pairs) of well names and whether or not the corresponding well objec...
Definition FlowGenericVanguard.hpp:328
SummaryState & summaryState()
Returns the summary state.
Definition FlowGenericVanguard.hpp:197
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
Definition FlowGenericVanguard.hpp:112
Definition FlowGenericVanguard.hpp:73
Definition FlowGenericVanguard.hpp:78
Definition FlowGenericVanguard.hpp:50
Definition FlowGenericVanguard.hpp:51
Definition FlowGenericVanguard.hpp:53
Definition FlowGenericVanguard.hpp:54
Definition FlowGenericVanguard.hpp:55
Definition FlowGenericVanguard.hpp:56
Definition FlowGenericVanguard.hpp:57
Definition FlowGenericVanguard.hpp:58
Definition FlowGenericVanguard.hpp:63
Definition FlowGenericVanguard.hpp:61
Definition FlowGenericVanguard.hpp:65
Definition FlowGenericVanguard.hpp:79
Definition FlowGenericVanguard.hpp:71
Definition FlowGenericVanguard.hpp:72
0: simple, 1: Zoltan, 2: METIS, 3: Zoltan with a all cells of a well represented by one vertex in the...
Definition FlowGenericVanguard.hpp:77
Definition FlowGenericVanguard.hpp:82
Definition FlowGenericVanguard.hpp:83
Definition FlowGenericVanguard.hpp:86
Definition FlowGenericVanguard.hpp:90
Definition FlowGenericVanguard.hpp:88