FlowGenericVanguard.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_FLOW_GENERIC_VANGUARD_HPP
28#define OPM_FLOW_GENERIC_VANGUARD_HPP
29
30#include <dune/common/parallel/communication.hh>
31
32#include <opm/common/OpmLog/OpmLog.hpp>
33
34#include <opm/grid/common/GridEnums.hpp>
35
36#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
37
39
40#include <cassert>
41#include <memory>
42#include <optional>
43#include <string>
44#include <unordered_map>
45#include <utility>
46#include <vector>
47
48namespace Opm::Parameters {
49
50struct AllowDistributedWells { static constexpr bool value = false; };
51struct AllowSplittingInactiveWells { static constexpr bool value = true; };
52
53struct EclOutputInterval { static constexpr int value = -1; };
54struct EdgeWeightsMethod { static constexpr auto value = "transmissibility"; };
55struct EnableDryRun { static constexpr auto value = "auto"; };
56struct EnableEclOutput { static constexpr auto value = true; };
57struct EnableOpmRstFile { static constexpr bool value = false; };
58struct ExternalPartition { static constexpr auto* value = ""; };
59
60template<class Scalar>
61struct ImbalanceTol { static constexpr Scalar value = 1.1; };
62
63struct IgnoreKeywords { static constexpr auto value = ""; };
64struct InputSkipMode { static constexpr auto value = "100"; };
65struct MetisParams { static constexpr auto value = "default"; };
66
67#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
68struct NumJacobiBlocks { static constexpr int value = 0; };
69#endif // HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
70
71struct OwnerCellsFirst { static constexpr bool value = true; };
72struct ParsingStrictness { static constexpr auto value = "normal"; };
73struct ActionParsingStrictness { static constexpr auto value = "normal"; };
74
77struct PartitionMethod { static constexpr auto value = "zoltanwell"; };
78struct AddCorners { static constexpr bool value = false; };
79struct NumOverlap { static constexpr int value = 1; };
80struct EdgeConformal { static constexpr bool value = false; };
81
82struct SchedRestart{ static constexpr bool value = false; };
83struct SerialPartitioning{ static constexpr bool value = false; };
84
85template<class Scalar>
86struct ZoltanImbalanceTol { static constexpr Scalar value = 1.1; };
87
88struct ZoltanPhgEdgeSizeThreshold { static constexpr auto value = 0.35; };
89
90struct ZoltanParams { static constexpr auto value = "graph"; };
91
92} // namespace Opm::Parameters
93
94namespace Opm {
95
96namespace Action { class State; }
97class Deck;
98class EclipseState;
99struct NumericalAquiferCell;
100class ParseContext;
101class Schedule;
102class Python;
103class SummaryConfig;
104class SummaryState;
105class UDQState;
106class WellTestState;
107
109public:
110 using ParallelWellStruct = std::vector<std::pair<std::string,bool>>;
111
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_;
121 };
122
124
131
137
139
146 static std::string canonicalDeckPath(const std::string& caseName);
147
151 double setupTime()
152 { return setupTime_; }
153
158 static void readDeck(const std::string& filename);
159
164
168 const EclipseState& eclState() const
169 { return *eclState_; }
170
171 EclipseState& eclState()
172 { return *eclState_; }
173
177 const Schedule& schedule() const
178 { return *eclSchedule_; }
179
180 Schedule& schedule()
181 { return *eclSchedule_; }
182
187 const SummaryConfig& summaryConfig() const
188 { return *eclSummaryConfig_; }
189
197 SummaryState& summaryState()
198 { return *summaryState_; }
199
200 const SummaryState& summaryState() const
201 { return *summaryState_; }
202
208 Action::State& actionState()
209 { return *actionState_; }
210
211 const Action::State& actionState() const
212 { return *actionState_; }
213
219 UDQState& udqState()
220 { return *udqState_; }
221
222 const UDQState& udqState() const
223 { return *udqState_; }
224
225 std::unique_ptr<WellTestState> transferWTestState() {
226 return std::move(this->wtestState_);
227 }
228
229
236 const std::string& caseName() const
237 { return caseName_; }
238
242 Dune::EdgeWeightMethod edgeWeightsMethod() const
243 { return edgeWeightsMethod_; }
244
248 int numJacobiBlocks() const
249 {
250#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
251 return numJacobiBlocks_;
252#else
253 return 0;
254#endif // HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
255 }
256
260 bool ownersFirst() const
261 { return ownersFirst_; }
262
263 bool edgeConformal() const
264 { return edgeConformal_; }
265
266#if HAVE_MPI
267 bool addCorners() const
268 { return addCorners_; }
269
270 int numOverlap() const
271 { return numOverlap_; }
272
276 Dune::PartitionMethod partitionMethod() const
277 { return partitionMethod_; }
278
284 { return serialPartitioning_; }
285
290 double imbalanceTol() const
291 {
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_;
297 } else {
298 return imbalanceTol_;
299 }
300 }
301
302 const std::string& externalPartitionFile() const
303 {
304 return this->externalPartitionFile_;
305 }
306#endif // HAVE_MPI
307
312 { return enableDistributedWells_; }
313
318 bool enableEclOutput() const
319 { return enableEclOutput_; }
320
329 { return parallelWells_; }
330
332 static void setCommunication(std::unique_ptr<Opm::Parallel::Communication> comm)
333 { comm_ = std::move(comm); }
334
337 {
338 assert(comm_ != nullptr);
339 return *comm_;
340 }
341
342 // Private to avoid pulling schedule in header.
343 // Static state is not serialized, only use for restart.
344 template<class Serializer>
345 void serializeOp(Serializer& serializer);
346
347 // Only compares dynamic state.
348 bool operator==(const FlowGenericVanguard& rhs) const;
349
350protected:
351 void updateOutputDir_(std::string outputDir,
352 bool enableEclCompatFile);
353
354 void updateNOSIM_(std::string_view enableDryRun);
355
356 bool drsdtconEnabled() const;
357
358 std::unordered_map<std::size_t, const NumericalAquiferCell*> allAquiferCells() const;
359
360 void init();
361
362 template<class Scalar>
363 static void registerParameters_();
364
366
367 // These variables may be owned by both Python and the simulator
368 static std::unique_ptr<Parallel::Communication> comm_;
369
370 std::string caseName_;
371 std::string fileName_;
372 Dune::EdgeWeightMethod edgeWeightsMethod_;
373
374#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
375 int numJacobiBlocks_{0};
376#endif // HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
377
380
381#if HAVE_MPI
384
385 Dune::PartitionMethod partitionMethod_;
388
392 std::string zoltanParams_;
393
394 std::string metisParams_;
395
397#endif // HAVE_MPI
398
402
403 std::string ignoredKeywords_;
404 std::optional<int> outputInterval_;
407
408 std::unique_ptr<SummaryState> summaryState_;
409 std::unique_ptr<UDQState> udqState_;
410 std::unique_ptr<Action::State> actionState_;
411
412 // Observe that this instance is handled differently from the other state
413 // variables, it will only be initialized for a restart run. While
414 // initializing a restarted run this instance is transferred to the WGState
415 // member in the well model.
416 std::unique_ptr<WellTestState> wtestState_;
417
418 // these attributes point either to the internal or to the external version of the
419 // parser objects.
420 std::shared_ptr<Python> python;
421 // These variables may be owned by both Python and the simulator
422 std::shared_ptr<EclipseState> eclState_;
423 std::shared_ptr<Schedule> eclSchedule_;
424 std::shared_ptr<SummaryConfig> eclSummaryConfig_;
425
432};
433
434} // namespace Opm
435
436#endif // OPM_FLOW_GENERIC_VANGUARD_HPP
Definition: FlowGenericVanguard.hpp:108
bool enableDistributedWells_
Definition: FlowGenericVanguard.hpp:399
bool useMultisegmentWell_
Definition: FlowGenericVanguard.hpp:405
static Parallel::Communication & comm()
Obtain global communicator.
Definition: FlowGenericVanguard.hpp:336
bool ownersFirst_
Definition: FlowGenericVanguard.hpp:378
int numJacobiBlocks() const
Number of blocks in the Block-Jacobi preconditioner.
Definition: FlowGenericVanguard.hpp:248
void updateNOSIM_(std::string_view enableDryRun)
const Schedule & schedule() const
Return a reference to the object that managages the ECL schedule.
Definition: FlowGenericVanguard.hpp:177
std::string zoltanParams_
Definition: FlowGenericVanguard.hpp:392
std::string caseName_
Definition: FlowGenericVanguard.hpp:370
std::unordered_map< std::size_t, const NumericalAquiferCell * > allAquiferCells() const
double setupTime()
Returns the wall time required to set up the simulator before it was born.
Definition: FlowGenericVanguard.hpp:151
bool addCorners() const
Definition: FlowGenericVanguard.hpp:267
std::string ignoredKeywords_
Definition: FlowGenericVanguard.hpp:403
Dune::PartitionMethod partitionMethod_
Definition: FlowGenericVanguard.hpp:385
const SummaryState & summaryState() const
Definition: FlowGenericVanguard.hpp:200
bool serialPartitioning() const
Parameter that decides if partitioning for parallel runs should be performed on a single process only...
Definition: FlowGenericVanguard.hpp:283
static void registerParameters_()
static SimulationModelParams modelParams_
Definition: FlowGenericVanguard.hpp:123
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
bool zoltanImbalanceTolSet_
Definition: FlowGenericVanguard.hpp:389
static SimulationModelParams serializationTestParams()
ParallelWellStruct parallelWells_
Information about wells in parallel.
Definition: FlowGenericVanguard.hpp:431
void serializeOp(Serializer &serializer)
bool enableExperiments_
Definition: FlowGenericVanguard.hpp:406
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
bool serialPartitioning_
Definition: FlowGenericVanguard.hpp:386
~FlowGenericVanguard()
Destructor.
Action::State & actionState()
Returns the action state.
Definition: FlowGenericVanguard.hpp:208
std::shared_ptr< SummaryConfig > eclSummaryConfig_
Definition: FlowGenericVanguard.hpp:424
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:168
double imbalanceTol() const
Parameter that sets the imbalance tolarance, depending on the chosen partition method.
Definition: FlowGenericVanguard.hpp:290
FlowGenericVanguard()
Constructor.
double zoltanPhgEdgeSizeThreshold_
Definition: FlowGenericVanguard.hpp:391
bool drsdtconEnabled() const
const std::string & externalPartitionFile() const
Definition: FlowGenericVanguard.hpp:302
Dune::EdgeWeightMethod edgeWeightsMethod_
Definition: FlowGenericVanguard.hpp:372
std::unique_ptr< WellTestState > transferWTestState()
Definition: FlowGenericVanguard.hpp:225
std::shared_ptr< Python > python
Definition: FlowGenericVanguard.hpp:420
const UDQState & udqState() const
Definition: FlowGenericVanguard.hpp:222
static std::string canonicalDeckPath(const std::string &caseName)
Returns the canonical path to a deck file.
void defineSimulationModel(SimulationModelParams &&params)
Set the simulation configuration objects.
bool enableEclOutput() const
Whether or not to emit result files that are compatible with a commercial reservoir simulator.
Definition: FlowGenericVanguard.hpp:318
bool enableEclOutput_
Definition: FlowGenericVanguard.hpp:400
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
std::unique_ptr< WellTestState > wtestState_
Definition: FlowGenericVanguard.hpp:416
UDQState & udqState()
Returns the udq state.
Definition: FlowGenericVanguard.hpp:219
Dune::PartitionMethod partitionMethod() const
Parameter deciding which partition method to use.
Definition: FlowGenericVanguard.hpp:276
bool addCorners_
Definition: FlowGenericVanguard.hpp:382
bool edgeConformal() const
Definition: FlowGenericVanguard.hpp:263
EclipseState & eclState()
Definition: FlowGenericVanguard.hpp:171
double zoltanImbalanceTol_
Definition: FlowGenericVanguard.hpp:390
static std::unique_ptr< Parallel::Communication > comm_
Definition: FlowGenericVanguard.hpp:368
bool operator==(const FlowGenericVanguard &rhs) const
static void readDeck(const std::string &filename)
Read a deck.
double imbalanceTol_
Definition: FlowGenericVanguard.hpp:387
std::unique_ptr< SummaryState > summaryState_
Definition: FlowGenericVanguard.hpp:408
std::shared_ptr< Schedule > eclSchedule_
Definition: FlowGenericVanguard.hpp:423
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
Schedule & schedule()
Definition: FlowGenericVanguard.hpp:180
std::string externalPartitionFile_
Definition: FlowGenericVanguard.hpp:396
FlowGenericVanguard(SimulationModelParams &&params)
std::vector< std::pair< std::string, bool > > ParallelWellStruct
Definition: FlowGenericVanguard.hpp:110
bool allow_splitting_inactive_wells_
Definition: FlowGenericVanguard.hpp:401
double setupTime_
Definition: FlowGenericVanguard.hpp:365
std::optional< int > outputInterval_
Definition: FlowGenericVanguard.hpp:404
const Action::State & actionState() const
Definition: FlowGenericVanguard.hpp:211
int numOverlap_
Definition: FlowGenericVanguard.hpp:383
std::unique_ptr< UDQState > udqState_
Definition: FlowGenericVanguard.hpp:409
std::string fileName_
Definition: FlowGenericVanguard.hpp:371
std::string metisParams_
Definition: FlowGenericVanguard.hpp:394
std::unique_ptr< Action::State > actionState_
Definition: FlowGenericVanguard.hpp:410
SummaryState & summaryState()
Returns the summary state.
Definition: FlowGenericVanguard.hpp:197
int numOverlap() const
Definition: FlowGenericVanguard.hpp:270
std::shared_ptr< EclipseState > eclState_
Definition: FlowGenericVanguard.hpp:422
void updateOutputDir_(std::string outputDir, bool enableEclCompatFile)
bool edgeConformal_
Definition: FlowGenericVanguard.hpp:379
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilnewtonmethodparams.hpp:31
Definition: blackoilboundaryratevector.hh:39
Definition: FlowGenericVanguard.hpp:112
std::unique_ptr< UDQState > udqState_
Definition: FlowGenericVanguard.hpp:114
std::shared_ptr< Schedule > eclSchedule_
Definition: FlowGenericVanguard.hpp:119
std::unique_ptr< SummaryState > summaryState_
Definition: FlowGenericVanguard.hpp:116
std::shared_ptr< EclipseState > eclState_
Definition: FlowGenericVanguard.hpp:118
std::shared_ptr< SummaryConfig > eclSummaryConfig_
Definition: FlowGenericVanguard.hpp:120
double setupTime_
Definition: FlowGenericVanguard.hpp:113
std::unique_ptr< WellTestState > wtestState_
Definition: FlowGenericVanguard.hpp:117
std::unique_ptr< Action::State > actionState_
Definition: FlowGenericVanguard.hpp:115
Definition: FlowGenericVanguard.hpp:73
static constexpr auto value
Definition: FlowGenericVanguard.hpp:73
Definition: FlowGenericVanguard.hpp:78
static constexpr bool value
Definition: FlowGenericVanguard.hpp:78
Definition: FlowGenericVanguard.hpp:50
static constexpr bool value
Definition: FlowGenericVanguard.hpp:50
Definition: FlowGenericVanguard.hpp:51
static constexpr bool value
Definition: FlowGenericVanguard.hpp:51
Definition: FlowGenericVanguard.hpp:53
static constexpr int value
Definition: FlowGenericVanguard.hpp:53
Definition: FlowGenericVanguard.hpp:80
static constexpr bool value
Definition: FlowGenericVanguard.hpp:80
Definition: FlowGenericVanguard.hpp:54
static constexpr auto value
Definition: FlowGenericVanguard.hpp:54
Definition: FlowGenericVanguard.hpp:55
static constexpr auto value
Definition: FlowGenericVanguard.hpp:55
Definition: FlowGenericVanguard.hpp:56
static constexpr auto value
Definition: FlowGenericVanguard.hpp:56
Definition: FlowGenericVanguard.hpp:57
static constexpr bool value
Definition: FlowGenericVanguard.hpp:57
Definition: FlowGenericVanguard.hpp:58
static constexpr auto * value
Definition: FlowGenericVanguard.hpp:58
Definition: FlowGenericVanguard.hpp:63
static constexpr auto value
Definition: FlowGenericVanguard.hpp:63
Definition: FlowGenericVanguard.hpp:61
static constexpr Scalar value
Definition: FlowGenericVanguard.hpp:61
Definition: FlowGenericVanguard.hpp:64
static constexpr auto value
Definition: FlowGenericVanguard.hpp:64
Definition: FlowGenericVanguard.hpp:65
static constexpr auto value
Definition: FlowGenericVanguard.hpp:65
Definition: FlowGenericVanguard.hpp:79
static constexpr int value
Definition: FlowGenericVanguard.hpp:79
Definition: FlowGenericVanguard.hpp:71
static constexpr bool value
Definition: FlowGenericVanguard.hpp:71
Definition: FlowGenericVanguard.hpp:72
static constexpr auto value
Definition: FlowGenericVanguard.hpp:72
Definition: FlowGenericVanguard.hpp:77
static constexpr auto value
Definition: FlowGenericVanguard.hpp:77
Definition: FlowGenericVanguard.hpp:82
static constexpr bool value
Definition: FlowGenericVanguard.hpp:82
Definition: FlowGenericVanguard.hpp:83
static constexpr bool value
Definition: FlowGenericVanguard.hpp:83
Definition: FlowGenericVanguard.hpp:86
static constexpr Scalar value
Definition: FlowGenericVanguard.hpp:86
Definition: FlowGenericVanguard.hpp:90
static constexpr auto value
Definition: FlowGenericVanguard.hpp:90
Definition: FlowGenericVanguard.hpp:88
static constexpr auto value
Definition: FlowGenericVanguard.hpp:88