opm-common
Loading...
Searching...
No Matches
Runspec.hpp
1/*
2 Copyright 2016 Statoil ASA.
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 it under the terms
7 of the GNU General Public License as published by the Free Software
8 Foundation, either version 3 of the License, or (at your option) any later
9 version.
10
11 OPM is distributed in the hope that it will be useful, but WITHOUT ANY
12 WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
13 A PARTICULAR PURPOSE. See the GNU General Public License for more details.
14
15 You should have received a copy of the GNU General Public License along with
16 OPM. If not, see <http://www.gnu.org/licenses/>.
17*/
18
19#ifndef OPM_RUNSPEC_HPP
20#define OPM_RUNSPEC_HPP
21
22#include <opm/common/OpmLog/KeywordLocation.hpp>
23
24#include <opm/input/eclipse/EclipseState/EndpointScaling.hpp>
25#include <opm/input/eclipse/EclipseState/Phase.hpp>
26#include <opm/input/eclipse/EclipseState/Tables/Regdims.hpp>
27#include <opm/input/eclipse/EclipseState/Tables/Tabdims.hpp>
28
29#include <opm/input/eclipse/Schedule/Action/Actdims.hpp>
30#include <opm/input/eclipse/Schedule/UDQ/UDQParams.hpp>
31
32#include <bitset>
33#include <cstddef>
34#include <ctime>
35#include <optional>
36
37namespace Opm {
38
39 class Deck;
40
41} // namespace Opm
42
43namespace Opm {
44
45class Phases
46{
47public:
48 Phases() noexcept = default;
49 Phases(bool oil, bool gas, bool wat,
50 bool solvent = false,
51 bool polymer = false,
52 bool energy = false,
53 bool polymw = false,
54 bool foam = false,
55 bool brine = false,
56 bool zfraction = false) noexcept;
57
58 static Phases serializationTestObject();
59
60 bool active( Phase ) const noexcept;
61 size_t size() const noexcept;
62
63 bool operator==(const Phases& data) const;
64
65 template<class Serializer>
66 void serializeOp(Serializer& serializer)
67 {
68 serializer(bits);
69 }
70
71private:
72 std::bitset<NUM_PHASES_IN_ENUM> bits;
73};
74
75class Welldims {
76public:
77 Welldims() = default;
78 explicit Welldims(const Deck& deck);
79
80 static Welldims serializationTestObject();
81
82 int maxConnPerWell() const
83 {
84 return this->nCWMax;
85 }
86
87 int maxWellsPerGroup() const
88 {
89 return this->nWGMax;
90 }
91
92 int maxGroupsInField() const
93 {
94 return this->nGMax;
95 }
96
97 int maxWellsInField() const
98 {
99 return this->nWMax;
100 }
101
102 int maxWellListsPrWell() const
103 {
104 return this->nWlistPrWellMax;
105 }
106
107 int maxDynamicWellLists() const
108 {
109 return this->nDynWlistMax;
110 }
111
112 const std::optional<KeywordLocation>& location() const
113 {
114 return this->m_location;
115 }
116
117 static bool rst_cmp(const Welldims& full_dims, const Welldims& rst_dims) {
118 return full_dims.maxConnPerWell() == rst_dims.maxConnPerWell() &&
119 full_dims.maxWellsPerGroup() == rst_dims.maxWellsPerGroup() &&
120 full_dims.maxGroupsInField() == rst_dims.maxGroupsInField() &&
121 full_dims.maxWellsInField() == rst_dims.maxWellsInField() &&
122 full_dims.maxWellListsPrWell() == rst_dims.maxWellListsPrWell() &&
123 full_dims.maxDynamicWellLists() == rst_dims.maxDynamicWellLists();
124 }
125
126 bool operator==(const Welldims& data) const {
127 return this->location() == data.location() &&
128 rst_cmp(*this, data);
129 }
130
131 template<class Serializer>
132 void serializeOp(Serializer& serializer)
133 {
134 serializer(nWMax);
135 serializer(nCWMax);
136 serializer(nWGMax);
137 serializer(nGMax);
138 serializer(nWlistPrWellMax);
139 serializer(nDynWlistMax);
140 serializer(m_location);
141 }
142
143private:
144 int nWMax { 0 };
145 int nCWMax { 0 };
146 int nWGMax { 0 };
147 int nGMax { 0 };
148 int nWlistPrWellMax { 1 };
149 int nDynWlistMax { 1 };
150 std::optional<KeywordLocation> m_location;
151};
152
153class WellSegmentDims {
154public:
155 WellSegmentDims();
156 explicit WellSegmentDims(const Deck& deck);
157
158 static WellSegmentDims serializationTestObject();
159
160 int maxSegmentedWells() const
161 {
162 return this->nSegWellMax;
163 }
164
165 int maxSegmentsPerWell() const
166 {
167 return this->nSegmentMax;
168 }
169
170 int maxLateralBranchesPerWell() const
171 {
172 return this->nLatBranchMax;
173 }
174
175 const std::optional<KeywordLocation>& location() const
176 {
177 return this->location_;
178 }
179
180 bool operator==(const WellSegmentDims& data) const;
181
182 template<class Serializer>
183 void serializeOp(Serializer& serializer)
184 {
185 serializer(nSegWellMax);
186 serializer(nSegmentMax);
187 serializer(nLatBranchMax);
188 serializer(location_);
189 }
190
191private:
192 int nSegWellMax;
193 int nSegmentMax;
194 int nLatBranchMax;
195 std::optional<KeywordLocation> location_;
196};
197
198class NetworkDims {
199public:
200 NetworkDims();
201 explicit NetworkDims(const Deck& deck);
202
203 static NetworkDims serializationTestObject();
204
205 int maxNONodes() const
206 {
207 return this->nMaxNoNodes;
208 }
209
210 int maxNoBranches() const
211 {
212 return this->nMaxNoBranches;
213 }
214
215 int maxNoBranchesConToNode() const
216 {
217 return this->nMaxNoBranchesConToNode;
218 }
219
220 bool extendedNetwork() const
221 {
222 return this->type_ == Type::Extended;
223 }
224
225 bool standardNetwork() const
226 {
227 return this->type_ == Type::Standard;
228 }
229
230 bool active() const
231 {
232 return this->extendedNetwork()
233 || this->standardNetwork();
234 }
235
236 bool operator==(const NetworkDims& data) const;
237
238 template<class Serializer>
239 void serializeOp(Serializer& serializer)
240 {
241 serializer(nMaxNoNodes);
242 serializer(nMaxNoBranches);
243 serializer(nMaxNoBranchesConToNode);
244 }
245
246private:
247 enum class Type { None, Extended, Standard, };
248
249 int nMaxNoNodes;
250 int nMaxNoBranches;
251 int nMaxNoBranchesConToNode;
252 Type type_{ Type::None };
253};
254
255class AquiferDimensions {
256public:
257 AquiferDimensions();
258 explicit AquiferDimensions(const Deck& deck);
259
260 static AquiferDimensions serializationTestObject();
261
262 int maxAnalyticAquifers() const
263 {
264 return this->maxNumAnalyticAquifers;
265 }
266
267 int maxAnalyticAquiferConnections() const
268 {
269 return this->maxNumAnalyticAquiferConn;
270 }
271
272 template <class Serializer>
273 void serializeOp(Serializer& serializer)
274 {
275 serializer(this->maxNumAnalyticAquifers);
276 serializer(this->maxNumAnalyticAquiferConn);
277 }
278
279private:
280 int maxNumAnalyticAquifers;
281 int maxNumAnalyticAquiferConn;
282};
283
284bool operator==(const AquiferDimensions& lhs, const AquiferDimensions& rhs);
285
286class EclHysterConfig
287{
288public:
289 EclHysterConfig() = default;
290 explicit EclHysterConfig(const Deck& deck);
291
292 static EclHysterConfig serializationTestObject();
293
297 //void setActive(bool yesno);
298
302 bool active() const;
303
310 int pcHysteresisModel() const;
311
318 int krHysteresisModel() const;
319
325 double modParamTrapped() const;
326
332 double curvatureCapPrs() const;
333
337 bool activeWag() const;
338
342 bool doPcScaling() const;
343
344 bool operator==(const EclHysterConfig& data) const;
345
346 template<class Serializer>
347 void serializeOp(Serializer& serializer)
348 {
349 serializer(activeHyst);
350 serializer(pcHystMod);
351 serializer(krHystMod);
352 serializer(modParamTrappedValue);
353 serializer(curvatureCapPrsValue);
354 serializer(activeWagHyst);
355 serializer(enablePcScaling);
356 }
357
358private:
359 // enable hysteresis at all
360 bool activeHyst { false };
361
362 // the capillary pressure and the relperm hysteresis models to be used
363 int pcHystMod { -1 };
364 int krHystMod { -1 };
365 // regularisation parameter used for Killough model
366 double modParamTrappedValue { 0.1 };
367 // curvature parameter for capillary pressure
368 double curvatureCapPrsValue { 0.1 };
369
370 // enable WAG hysteresis
371 bool activeWagHyst { false };
372
373 // flag to enable Pc scaling
374 bool enablePcScaling { false };
375};
376
377class SatFuncControls {
378public:
379 enum class ThreePhaseOilKrModel {
380 Default,
381 Stone1,
382 Stone2
383 };
384
385 enum class KeywordFamily {
386 Family_I, // SGOF, SWOF, SLGOF
387 Family_II, // SGFN, SOF{2,3}, SWFN, SGWFN
388 Family_III, // GSF, WSF
389
390 Undefined,
391 };
392
393 SatFuncControls();
394 explicit SatFuncControls(const Deck& deck);
395 explicit SatFuncControls(const double tolcritArg,
396 const ThreePhaseOilKrModel model,
397 const KeywordFamily family);
398
399 static SatFuncControls serializationTestObject();
400
401 double minimumRelpermMobilityThreshold() const
402 {
403 return this->tolcrit;
404 }
405
406 ThreePhaseOilKrModel krModel() const
407 {
408 return this->krmodel;
409 }
410
411 KeywordFamily family() const
412 {
413 return this->satfunc_family;
414 }
415
416 bool operator==(const SatFuncControls& rhs) const;
417
418 template<class Serializer>
419 void serializeOp(Serializer& serializer)
420 {
421 serializer(tolcrit);
422 serializer(krmodel);
423 serializer(satfunc_family);
424 }
425
426private:
427 double tolcrit;
428 ThreePhaseOilKrModel krmodel = ThreePhaseOilKrModel::Default;
429 KeywordFamily satfunc_family = KeywordFamily::Undefined;
430};
431
432
433class Nupcol {
434public:
435 Nupcol();
436 explicit Nupcol(int min_value);
437 void update(int value);
438 int value() const;
439
440 static Nupcol serializationTestObject();
441 bool operator==(const Nupcol& data) const;
442
443 template<class Serializer>
444 void serializeOp(Serializer& serializer) {
445 serializer(this->nupcol_value);
446 serializer(this->min_nupcol);
447 }
448
449private:
450 int min_nupcol;
451 int nupcol_value;
452};
453
454
455class Tracers {
456public:
457 Tracers() = default;
458
459 explicit Tracers(const Deck& );
460 int water_tracers() const;
461
462 template<class Serializer>
463 void serializeOp(Serializer& serializer) {
464 serializer(this->m_oil_tracers);
465 serializer(this->m_water_tracers);
466 serializer(this->m_gas_tracers);
467 serializer(this->m_env_tracers);
468 serializer(this->diffusion_control);
469 serializer(this->max_iter);
470 serializer(this->min_iter);
471 }
472
473 static Tracers serializationTestObject();
474 bool operator==(const Tracers& data) const;
475
476private:
477 int m_oil_tracers{};
478 int m_water_tracers{};
479 int m_gas_tracers{};
480 int m_env_tracers{};
481 bool diffusion_control{false};
482 int max_iter{};
483 int min_iter{};
484 // The TRACERS keyword has some additional options which seem quite arcane,
485 // for now not included here.
486};
487
488
489class Runspec {
490public:
491 Runspec() = default;
492 explicit Runspec( const Deck& );
493
494 static Runspec serializationTestObject();
495
496 std::time_t start_time() const noexcept;
497 const UDQParams& udqParams() const noexcept;
498 const Phases& phases() const noexcept;
499 const Tabdims& tabdims() const noexcept;
500 const Regdims& regdims() const noexcept;
501 const EndpointScaling& endpointScaling() const noexcept;
502 const Welldims& wellDimensions() const noexcept;
503 const WellSegmentDims& wellSegmentDimensions() const noexcept;
504 const NetworkDims& networkDimensions() const noexcept;
505 const AquiferDimensions& aquiferDimensions() const noexcept;
506 int eclPhaseMask( ) const noexcept;
507 const EclHysterConfig& hysterPar() const noexcept;
508 const Actdims& actdims() const noexcept;
509 const SatFuncControls& saturationFunctionControls() const noexcept;
510 const Nupcol& nupcol() const noexcept;
511 const Tracers& tracers() const;
512 bool compositionalMode() const;
513 size_t numComps() const;
514 bool co2Storage() const noexcept;
515 bool co2Sol() const noexcept;
516 bool h2Sol() const noexcept;
517 bool h2Storage() const noexcept;
518 bool micp() const noexcept;
519 bool mech() const noexcept;
520 bool frac() const noexcept;
521 bool temp() const noexcept;
522 bool compositional() const noexcept;
523 bool biof() const noexcept;
524
525 bool operator==(const Runspec& data) const;
526 static bool rst_cmp(const Runspec& full_state, const Runspec& rst_state);
527
528 template<class Serializer>
529 void serializeOp(Serializer& serializer)
530 {
531 serializer(this->m_start_time);
532 serializer(active_phases);
533 serializer(m_tabdims);
534 serializer(m_regdims);
535 serializer(endscale);
536 serializer(welldims);
537 serializer(wsegdims);
538 serializer(netwrkdims);
539 serializer(aquiferdims);
540 serializer(udq_params);
541 serializer(hystpar);
542 serializer(m_actdims);
543 serializer(m_sfuncctrl);
544 serializer(m_nupcol);
545 serializer(m_tracers);
546 serializer(m_comps);
547 serializer(m_co2storage);
548 serializer(m_co2sol);
549 serializer(m_h2sol);
550 serializer(m_h2storage);
551 serializer(m_micp);
552 serializer(m_mech);
553 serializer(m_frac);
554 serializer(m_temp);
555 serializer(m_biof);
556 }
557
558private:
559 std::time_t m_start_time{};
560 Phases active_phases{};
561 Tabdims m_tabdims{};
562 Regdims m_regdims{};
563 EndpointScaling endscale{};
564 Welldims welldims{};
565 WellSegmentDims wsegdims{};
566 NetworkDims netwrkdims{};
567 AquiferDimensions aquiferdims{};
568 UDQParams udq_params{};
569 EclHysterConfig hystpar{};
570 Actdims m_actdims{};
571 SatFuncControls m_sfuncctrl{};
572 Nupcol m_nupcol{};
573 Tracers m_tracers{};
574 size_t m_comps = 0;
575 bool m_co2storage{false};
576 bool m_co2sol{false};
577 bool m_h2sol{false};
578 bool m_h2storage{false};
579 bool m_micp{false};
580 bool m_mech{false};
581 bool m_frac{false};
582 bool m_temp{false};
583 bool m_biof{false};
584};
585
586std::size_t declaredMaxRegionID(const Runspec& rspec);
587
588} // namespace Opm
589
590#endif // OPM_RUNSPEC_HPP
Definition Actdims.hpp:30
Definition Runspec.hpp:255
Definition Deck.hpp:46
Definition Runspec.hpp:287
int pcHysteresisModel() const
Return the type of the hysteresis model which is used for capillary pressure.
Definition Runspec.cpp:478
double modParamTrapped() const
Regularisation parameter used for Killough model.
Definition Runspec.cpp:484
double curvatureCapPrs() const
Curvature parameter used for capillary pressure hysteresis.
Definition Runspec.cpp:487
bool activeWag() const
Wag hysteresis.
Definition Runspec.cpp:490
bool doPcScaling() const
Do Pc scaling for scanning curves.
Definition Runspec.cpp:493
bool active() const
Specify whether hysteresis is enabled or not.
Definition Runspec.cpp:475
int krHysteresisModel() const
Return the type of the hysteresis model which is used for relative permeability.
Definition Runspec.cpp:481
Definition EndpointScaling.hpp:28
Definition Runspec.hpp:198
Definition Runspec.hpp:433
Definition Runspec.hpp:46
Definition Regdims.hpp:36
Definition Runspec.hpp:489
Definition Runspec.hpp:377
Class for (de-)serializing.
Definition Serializer.hpp:94
Definition Tabdims.hpp:36
Definition Runspec.hpp:455
Definition UDQParams.hpp:31
Definition Runspec.hpp:153
Definition Runspec.hpp:75
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30