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
99template<class Scalar>
101{
102public:
103 virtual ~RsFunction() = default;
104
120 virtual Scalar operator()(const Scalar depth,
121 const Scalar press,
122 const Scalar temp,
123 const Scalar sat = 0.0) const = 0;
124};
125
126
130template<class Scalar>
131class NoMixing : public RsFunction<Scalar>
132{
133public:
134 virtual ~NoMixing() = default;
135
152 Scalar
153 operator()(const Scalar /* depth */,
154 const Scalar /* press */,
155 const Scalar /* temp */,
156 const Scalar /* sat */ = 0.0) const override
157 {
158 return 0.0;
159 }
160};
161
162
168template <class FluidSystem>
169class RsVD : public RsFunction<typename FluidSystem::Scalar>
170{
171public:
172 using Scalar = typename FluidSystem::Scalar;
180 RsVD(const int pvtRegionIdx,
181 const std::vector<Scalar>& depth,
182 const std::vector<Scalar>& rs);
183
184 virtual ~RsVD() = default;
185
201 Scalar operator()(const Scalar depth,
202 const Scalar press,
203 const Scalar temp,
204 const Scalar satGas = 0.0) const override;
205
206private:
207 using RsVsDepthFunc = Tabulated1DFunction<Scalar>;
208
209 const int pvtRegionIdx_;
210 RsVsDepthFunc rsVsDepth_;
211
212 Scalar satRs(const Scalar press, const Scalar temp) const;
213};
214
215
221template <class FluidSystem>
222class PBVD : public RsFunction<typename FluidSystem::Scalar>
223{
224public:
225 using Scalar = typename FluidSystem::Scalar;
233 PBVD(const int pvtRegionIdx,
234 const std::vector<Scalar>& depth,
235 const std::vector<Scalar>& pbub);
236
237 virtual ~PBVD() = default;
238
253 Scalar operator()(const Scalar depth,
254 const Scalar cellPress,
255 const Scalar temp,
256 const Scalar satGas = 0.0) const override;
257
258private:
259 using PbubVsDepthFunc = Tabulated1DFunction<Scalar>;
260
261 const int pvtRegionIdx_;
262 PbubVsDepthFunc pbubVsDepth_;
263
264 Scalar satRs(const Scalar press, const Scalar temp) const;
265};
266
267
273template <class FluidSystem>
274class PDVD : public RsFunction<typename FluidSystem::Scalar>
275{
276 using Scalar = typename FluidSystem::Scalar;
277public:
285 PDVD(const int pvtRegionIdx,
286 const std::vector<Scalar>& depth,
287 const std::vector<Scalar>& pdew);
288
289 virtual ~PDVD() = default;
290
305 Scalar operator()(const Scalar depth,
306 const Scalar cellPress,
307 const Scalar temp,
308 const Scalar satOil = 0.0) const override;
309
310private:
311 using PdewVsDepthFunc = Tabulated1DFunction<Scalar>;
312
313 const int pvtRegionIdx_;
314 PdewVsDepthFunc pdewVsDepth_;
315
316 Scalar satRv(const Scalar press, const Scalar temp) const;
317};
318
319
325template <class FluidSystem>
326class RvVD : public RsFunction<typename FluidSystem::Scalar>
327{
328 using Scalar = typename FluidSystem::Scalar;
329public:
337 RvVD(const int pvtRegionIdx,
338 const std::vector<Scalar>& depth,
339 const std::vector<Scalar>& rv);
340
356 Scalar operator()(const Scalar depth,
357 const Scalar press,
358 const Scalar temp,
359 const Scalar satOil = 0.0) const override;
360
361private:
362 using RvVsDepthFunc = Tabulated1DFunction<Scalar>;
363
364 const int pvtRegionIdx_;
365 RvVsDepthFunc rvVsDepth_;
366
367 Scalar satRv(const Scalar press, const Scalar temp) const;
368};
369
370
376template <class FluidSystem>
377class RvwVD : public RsFunction<typename FluidSystem::Scalar>
378{
379 using Scalar = typename FluidSystem::Scalar;
380public:
388 RvwVD(const int pvtRegionIdx,
389 const std::vector<Scalar>& depth,
390 const std::vector<Scalar>& rvw);
391
407 Scalar operator()(const Scalar depth,
408 const Scalar press,
409 const Scalar temp,
410 const Scalar satWat = 0.0) const override;
411
412private:
413 using RvwVsDepthFunc = Tabulated1DFunction<Scalar>;
414
415 const int pvtRegionIdx_;
416 RvwVsDepthFunc rvwVsDepth_;
417
418 Scalar satRvw(const Scalar press, const Scalar temp) const;
419};
420
421
436template <class FluidSystem>
437class RsSatAtContact : public RsFunction<typename FluidSystem::Scalar>
438{
439 using Scalar = typename FluidSystem::Scalar;
440public:
448 RsSatAtContact(const int pvtRegionIdx,
449 const Scalar pContact,
450 const Scalar T_contact);
451
467 Scalar operator()(const Scalar /* depth */,
468 const Scalar press,
469 const Scalar temp,
470 const Scalar satGas = 0.0) const override;
471
472private:
473 const int pvtRegionIdx_;
474 Scalar rsSatContact_;
475
476 Scalar satRs(const Scalar press, const Scalar temp) const;
477};
478
479
494template <class FluidSystem>
495class RvSatAtContact : public RsFunction<typename FluidSystem::Scalar>
496{
497 using Scalar = typename FluidSystem::Scalar;
498public:
506 RvSatAtContact(const int pvtRegionIdx,
507 const Scalar pContact,
508 const Scalar T_contact);
509
525 Scalar operator()(const Scalar /*depth*/,
526 const Scalar press,
527 const Scalar temp,
528 const Scalar satOil = 0.0) const override;
529
530private:
531 const int pvtRegionIdx_;
532 Scalar rvSatContact_;
533
534 Scalar satRv(const Scalar press, const Scalar temp) const;
535};
536
551template <class FluidSystem>
552class RvwSatAtContact : public RsFunction<typename FluidSystem::Scalar>
553{
554 using Scalar = typename FluidSystem::Scalar;
555public:
563 RvwSatAtContact(const int pvtRegionIdx,
564 const Scalar pContact,
565 const Scalar T_contact);
566
582 Scalar operator()(const Scalar /*depth*/,
583 const Scalar press,
584 const Scalar temp,
585 const Scalar satWat = 0.0) const override;
586
587private:
588 const int pvtRegionIdx_;
589 Scalar rvwSatContact_;
590
591 Scalar satRvw(const Scalar press, const Scalar temp) const;
592};
593
594} // namespace Miscibility
595
615template<class Scalar>
617{
618 using TabulatedFunction = Tabulated1DFunction<Scalar>;
619
620public:
630 EquilReg(const EquilRecord& rec,
631 std::shared_ptr<Miscibility::RsFunction<Scalar>> rs,
632 std::shared_ptr<Miscibility::RsFunction<Scalar>> rv,
633 std::shared_ptr<Miscibility::RsFunction<Scalar>> rvw,
634 const TabulatedFunction& tempVdTable,
635 const TabulatedFunction& saltVdTable,
636 const int pvtIdx);
637
642
647
652
653
657 Scalar datum() const;
658
662 Scalar pressure() const;
663
667 Scalar zwoc() const;
668
674 Scalar pcowWoc() const;
675
679 Scalar zgoc() const;
680
686 Scalar pcgoGoc() const;
687
695 int equilibrationAccuracy() const;
696
702
708
714
715 const TabulatedFunction& saltVdTable() const;
716 const TabulatedFunction& tempVdTable() const;
717
721 int pvtIdx() const;
722
723private:
724 EquilRecord rec_;
725 std::shared_ptr<Miscibility::RsFunction<Scalar>> rs_;
726 std::shared_ptr<Miscibility::RsFunction<Scalar>> rv_;
727 std::shared_ptr<Miscibility::RsFunction<Scalar>> rvw_;
728 const TabulatedFunction& tempVdTable_;
729 const TabulatedFunction& saltVdTable_;
730 const int pvtIdx_;
731};
732
733
734
738template <class FluidSystem, class MaterialLawManager>
739struct PcEq
740{
741 using Scalar = typename FluidSystem::Scalar;
742 PcEq(const MaterialLawManager& materialLawManager,
743 const int phase,
744 const int cell,
745 const Scalar targetPc);
746
747 Scalar operator()(Scalar s) const;
748
749private:
750 const MaterialLawManager& materialLawManager_;
751 const int phase_;
752 const int cell_;
753 const Scalar targetPc_;
754};
755
756template <class FluidSystem, class MaterialLawManager>
757typename FluidSystem::Scalar
758minSaturations(const MaterialLawManager& materialLawManager,
759 const int phase, const int cell);
760
761template <class FluidSystem, class MaterialLawManager>
762typename FluidSystem::Scalar
763maxSaturations(const MaterialLawManager& materialLawManager,
764 const int phase, const int cell);
765
768template <class FluidSystem, class MaterialLawManager>
769typename FluidSystem::Scalar
770satFromPc(const MaterialLawManager& materialLawManager,
771 const int phase,
772 const int cell,
773 const typename FluidSystem::Scalar targetPc,
774 const bool increasing = false);
775
779template <class FluidSystem, class MaterialLawManager>
781{
782 using Scalar = typename FluidSystem::Scalar;
783 PcEqSum(const MaterialLawManager& materialLawManager,
784 const int phase1,
785 const int phase2,
786 const int cell,
787 const Scalar targetPc);
788
789 Scalar operator()(Scalar s) const;
790
791private:
792 const MaterialLawManager& materialLawManager_;
793 const int phase1_;
794 const int phase2_;
795 const int cell_;
796 const Scalar targetPc_;
797};
798
802template <class FluidSystem, class MaterialLawManager>
803typename FluidSystem::Scalar
804satFromSumOfPcs(const MaterialLawManager& materialLawManager,
805 const int phase1,
806 const int phase2,
807 const int cell,
808 const typename FluidSystem::Scalar targetPc);
809
811template <class FluidSystem, class MaterialLawManager>
812typename FluidSystem::Scalar
813satFromDepth(const MaterialLawManager& materialLawManager,
814 const typename FluidSystem::Scalar cellDepth,
815 const typename FluidSystem::Scalar contactDepth,
816 const int phase,
817 const int cell,
818 const bool increasing = false);
819
821template <class FluidSystem, class MaterialLawManager>
822bool isConstPc(const MaterialLawManager& materialLawManager,
823 const int phase,
824 const int cell);
825
826} // namespace Equil
827} // namespace Opm
828
829#endif // OPM_EQUILIBRATION_HELPERS_HPP
Definition: EquilibrationHelpers.hpp:617
Scalar datum() const
Definition: EquilibrationHelpers_impl.hpp:373
int pvtIdx() const
Definition: EquilibrationHelpers_impl.hpp:450
EquilReg(const EquilRecord &rec, std::shared_ptr< Miscibility::RsFunction< Scalar > > rs, std::shared_ptr< Miscibility::RsFunction< Scalar > > rv, std::shared_ptr< Miscibility::RsFunction< Scalar > > rvw, const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtIdx)
Definition: EquilibrationHelpers_impl.hpp:355
Scalar pcgoGoc() const
Definition: EquilibrationHelpers_impl.hpp:403
int equilibrationAccuracy() const
Definition: EquilibrationHelpers_impl.hpp:409
const CalcEvaporation & evaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:423
Scalar zgoc() const
Definition: EquilibrationHelpers_impl.hpp:397
Scalar pressure() const
Definition: EquilibrationHelpers_impl.hpp:379
const TabulatedFunction & saltVdTable() const
Definition: EquilibrationHelpers_impl.hpp:437
Scalar pcowWoc() const
Definition: EquilibrationHelpers_impl.hpp:391
const TabulatedFunction & tempVdTable() const
Definition: EquilibrationHelpers_impl.hpp:444
Scalar zwoc() const
Definition: EquilibrationHelpers_impl.hpp:385
const CalcDissolution & dissolutionCalculator() const
Definition: EquilibrationHelpers_impl.hpp:416
const CalcWaterEvaporation & waterEvaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:430
Definition: EquilibrationHelpers.hpp:132
Scalar operator()(const Scalar, const Scalar, const Scalar, const Scalar=0.0) const override
Definition: EquilibrationHelpers.hpp:153
Definition: EquilibrationHelpers.hpp:223
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:107
PBVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pbub)
Definition: EquilibrationHelpers_impl.hpp:96
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:225
Definition: EquilibrationHelpers.hpp:275
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:144
PDVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pdew)
Definition: EquilibrationHelpers_impl.hpp:133
Definition: EquilibrationHelpers.hpp:101
virtual Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar sat=0.0) const =0
Definition: EquilibrationHelpers.hpp:438
RsSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:258
Scalar operator()(const Scalar, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:267
Definition: EquilibrationHelpers.hpp:170
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:69
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:172
RsVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rs)
Definition: EquilibrationHelpers_impl.hpp:58
Definition: EquilibrationHelpers.hpp:496
Scalar operator()(const Scalar, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:299
RvSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:290
Definition: EquilibrationHelpers.hpp:327
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:181
RvVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rv)
Definition: EquilibrationHelpers_impl.hpp:170
Definition: EquilibrationHelpers.hpp:553
Scalar operator()(const Scalar, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:331
RvwSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:322
Definition: EquilibrationHelpers.hpp:378
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:224
RvwVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rvw)
Definition: EquilibrationHelpers_impl.hpp:213
FluidSystem::Scalar satFromDepth(const MaterialLawManager &materialLawManager, const typename FluidSystem::Scalar cellDepth, const typename FluidSystem::Scalar 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:650
FluidSystem::Scalar satFromPc(const MaterialLawManager &materialLawManager, const int phase, const int cell, const typename FluidSystem::Scalar targetPc, const bool increasing=false)
Definition: EquilibrationHelpers_impl.hpp:586
FluidSystem::Scalar minSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:536
FluidSystem::Scalar maxSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:561
bool isConstPc(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Return true if capillary pressure function is constant.
Definition: EquilibrationHelpers_impl.hpp:670
FluidSystem::Scalar satFromSumOfPcs(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const typename FluidSystem::Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:619
Definition: blackoilboundaryratevector.hh:37
Definition: EquilibrationHelpers.hpp:740
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:471
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:741
PcEq(const MaterialLawManager &materialLawManager, const int phase, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:457
Definition: EquilibrationHelpers.hpp:781
PcEqSum(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:497
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:513
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:782