EquilibrationHelpers_impl.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*/
23
24#include <opm/common/TimingMacros.hpp>
25
26#include <opm/common/utility/numeric/RootFinders.hpp>
27
28#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
29#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
30#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
31
33
34#include <fmt/format.h>
35
36namespace Opm {
37namespace EQUIL {
38
39using FluidSystemSimple = BlackOilFluidSystem<double>;
40
41// Adjust oil pressure according to gas saturation and cap pressure
42using SatOnlyFluidState = SimpleModularFluidState<double,
43 /*numPhases=*/3,
44 /*numComponents=*/3,
46 /*storePressure=*/false,
47 /*storeTemperature=*/false,
48 /*storeComposition=*/false,
49 /*storeFugacity=*/false,
50 /*storeSaturation=*/true,
51 /*storeDensity=*/false,
52 /*storeViscosity=*/false,
53 /*storeEnthalpy=*/false>;
54
55namespace Miscibility {
56
57template<class FluidSystem>
58RsVD<FluidSystem>::RsVD(const int pvtRegionIdx,
59 const std::vector<Scalar>& depth,
60 const std::vector<Scalar>& rs)
61 : pvtRegionIdx_(pvtRegionIdx)
62 , rsVsDepth_(depth, rs)
63{
64}
65
66template<class FluidSystem>
69operator()(const Scalar depth,
70 const Scalar press,
71 const Scalar temp,
72 const Scalar satGas) const
73{
74 const auto sat_rs = satRs(press, temp);
75 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
76 return sat_rs;
77 }
78 else {
79 if (rsVsDepth_.xMin() > depth)
80 return std::min(sat_rs, rsVsDepth_.valueAt(0));
81 else if (rsVsDepth_.xMax() < depth)
82 return std::min(sat_rs, rsVsDepth_.valueAt(rsVsDepth_.numSamples() - 1));
83 return std::min(sat_rs, rsVsDepth_.eval(depth, /*extrapolate=*/false));
84 }
85}
86
87template<class FluidSystem>
90satRs(const Scalar press, const Scalar temp) const
91{
92 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
93}
94
95template<class FluidSystem>
96PBVD<FluidSystem>::PBVD(const int pvtRegionIdx,
97 const std::vector<Scalar>& depth,
98 const std::vector<Scalar>& pbub)
99 : pvtRegionIdx_(pvtRegionIdx)
100 , pbubVsDepth_(depth, pbub)
101{
102}
103
104template<class FluidSystem>
107operator()(const Scalar depth,
108 const Scalar cellPress,
109 const Scalar temp,
110 const Scalar satGas) const
111{
112 Scalar press = cellPress;
113 if (satGas <= 0.0) {
114 if (pbubVsDepth_.xMin() > depth)
115 press = pbubVsDepth_.valueAt(0);
116 else if (pbubVsDepth_.xMax() < depth)
117 press = pbubVsDepth_.valueAt(pbubVsDepth_.numSamples() - 1);
118 else
119 press = pbubVsDepth_.eval(depth, /*extrapolate=*/false);
120 }
121 return satRs(std::min(press, cellPress), temp);
122}
123
124template<class FluidSystem>
127satRs(const Scalar press, const Scalar temp) const
128{
129 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
130}
131
132template<class FluidSystem>
133PDVD<FluidSystem>::PDVD(const int pvtRegionIdx,
134 const std::vector<Scalar>& depth,
135 const std::vector<Scalar>& pdew)
136 : pvtRegionIdx_(pvtRegionIdx)
137 , pdewVsDepth_(depth, pdew)
138{
139}
140
141template<class FluidSystem>
142typename PDVD<FluidSystem>::Scalar
144operator()(const Scalar depth,
145 const Scalar cellPress,
146 const Scalar temp,
147 const Scalar satOil) const
148{
149 Scalar press = cellPress;
150 if (satOil <= 0.0) {
151 if (pdewVsDepth_.xMin() > depth)
152 press = pdewVsDepth_.valueAt(0);
153 else if (pdewVsDepth_.xMax() < depth)
154 press = pdewVsDepth_.valueAt(pdewVsDepth_.numSamples() - 1);
155 else
156 press = pdewVsDepth_.eval(depth, /*extrapolate=*/false);
157 }
158 return satRv(std::min(press, cellPress), temp);
159}
160
161template<class FluidSystem>
162typename PDVD<FluidSystem>::Scalar
164satRv(const Scalar press, const Scalar temp) const
165{
166 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
167}
168
169template<class FluidSystem>
170RvVD<FluidSystem>::RvVD(const int pvtRegionIdx,
171 const std::vector<Scalar>& depth,
172 const std::vector<Scalar>& rv)
173 : pvtRegionIdx_(pvtRegionIdx)
174 , rvVsDepth_(depth, rv)
175{
176}
177
178template<class FluidSystem>
179typename RvVD<FluidSystem>::Scalar
181operator()(const Scalar depth,
182 const Scalar press,
183 const Scalar temp,
184 const Scalar satOil) const
185{
186 if (satOil < - std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
187 throw std::logic_error {
188 "Must not pass negative oil saturation"
189 };
190 }
191 const auto sat_rv = satRv(press, temp);
192 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
193 return sat_rv;
194 }
195 else {
196 if (rvVsDepth_.xMin() > depth)
197 return std::min(sat_rv, rvVsDepth_.valueAt(0));
198 else if (rvVsDepth_.xMax() < depth)
199 return std::min(sat_rv, rvVsDepth_.valueAt(rvVsDepth_.numSamples() - 1));
200 return std::min(sat_rv, rvVsDepth_.eval(depth, /*extrapolate=*/false));
201 }
202}
203
204template<class FluidSystem>
205typename RvVD<FluidSystem>::Scalar
207satRv(const Scalar press, const Scalar temp) const
208{
209 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
210}
211
212template<class FluidSystem>
213RvwVD<FluidSystem>::RvwVD(const int pvtRegionIdx,
214 const std::vector<Scalar>& depth,
215 const std::vector<Scalar>& rvw)
216 : pvtRegionIdx_(pvtRegionIdx)
217 , rvwVsDepth_(depth, rvw)
218{
219}
220
221template<class FluidSystem>
222typename RvwVD<FluidSystem>::Scalar
224operator()(const Scalar depth,
225 const Scalar press,
226 const Scalar temp,
227 const Scalar satWat) const
228{
229 if (satWat < - std::sqrt(std::numeric_limits<double>::epsilon())) {
230 throw std::logic_error {
231 "Must not pass negative water saturation"
232 };
233 }
234
235 const auto sat_rvw = satRvw(press, temp);
236 if (satWat > std::sqrt(std::numeric_limits<double>::epsilon())) {
237 return sat_rvw; //saturated Rvw
238 }
239 else {
240 if (rvwVsDepth_.xMin() > depth)
241 return std::min(sat_rvw,rvwVsDepth_.valueAt(0));
242 else if (rvwVsDepth_.xMax() < depth)
243 return std::min(sat_rvw, rvwVsDepth_.valueAt(rvwVsDepth_.numSamples() - 1));
244 return std::min(sat_rvw, rvwVsDepth_.eval(depth, /*extrapolate=*/false));
245 }
246}
247
248template<class FluidSystem>
249typename RvwVD<FluidSystem>::Scalar
251satRvw(const Scalar press, const Scalar temp) const
252{
253 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);
254}
255
256template<class FluidSystem>
258RsSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
259 : pvtRegionIdx_(pvtRegionIdx)
260{
261 rsSatContact_ = satRs(pContact, T_contact);
262}
263
264template<class FluidSystem>
265typename RsSatAtContact<FluidSystem>::Scalar
267operator()(const Scalar /* depth */,
268 const Scalar press,
269 const Scalar temp,
270 const Scalar satGas) const
271{
272 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
273 return satRs(press, temp);
274 }
275 else {
276 return std::min(satRs(press, temp), rsSatContact_);
277 }
278}
279
280template<class FluidSystem>
281typename RsSatAtContact<FluidSystem>::Scalar
283satRs(const Scalar press, const Scalar temp) const
284{
285 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
286}
287
288template<class FluidSystem>
290RvSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
291 : pvtRegionIdx_(pvtRegionIdx)
292{
293 rvSatContact_ = satRv(pContact, T_contact);
294}
295
296template<class FluidSystem>
297typename RvSatAtContact<FluidSystem>::Scalar
299operator()(const Scalar /*depth*/,
300 const Scalar press,
301 const Scalar temp,
302 const Scalar satOil) const
303{
304 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
305 return satRv(press, temp);
306 }
307 else {
308 return std::min(satRv(press, temp), rvSatContact_);
309 }
310}
311
312template<class FluidSystem>
313typename RvSatAtContact<FluidSystem>::Scalar
315satRv(const Scalar press, const Scalar temp) const
316{
317 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);;
318}
319
320template<class FluidSystem>
322RvwSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
323 : pvtRegionIdx_(pvtRegionIdx)
324{
325 rvwSatContact_ = satRvw(pContact, T_contact);
326}
327
328template<class FluidSystem>
329typename RvwSatAtContact<FluidSystem>::Scalar
331operator()(const Scalar /*depth*/,
332 const Scalar press,
333 const Scalar temp,
334 const Scalar satWat) const
335{
336 if (satWat > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
337 return satRvw(press, temp);
338 }
339 else {
340 return std::min(satRvw(press, temp), rvwSatContact_);
341 }
342}
343
344template<class FluidSystem>
345typename RvwSatAtContact<FluidSystem>::Scalar
347satRvw(const Scalar press, const Scalar temp) const
348{
349 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);;
350}
351
352} // namespace Miscibility
353
354template<class Scalar>
355EquilReg<Scalar>::EquilReg(const EquilRecord& rec,
356 std::shared_ptr<Miscibility::RsFunction<Scalar>> rs,
357 std::shared_ptr<Miscibility::RsFunction<Scalar>> rv,
358 std::shared_ptr<Miscibility::RsFunction<Scalar>> rvw,
359 const TabulatedFunction& tempVdTable,
360 const TabulatedFunction& saltVdTable,
361 const int pvtIdx)
362 : rec_ (rec)
363 , rs_ (rs)
364 , rv_ (rv)
365 , rvw_ (rvw)
366 , tempVdTable_ (tempVdTable)
367 , saltVdTable_ (saltVdTable)
368 , pvtIdx_ (pvtIdx)
369{
370}
371
372template<class Scalar>
374{
375 return this->rec_.datumDepth();
376}
377
378template<class Scalar>
380{
381 return this->rec_.datumDepthPressure();
382}
383
384template<class Scalar>
386{
387 return this->rec_.waterOilContactDepth();
388}
389
390template<class Scalar>
392{
393 return this->rec_.waterOilContactCapillaryPressure();
394}
395
396template<class Scalar>
398{
399 return this->rec_.gasOilContactDepth();
400}
401
402template<class Scalar>
404{
405 return this->rec_.gasOilContactCapillaryPressure();
406}
407
408template<class Scalar>
410{
411 return this->rec_.initializationTargetAccuracy();
412}
413
414template<class Scalar>
417{
418 return *this->rs_;
419}
420
421template<class Scalar>
424{
425 return *this->rv_;
426}
427
428template<class Scalar>
431{
432 return *this->rvw_;
433}
434
435template<class Scalar>
436const typename EquilReg<Scalar>::TabulatedFunction&
438{
439 return saltVdTable_;
440}
441
442template<class Scalar>
443const typename EquilReg<Scalar>::TabulatedFunction&
445{
446 return tempVdTable_;
447}
448
449template<class Scalar>
451{
452 return this->pvtIdx_;
453}
454
455template<class FluidSystem, class MaterialLawManager>
457PcEq(const MaterialLawManager& materialLawManager,
458 const int phase,
459 const int cell,
460 const Scalar targetPc)
461 : materialLawManager_(materialLawManager)
462 , phase_(phase)
463 , cell_(cell)
464 , targetPc_(targetPc)
465{
466}
467
468template<class FluidSystem, class MaterialLawManager>
471operator()(Scalar s) const
472{
473 const auto& matParams = materialLawManager_.materialLawParams(cell_);
474 SatOnlyFluidState fluidState;
475 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
476 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
477 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
478 fluidState.setSaturation(phase_, s);
479
480 // This code should only be called for water or gas phase
481 assert(phase_ != FluidSystem::oilPhaseIdx);
482 // The capillaryPressure code only uses the wetting phase saturation
483 if (phase_ == FluidSystem::gasPhaseIdx) {
484 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0 - s);
485 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0 - s);
486 }
487 std::array<Scalar, FluidSystem::numPhases> pc{0.0};
488 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
489 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
490 Scalar sign = (phase_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
491 Scalar pcPhase = pc[FluidSystem::oilPhaseIdx] + sign * pc[phase_];
492 return pcPhase - targetPc_;
493}
494
495template<class FluidSystem, class MaterialLawManager>
497PcEqSum(const MaterialLawManager& materialLawManager,
498 const int phase1,
499 const int phase2,
500 const int cell,
501 const Scalar targetPc)
502 : materialLawManager_(materialLawManager)
503 , phase1_(phase1)
504 , phase2_(phase2)
505 , cell_(cell)
506 , targetPc_(targetPc)
507{
508}
509
510template<class FluidSystem, class MaterialLawManager>
513operator()(Scalar s) const
514{
515 const auto& matParams = materialLawManager_.materialLawParams(cell_);
516 SatOnlyFluidState fluidState;
517 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
518 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
519 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
520 fluidState.setSaturation(phase1_, s);
521 fluidState.setSaturation(phase2_, 1.0 - s);
522
523 std::array<Scalar, FluidSystem::numPhases> pc {0.0};
524
525 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
526 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
527 Scalar sign1 = (phase1_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
528 Scalar pc1 = pc[FluidSystem::oilPhaseIdx] + sign1 * pc[phase1_];
529 Scalar sign2 = (phase2_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
530 Scalar pc2 = pc[FluidSystem::oilPhaseIdx] + sign2 * pc[phase2_];
531 return pc1 + pc2 - targetPc_;
532}
533
534template <class FluidSystem, class MaterialLawManager>
535typename FluidSystem::Scalar
536minSaturations(const MaterialLawManager& materialLawManager,
537 const int phase, const int cell)
538{
539 const auto& scaledDrainageInfo =
540 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
541
542 // Find minimum and maximum saturations.
543 switch(phase) {
544 case FluidSystem::waterPhaseIdx:
545 return scaledDrainageInfo.Swl;
546
547 case FluidSystem::gasPhaseIdx:
548 return scaledDrainageInfo.Sgl;
549
550 case FluidSystem::oilPhaseIdx:
551 throw std::runtime_error("Min saturation not implemented for oil phase.");
552
553 default:
554 throw std::runtime_error("Unknown phaseIdx .");
555 }
556 return -1.0;
557}
558
559template <class FluidSystem, class MaterialLawManager>
560typename FluidSystem::Scalar
561maxSaturations(const MaterialLawManager& materialLawManager,
562 const int phase, const int cell)
563{
564 const auto& scaledDrainageInfo =
565 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
566
567 // Find minimum and maximum saturations.
568 switch(phase) {
569 case FluidSystem::waterPhaseIdx:
570 return scaledDrainageInfo.Swu;
571
572 case FluidSystem::gasPhaseIdx:
573 return scaledDrainageInfo.Sgu;
574
575 case FluidSystem::oilPhaseIdx:
576 throw std::runtime_error("Max saturation not implemented for oil phase.");
577
578 default:
579 throw std::runtime_error("Unknown phaseIdx .");
580 }
581 return -1.0;
582}
583
584template <class FluidSystem, class MaterialLawManager>
585typename FluidSystem::Scalar
586satFromPc(const MaterialLawManager& materialLawManager,
587 const int phase,
588 const int cell,
589 const typename FluidSystem::Scalar targetPc,
590 const bool increasing)
591{
592 using Scalar = typename FluidSystem::Scalar;
593 // Find minimum and maximum saturations.
594 Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
595 Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
596
597 // Create the equation f(s) = pc(s) - targetPc
598 const PcEq<FluidSystem, MaterialLawManager> f(materialLawManager, phase, cell, targetPc);
599 Scalar f0 = f(s0);
600 Scalar f1 = f(s1);
601 if (!std::isfinite(f0 + f1))
602 throw std::logic_error(fmt::format("The capillary pressure values {} and {} are not finite", f0, f1));
603
604 if (f0 <= 0.0)
605 return s0;
606 else if (f1 >= 0.0)
607 return s1;
608
609 const Scalar tol = 1e-10;
610 // should at least converge in 2 times bisection but some safety here:
611 const int maxIter = -2*static_cast<int>(std::log2(tol)) + 10;
612 int usedIterations = -1;
613 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
614 return root;
615}
616
617template<class FluidSystem, class MaterialLawManager>
618typename FluidSystem::Scalar
619satFromSumOfPcs(const MaterialLawManager& materialLawManager,
620 const int phase1,
621 const int phase2,
622 const int cell,
623 const typename FluidSystem::Scalar targetPc)
624{
625 using Scalar = typename FluidSystem::Scalar;
626 // Find minimum and maximum saturations.
627 Scalar s0 = minSaturations<FluidSystem>(materialLawManager, phase1, cell);
628 Scalar s1 = maxSaturations<FluidSystem>(materialLawManager, phase1, cell);
629
630 // Create the equation f(s) = pc1(s) + pc2(1-s) - targetPc
631 const PcEqSum<FluidSystem, MaterialLawManager> f(materialLawManager, phase1, phase2, cell, targetPc);
632 Scalar f0 = f(s0);
633 Scalar f1 = f(s1);
634 if (f0 <= 0.0)
635 return s0;
636 else if (f1 >= 0.0)
637 return s1;
638
639 assert(f0 > 0.0 && f1 < 0.0);
640 const Scalar tol = 1e-10;
641 // should at least converge in 2 times bisection but some safety here:
642 const int maxIter = -2*static_cast<int>(std::log2(tol)) + 10;
643 int usedIterations = -1;
644 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
645 return root;
646}
647
648template<class FluidSystem, class MaterialLawManager>
649typename FluidSystem::Scalar
650satFromDepth(const MaterialLawManager& materialLawManager,
651 const typename FluidSystem::Scalar cellDepth,
652 const typename FluidSystem::Scalar contactDepth,
653 const int phase,
654 const int cell,
655 const bool increasing)
656{
657 using Scalar = typename FluidSystem::Scalar;
658 const Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
659 const Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
660
661 if (cellDepth < contactDepth) {
662 return s0;
663 }
664 else {
665 return s1;
666 }
667}
668
669template<class FluidSystem, class MaterialLawManager>
670bool isConstPc(const MaterialLawManager& materialLawManager,
671 const int phase,
672 const int cell)
673{
674 using Scalar = typename FluidSystem::Scalar;
675 // Create the equation f(s) = pc(s);
676 const PcEq<FluidSystem, MaterialLawManager> f(materialLawManager, phase, cell, 0);
677 const Scalar f0 = f(minSaturations<FluidSystem>(materialLawManager, phase, cell));
678 const Scalar f1 = f(maxSaturations<FluidSystem>(materialLawManager, phase, cell));
679 return std::abs(f0 - f1) < std::numeric_limits<Scalar>::epsilon();
680}
681
682}
683}
Auxiliary routines that to solve the ODEs that emerge from the hydrostatic equilibrium problem.
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: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
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
BlackOilFluidSystem< double > FluidSystemSimple
Definition: EquilibrationHelpers_impl.hpp:39
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
SimpleModularFluidState< double, 3, 3, FluidSystemSimple, false, false, false, false, true, false, false, false > SatOnlyFluidState
Definition: EquilibrationHelpers_impl.hpp:53
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