EquilibrationHelpers.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 3 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*/
29#ifndef OPM_EQUILIBRATION_HELPERS_HPP
30#define OPM_EQUILIBRATION_HELPERS_HPP
31
32#include <opm/material/common/Tabulated1DFunction.hpp>
33
34#include <opm/input/eclipse/EclipseState/InitConfig/Equil.hpp>
35
36#include <memory>
37#include <vector>
38
39/*
40 ---- synopsis of EquilibrationHelpers.hpp ----
41
42 namespace Opm
43 {
44 namespace EQUIL {
45
46 namespace Miscibility {
47 class RsFunction;
48 class NoMixing;
49 template <class FluidSystem>
50 class RsVD;
51 template <class FluidSystem>
52 class RsSatAtContact;
53 }
54
55 class EquilReg;
56
57
58 template <class FluidSystem, class MaterialLawManager>
59 struct PcEq;
60
61 template <class FluidSystem, class MaterialLawManager>
62 double satFromPc(const MaterialLawManager& materialLawManager,
63 const int phase,
64 const int cell,
65 const double targetPc,
66 const bool increasing = false)
67 template <class FluidSystem, class MaterialLawManager>
68 double satFromSumOfPcs(const MaterialLawManager& materialLawManager,
69 const int phase1,
70 const int phase2,
71 const int cell,
72 const double targetPc)
73 } // namespace Equil
74 } // namespace Opm
75
76 ---- end of synopsis of EquilibrationHelpers.hpp ----
77*/
78
79namespace Opm {
80
88namespace EQUIL {
89
94namespace Miscibility {
95
100{
101public:
102 virtual ~RsFunction() = default;
103
119 virtual double operator()(const double depth,
120 const double press,
121 const double temp,
122 const double sat = 0.0) const = 0;
123};
124
125
129class NoMixing : public RsFunction
130{
131public:
132 virtual ~NoMixing() = default;
133
150 double
151 operator()(const double /* depth */,
152 const double /* press */,
153 const double /* temp */,
154 const double /* sat */ = 0.0) const
155 {
156 return 0.0;
157 }
158};
159
160
166template <class FluidSystem>
167class RsVD : public RsFunction
168{
169public:
177 RsVD(const int pvtRegionIdx,
178 const std::vector<double>& depth,
179 const std::vector<double>& rs);
180
181 virtual ~RsVD() = default;
182
198 double operator()(const double depth,
199 const double press,
200 const double temp,
201 const double satGas = 0.0) const;
202
203private:
204 using RsVsDepthFunc = Tabulated1DFunction<double>;
205
206 const int pvtRegionIdx_;
207 RsVsDepthFunc rsVsDepth_;
208
209 double satRs(const double press, const double temp) const;
210};
211
212
218template <class FluidSystem>
219class PBVD : public RsFunction
220{
221public:
229 PBVD(const int pvtRegionIdx,
230 const std::vector<double>& depth,
231 const std::vector<double>& pbub);
232
233 virtual ~PBVD() = default;
234
249 double operator()(const double depth,
250 const double cellPress,
251 const double temp,
252 const double satGas = 0.0) const;
253
254private:
255 using PbubVsDepthFunc = Tabulated1DFunction<double>;
256
257 const int pvtRegionIdx_;
258 PbubVsDepthFunc pbubVsDepth_;
259
260 double satRs(const double press, const double temp) const;
261};
262
263
269template <class FluidSystem>
270class PDVD : public RsFunction
271{
272public:
280 PDVD(const int pvtRegionIdx,
281 const std::vector<double>& depth,
282 const std::vector<double>& pdew);
283
284 virtual ~PDVD() = default;
285
300 double operator()(const double depth,
301 const double cellPress,
302 const double temp,
303 const double satOil = 0.0) const;
304
305private:
306 using PdewVsDepthFunc = Tabulated1DFunction<double>;
307
308 const int pvtRegionIdx_;
309 PdewVsDepthFunc pdewVsDepth_;
310
311 double satRv(const double press, const double temp) const;
312};
313
314
320template <class FluidSystem>
321class RvVD : public RsFunction
322{
323public:
331 RvVD(const int pvtRegionIdx,
332 const std::vector<double>& depth,
333 const std::vector<double>& rv);
334
350 double operator()(const double depth,
351 const double press,
352 const double temp,
353 const double satOil = 0.0) const;
354
355private:
356 using RvVsDepthFunc = Tabulated1DFunction<double>;
357
358 const int pvtRegionIdx_;
359 RvVsDepthFunc rvVsDepth_;
360
361 double satRv(const double press, const double temp) const;
362};
363
364
370template <class FluidSystem>
371class RvwVD : public RsFunction
372{
373public:
381 RvwVD(const int pvtRegionIdx,
382 const std::vector<double>& depth,
383 const std::vector<double>& rvw);
384
400 double operator()(const double depth,
401 const double press,
402 const double temp,
403 const double satWat = 0.0) const;
404
405private:
406 using RvwVsDepthFunc = Tabulated1DFunction<double>;
407
408 const int pvtRegionIdx_;
409 RvwVsDepthFunc rvwVsDepth_;
410
411 double satRvw(const double press, const double temp) const;
412};
413
414
429template <class FluidSystem>
431{
432public:
440 RsSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact);
441
457 double operator()(const double /* depth */,
458 const double press,
459 const double temp,
460 const double satGas = 0.0) const;
461
462private:
463 const int pvtRegionIdx_;
464 double rsSatContact_;
465
466 double satRs(const double press, const double temp) const;
467};
468
469
484template <class FluidSystem>
486{
487public:
495 RvSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact);
496
512 double operator()(const double /*depth*/,
513 const double press,
514 const double temp,
515 const double satOil = 0.0) const;
516
517private:
518 const int pvtRegionIdx_;
519 double rvSatContact_;
520
521 double satRv(const double press, const double temp) const;
522};
523
538template <class FluidSystem>
540{
541public:
549 RvwSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact);
550
566 double operator()(const double /*depth*/,
567 const double press,
568 const double temp,
569 const double satWat = 0.0) const;
570
571private:
572 const int pvtRegionIdx_;
573 double rvwSatContact_;
574
575 double satRvw(const double press, const double temp) const;
576};
577
578} // namespace Miscibility
579
600{
601 using TabulatedFunction = Tabulated1DFunction<double>;
602
603public:
613 EquilReg(const EquilRecord& rec,
614 std::shared_ptr<Miscibility::RsFunction> rs,
615 std::shared_ptr<Miscibility::RsFunction> rv,
616 std::shared_ptr<Miscibility::RsFunction> rvw,
617 const TabulatedFunction& tempVdTable,
618 const TabulatedFunction& saltVdTable,
619 const int pvtIdx);
620
625
630
635
636
640 double datum() const;
641
645 double pressure() const;
646
650 double zwoc() const;
651
657 double pcowWoc() const;
658
662 double zgoc() const;
663
669 double pcgoGoc() const;
670
678 int equilibrationAccuracy() const;
679
685
691
697
698 const TabulatedFunction& saltVdTable() const;
699 const TabulatedFunction& tempVdTable() const;
700
704 int pvtIdx() const;
705
706private:
707 EquilRecord rec_;
708 std::shared_ptr<Miscibility::RsFunction> rs_;
709 std::shared_ptr<Miscibility::RsFunction> rv_;
710 std::shared_ptr<Miscibility::RsFunction> rvw_;
711 const TabulatedFunction& tempVdTable_;
712 const TabulatedFunction& saltVdTable_;
713 const int pvtIdx_;
714};
715
716
717
721template <class FluidSystem, class MaterialLawManager>
722struct PcEq
723{
724 PcEq(const MaterialLawManager& materialLawManager,
725 const int phase,
726 const int cell,
727 const double targetPc);
728
729 double operator()(double s) const;
730
731private:
732 const MaterialLawManager& materialLawManager_;
733 const int phase_;
734 const int cell_;
735 const double targetPc_;
736};
737
738template <class FluidSystem, class MaterialLawManager>
739double minSaturations(const MaterialLawManager& materialLawManager,
740 const int phase, const int cell);
741
742template <class FluidSystem, class MaterialLawManager>
743double maxSaturations(const MaterialLawManager& materialLawManager,
744 const int phase, const int cell);
745
748template <class FluidSystem, class MaterialLawManager>
749double satFromPc(const MaterialLawManager& materialLawManager,
750 const int phase,
751 const int cell,
752 const double targetPc,
753 const bool increasing = false);
754
758template <class FluidSystem, class MaterialLawManager>
760{
761 PcEqSum(const MaterialLawManager& materialLawManager,
762 const int phase1,
763 const int phase2,
764 const int cell,
765 const double targetPc);
766
767 double operator()(double s) const;
768
769private:
770 const MaterialLawManager& materialLawManager_;
771 const int phase1_;
772 const int phase2_;
773 const int cell_;
774 const double targetPc_;
775};
776
780template <class FluidSystem, class MaterialLawManager>
781double satFromSumOfPcs(const MaterialLawManager& materialLawManager,
782 const int phase1,
783 const int phase2,
784 const int cell,
785 const double targetPc);
786
788template <class FluidSystem, class MaterialLawManager>
789double satFromDepth(const MaterialLawManager& materialLawManager,
790 const double cellDepth,
791 const double contactDepth,
792 const int phase,
793 const int cell,
794 const bool increasing = false);
795
797template <class FluidSystem, class MaterialLawManager>
798bool isConstPc(const MaterialLawManager& materialLawManager,
799 const int phase,
800 const int cell);
801
802} // namespace Equil
803} // namespace Opm
804
805#endif // OPM_EQUILIBRATION_HELPERS_HPP
Definition: EquilibrationHelpers.hpp:600
const CalcDissolution & dissolutionCalculator() const
Definition: EquilibrationHelpers_impl.hpp:393
double datum() const
Definition: EquilibrationHelpers_impl.hpp:357
double pcgoGoc() const
Definition: EquilibrationHelpers_impl.hpp:382
int equilibrationAccuracy() const
Definition: EquilibrationHelpers_impl.hpp:387
double zgoc() const
Definition: EquilibrationHelpers_impl.hpp:377
const CalcWaterEvaporation & waterEvaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:405
const TabulatedFunction & tempVdTable() const
Definition: EquilibrationHelpers_impl.hpp:417
const CalcEvaporation & evaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:399
double pcowWoc() const
Definition: EquilibrationHelpers_impl.hpp:372
double zwoc() const
Definition: EquilibrationHelpers_impl.hpp:367
double pressure() const
Definition: EquilibrationHelpers_impl.hpp:362
int pvtIdx() const
Definition: EquilibrationHelpers_impl.hpp:422
EquilReg(const EquilRecord &rec, std::shared_ptr< Miscibility::RsFunction > rs, std::shared_ptr< Miscibility::RsFunction > rv, std::shared_ptr< Miscibility::RsFunction > rvw, const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtIdx)
Definition: EquilibrationHelpers_impl.hpp:340
const TabulatedFunction & saltVdTable() const
Definition: EquilibrationHelpers_impl.hpp:411
Definition: EquilibrationHelpers.hpp:130
double operator()(const double, const double, const double, const double=0.0) const
Definition: EquilibrationHelpers.hpp:151
Definition: EquilibrationHelpers.hpp:220
double operator()(const double depth, const double cellPress, const double temp, const double satGas=0.0) const
Definition: EquilibrationHelpers_impl.hpp:103
PBVD(const int pvtRegionIdx, const std::vector< double > &depth, const std::vector< double > &pbub)
Definition: EquilibrationHelpers_impl.hpp:93
Definition: EquilibrationHelpers.hpp:271
PDVD(const int pvtRegionIdx, const std::vector< double > &depth, const std::vector< double > &pdew)
Definition: EquilibrationHelpers_impl.hpp:128
double operator()(const double depth, const double cellPress, const double temp, const double satOil=0.0) const
Definition: EquilibrationHelpers_impl.hpp:138
Definition: EquilibrationHelpers.hpp:100
virtual double operator()(const double depth, const double press, const double temp, const double sat=0.0) const =0
Definition: EquilibrationHelpers.hpp:431
RsSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact)
Definition: EquilibrationHelpers_impl.hpp:250
double operator()(const double, const double press, const double temp, const double satGas=0.0) const
Definition: EquilibrationHelpers_impl.hpp:258
Definition: EquilibrationHelpers.hpp:168
RsVD(const int pvtRegionIdx, const std::vector< double > &depth, const std::vector< double > &rs)
Definition: EquilibrationHelpers_impl.hpp:58
double operator()(const double depth, const double press, const double temp, const double satGas=0.0) const
Definition: EquilibrationHelpers_impl.hpp:68
Definition: EquilibrationHelpers.hpp:486
RvSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact)
Definition: EquilibrationHelpers_impl.hpp:280
double operator()(const double, const double press, const double temp, const double satOil=0.0) const
Definition: EquilibrationHelpers_impl.hpp:288
Definition: EquilibrationHelpers.hpp:322
double operator()(const double depth, const double press, const double temp, const double satOil=0.0) const
Definition: EquilibrationHelpers_impl.hpp:174
RvVD(const int pvtRegionIdx, const std::vector< double > &depth, const std::vector< double > &rv)
Definition: EquilibrationHelpers_impl.hpp:164
Definition: EquilibrationHelpers.hpp:540
RvwSatAtContact(const int pvtRegionIdx, const double pContact, const double T_contact)
Definition: EquilibrationHelpers_impl.hpp:310
double operator()(const double, const double press, const double temp, const double satWat=0.0) const
Definition: EquilibrationHelpers_impl.hpp:318
Definition: EquilibrationHelpers.hpp:372
double operator()(const double depth, const double press, const double temp, const double satWat=0.0) const
Definition: EquilibrationHelpers_impl.hpp:216
RvwVD(const int pvtRegionIdx, const std::vector< double > &depth, const std::vector< double > &rvw)
Definition: EquilibrationHelpers_impl.hpp:206
double satFromDepth(const MaterialLawManager &materialLawManager, const double cellDepth, const double contactDepth, const int phase, const int cell, const bool increasing=false)
Compute saturation from depth. Used for constant capillary pressure function.
Definition: EquilibrationHelpers_impl.hpp:613
double satFromSumOfPcs(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const double targetPc)
Definition: EquilibrationHelpers_impl.hpp:584
double satFromPc(const MaterialLawManager &materialLawManager, const int phase, const int cell, const double targetPc, const bool increasing=false)
Definition: EquilibrationHelpers_impl.hpp:553
bool isConstPc(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Return true if capillary pressure function is constant.
Definition: EquilibrationHelpers_impl.hpp:632
double maxSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:529
double minSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:505
Definition: BlackoilPhases.hpp:27
Definition: EquilibrationHelpers.hpp:723
double operator()(double s) const
Definition: EquilibrationHelpers_impl.hpp:442
PcEq(const MaterialLawManager &materialLawManager, const int phase, const int cell, const double targetPc)
Definition: EquilibrationHelpers_impl.hpp:429
Definition: EquilibrationHelpers.hpp:760
PcEqSum(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const double targetPc)
Definition: EquilibrationHelpers_impl.hpp:468
double operator()(double s) const
Definition: EquilibrationHelpers_impl.hpp:483