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 <algorithm>
35#include <array>
36#include <cmath>
37#include <limits>
38#include <stdexcept>
39
40#include <fmt/format.h>
41
42namespace Opm {
43namespace EQUIL {
44
45using FluidSystemSimple = BlackOilFluidSystem<double>;
46
47// Adjust oil pressure according to gas saturation and cap pressure
48using SatOnlyFluidState = SimpleModularFluidState<double,
49 /*numPhases=*/3,
50 /*numComponents=*/3,
52 /*storePressure=*/false,
53 /*storeTemperature=*/false,
54 /*storeComposition=*/false,
55 /*storeFugacity=*/false,
56 /*storeSaturation=*/true,
57 /*storeDensity=*/false,
58 /*storeViscosity=*/false,
59 /*storeEnthalpy=*/false>;
60
61namespace Miscibility {
62
63template<class FluidSystem>
64RsVD<FluidSystem>::RsVD(const int pvtRegionIdx,
65 const std::vector<Scalar>& depth,
66 const std::vector<Scalar>& rs)
67 : pvtRegionIdx_(pvtRegionIdx)
68 , rsVsDepth_(depth, rs)
69{
70}
71
72template<class FluidSystem>
75operator()(const Scalar depth,
76 const Scalar press,
77 const Scalar temp,
78 const Scalar satGas) const
79{
80 const auto sat_rs = satRs(press, temp);
81 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
82 return sat_rs;
83 }
84 else {
85 if (rsVsDepth_.xMin() > depth)
86 return std::min(sat_rs, rsVsDepth_.valueAt(0));
87 else if (rsVsDepth_.xMax() < depth)
88 return std::min(sat_rs, rsVsDepth_.valueAt(rsVsDepth_.numSamples() - 1));
89 return std::min(sat_rs, rsVsDepth_.eval(depth, /*extrapolate=*/false));
90 }
91}
92
93template<class FluidSystem>
96satRs(const Scalar press, const Scalar temp) const
97{
98 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
99}
100
101template<class FluidSystem>
102PBVD<FluidSystem>::PBVD(const int pvtRegionIdx,
103 const std::vector<Scalar>& depth,
104 const std::vector<Scalar>& pbub)
105 : pvtRegionIdx_(pvtRegionIdx)
106 , pbubVsDepth_(depth, pbub)
107{
108}
109
110template<class FluidSystem>
113operator()(const Scalar depth,
114 const Scalar cellPress,
115 const Scalar temp,
116 const Scalar satGas) const
117{
118 Scalar press = cellPress;
119 if (satGas <= 0.0) {
120 if (pbubVsDepth_.xMin() > depth)
121 press = pbubVsDepth_.valueAt(0);
122 else if (pbubVsDepth_.xMax() < depth)
123 press = pbubVsDepth_.valueAt(pbubVsDepth_.numSamples() - 1);
124 else
125 press = pbubVsDepth_.eval(depth, /*extrapolate=*/false);
126 }
127 return satRs(std::min(press, cellPress), temp);
128}
129
130template<class FluidSystem>
133satRs(const Scalar press, const Scalar temp) const
134{
135 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
136}
137
138template<class FluidSystem>
139PDVD<FluidSystem>::PDVD(const int pvtRegionIdx,
140 const std::vector<Scalar>& depth,
141 const std::vector<Scalar>& pdew)
142 : pvtRegionIdx_(pvtRegionIdx)
143 , pdewVsDepth_(depth, pdew)
144{
145}
146
147template<class FluidSystem>
148typename PDVD<FluidSystem>::Scalar
150operator()(const Scalar depth,
151 const Scalar cellPress,
152 const Scalar temp,
153 const Scalar satOil) const
154{
155 Scalar press = cellPress;
156 if (satOil <= 0.0) {
157 if (pdewVsDepth_.xMin() > depth)
158 press = pdewVsDepth_.valueAt(0);
159 else if (pdewVsDepth_.xMax() < depth)
160 press = pdewVsDepth_.valueAt(pdewVsDepth_.numSamples() - 1);
161 else
162 press = pdewVsDepth_.eval(depth, /*extrapolate=*/false);
163 }
164 return satRv(std::min(press, cellPress), temp);
165}
166
167template<class FluidSystem>
168typename PDVD<FluidSystem>::Scalar
170satRv(const Scalar press, const Scalar temp) const
171{
172 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
173}
174
175template<class FluidSystem>
176RvVD<FluidSystem>::RvVD(const int pvtRegionIdx,
177 const std::vector<Scalar>& depth,
178 const std::vector<Scalar>& rv)
179 : pvtRegionIdx_(pvtRegionIdx)
180 , rvVsDepth_(depth, rv)
181{
182}
183
184template<class FluidSystem>
185typename RvVD<FluidSystem>::Scalar
187operator()(const Scalar depth,
188 const Scalar press,
189 const Scalar temp,
190 const Scalar satOil) const
191{
192 if (satOil < - std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
193 throw std::logic_error {
194 "Must not pass negative oil saturation"
195 };
196 }
197 const auto sat_rv = satRv(press, temp);
198 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
199 return sat_rv;
200 }
201 else {
202 if (rvVsDepth_.xMin() > depth)
203 return std::min(sat_rv, rvVsDepth_.valueAt(0));
204 else if (rvVsDepth_.xMax() < depth)
205 return std::min(sat_rv, rvVsDepth_.valueAt(rvVsDepth_.numSamples() - 1));
206 return std::min(sat_rv, rvVsDepth_.eval(depth, /*extrapolate=*/false));
207 }
208}
209
210template<class FluidSystem>
211typename RvVD<FluidSystem>::Scalar
213satRv(const Scalar press, const Scalar temp) const
214{
215 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
216}
217
218template<class FluidSystem>
219RvwVD<FluidSystem>::RvwVD(const int pvtRegionIdx,
220 const std::vector<Scalar>& depth,
221 const std::vector<Scalar>& rvw)
222 : pvtRegionIdx_(pvtRegionIdx)
223 , rvwVsDepth_(depth, rvw)
224{
225}
226
227template<class FluidSystem>
228typename RvwVD<FluidSystem>::Scalar
230operator()(const Scalar depth,
231 const Scalar press,
232 const Scalar temp,
233 const Scalar satWat) const
234{
235 if (satWat < - std::sqrt(std::numeric_limits<double>::epsilon())) {
236 throw std::logic_error {
237 "Must not pass negative water saturation"
238 };
239 }
240
241 const auto sat_rvw = satRvw(press, temp);
242 if (satWat > std::sqrt(std::numeric_limits<double>::epsilon())) {
243 return sat_rvw; //saturated Rvw
244 }
245 else {
246 if (rvwVsDepth_.xMin() > depth)
247 return std::min(sat_rvw,rvwVsDepth_.valueAt(0));
248 else if (rvwVsDepth_.xMax() < depth)
249 return std::min(sat_rvw, rvwVsDepth_.valueAt(rvwVsDepth_.numSamples() - 1));
250 return std::min(sat_rvw, rvwVsDepth_.eval(depth, /*extrapolate=*/false));
251 }
252}
253
254template<class FluidSystem>
255typename RvwVD<FluidSystem>::Scalar
257satRvw(const Scalar press, const Scalar temp) const
258{
259 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);
260}
261
262template<class FluidSystem>
264RsSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
265 : pvtRegionIdx_(pvtRegionIdx)
266{
267 rsSatContact_ = satRs(pContact, T_contact);
268}
269
270template<class FluidSystem>
271typename RsSatAtContact<FluidSystem>::Scalar
273operator()(const Scalar /* depth */,
274 const Scalar press,
275 const Scalar temp,
276 const Scalar satGas) const
277{
278 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
279 return satRs(press, temp);
280 }
281 else {
282 return std::min(satRs(press, temp), rsSatContact_);
283 }
284}
285
286template<class FluidSystem>
287typename RsSatAtContact<FluidSystem>::Scalar
289satRs(const Scalar press, const Scalar temp) const
290{
291 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
292}
293
294template<class FluidSystem>
296RvSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
297 : pvtRegionIdx_(pvtRegionIdx)
298{
299 rvSatContact_ = satRv(pContact, T_contact);
300}
301
302template<class FluidSystem>
303typename RvSatAtContact<FluidSystem>::Scalar
305operator()(const Scalar /*depth*/,
306 const Scalar press,
307 const Scalar temp,
308 const Scalar satOil) const
309{
310 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
311 return satRv(press, temp);
312 }
313 else {
314 return std::min(satRv(press, temp), rvSatContact_);
315 }
316}
317
318template<class FluidSystem>
319typename RvSatAtContact<FluidSystem>::Scalar
321satRv(const Scalar press, const Scalar temp) const
322{
323 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);;
324}
325
326template<class FluidSystem>
328RvwSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
329 : pvtRegionIdx_(pvtRegionIdx)
330{
331 rvwSatContact_ = satRvw(pContact, T_contact);
332}
333
334template<class FluidSystem>
335typename RvwSatAtContact<FluidSystem>::Scalar
337operator()(const Scalar /*depth*/,
338 const Scalar press,
339 const Scalar temp,
340 const Scalar satWat) const
341{
342 if (satWat > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
343 return satRvw(press, temp);
344 }
345 else {
346 return std::min(satRvw(press, temp), rvwSatContact_);
347 }
348}
349
350template<class FluidSystem>
351typename RvwSatAtContact<FluidSystem>::Scalar
353satRvw(const Scalar press, const Scalar temp) const
354{
355 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);;
356}
357
358template <class FluidSystem>
360 const Scalar pb_constant)
361 : rs_constant_(rs_constant)
362 , pb_constant_(pb_constant)
363{
364 // Validate that rs_constant is reasonable
365 if (this->rs_constant_ < 0.0) {
366 throw std::invalid_argument("RSCONST RS value cannot be negative");
367 }
368
369 if (this->pb_constant_ <= 0.0) {
370 throw std::invalid_argument("RSCONST bubble point pressure must be positive");
371 }
372}
373
374template <class FluidSystem>
375typename FluidSystem::Scalar
376RsConst<FluidSystem>::satRs([[maybe_unused]] const Scalar press,
377 [[maybe_unused]] const Scalar temp) const
378{
379 // For RSCONST, the saturated Rs is the constant value.
380 return this->rs_constant_;
381}
382
383template <class FluidSystem>
384typename FluidSystem::Scalar
385RsConst<FluidSystem>::operator()([[maybe_unused]] const Scalar depth,
386 [[maybe_unused]] const Scalar press,
387 [[maybe_unused]] const Scalar temp,
388 [[maybe_unused]] const Scalar satGas) const
389{
390 // For RSCONST with PVDO (dead oil):
391 // 1. Oil is always treated as undersaturated with constant dissolved gas.
392 // 2. Simulation should terminate if pressure < pb_constant_.
393 // 3. Always ensure we use undersaturated PVT calculation.
394
395 if (press < this->pb_constant_) {
396 throw std::invalid_argument {
397 fmt::format("Pressure {} is below bubble point pressure {}",
398 press, this->pb_constant_)
399 };
400 }
401
402 return this->rs_constant_;
403}
404
405} // namespace Miscibility
406
407template<class Scalar>
408EquilReg<Scalar>::EquilReg(const EquilRecord& rec,
409 std::shared_ptr<Miscibility::RsFunction<Scalar>> rs,
410 std::shared_ptr<Miscibility::RsFunction<Scalar>> rv,
411 std::shared_ptr<Miscibility::RsFunction<Scalar>> rvw,
412 const TabulatedFunction& tempVdTable,
413 const TabulatedFunction& saltVdTable,
414 const int pvtIdx)
415 : rec_ (rec)
416 , rs_ (rs)
417 , rv_ (rv)
418 , rvw_ (rvw)
419 , tempVdTable_ (tempVdTable)
420 , saltVdTable_ (saltVdTable)
421 , pvtIdx_ (pvtIdx)
422{
423}
424
425template<class Scalar>
427{
428 return this->rec_.datumDepth();
429}
430
431template<class Scalar>
433{
434 return this->rec_.datumDepthPressure();
435}
436
437template<class Scalar>
439{
440 return this->rec_.waterOilContactDepth();
441}
442
443template<class Scalar>
445{
446 return this->rec_.waterOilContactCapillaryPressure();
447}
448
449template<class Scalar>
451{
452 return this->rec_.gasOilContactDepth();
453}
454
455template<class Scalar>
457{
458 return this->rec_.gasOilContactCapillaryPressure();
459}
460
461template<class Scalar>
463{
464 return this->rec_.initializationTargetAccuracy();
465}
466
467template<class Scalar>
470{
471 return *this->rs_;
472}
473
474template<class Scalar>
477{
478 return *this->rv_;
479}
480
481template<class Scalar>
484{
485 return *this->rvw_;
486}
487
488template<class Scalar>
489const typename EquilReg<Scalar>::TabulatedFunction&
491{
492 return saltVdTable_;
493}
494
495template<class Scalar>
496const typename EquilReg<Scalar>::TabulatedFunction&
498{
499 return tempVdTable_;
500}
501
502template<class Scalar>
504{
505 return this->pvtIdx_;
506}
507
508template<class FluidSystem, class MaterialLawManager>
510PcEq(const MaterialLawManager& materialLawManager,
511 const int phase,
512 const int cell,
513 const Scalar targetPc)
514 : materialLawManager_(materialLawManager)
515 , phase_(phase)
516 , cell_(cell)
517 , targetPc_(targetPc)
518{
519}
520
521template<class FluidSystem, class MaterialLawManager>
524operator()(Scalar s) const
525{
526 const auto& matParams = materialLawManager_.materialLawParams(cell_);
527 SatOnlyFluidState fluidState;
528 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
529 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
530 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
531 fluidState.setSaturation(phase_, s);
532
533 // This code should only be called for water or gas phase
534 assert(phase_ != FluidSystem::oilPhaseIdx);
535 // The capillaryPressure code only uses the wetting phase saturation
536 if (phase_ == FluidSystem::gasPhaseIdx) {
537 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0 - s);
538 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0 - s);
539 }
540 std::array<Scalar, FluidSystem::numPhases> pc{0.0};
541 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
542 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
543 Scalar sign = (phase_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
544 Scalar pcPhase = pc[FluidSystem::oilPhaseIdx] + sign * pc[phase_];
545 return pcPhase - targetPc_;
546}
547
548template<class FluidSystem, class MaterialLawManager>
550PcEqSum(const MaterialLawManager& materialLawManager,
551 const int phase1,
552 const int phase2,
553 const int cell,
554 const Scalar targetPc)
555 : materialLawManager_(materialLawManager)
556 , phase1_(phase1)
557 , phase2_(phase2)
558 , cell_(cell)
559 , targetPc_(targetPc)
560{
561}
562
563template<class FluidSystem, class MaterialLawManager>
566operator()(Scalar s) const
567{
568 const auto& matParams = materialLawManager_.materialLawParams(cell_);
569 SatOnlyFluidState fluidState;
570 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
571 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
572 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
573 fluidState.setSaturation(phase1_, s);
574 fluidState.setSaturation(phase2_, 1.0 - s);
575
576 std::array<Scalar, FluidSystem::numPhases> pc {0.0};
577
578 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
579 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
580 Scalar sign1 = (phase1_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
581 Scalar pc1 = pc[FluidSystem::oilPhaseIdx] + sign1 * pc[phase1_];
582 Scalar sign2 = (phase2_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
583 Scalar pc2 = pc[FluidSystem::oilPhaseIdx] + sign2 * pc[phase2_];
584 return pc1 + pc2 - targetPc_;
585}
586
587template <class FluidSystem, class MaterialLawManager>
588typename FluidSystem::Scalar
589minSaturations(const MaterialLawManager& materialLawManager,
590 const int phase, const int cell)
591{
592 const auto& scaledDrainageInfo =
593 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
594
595 // Find minimum and maximum saturations.
596 switch(phase) {
597 case FluidSystem::waterPhaseIdx:
598 return scaledDrainageInfo.Swl;
599
600 case FluidSystem::gasPhaseIdx:
601 return scaledDrainageInfo.Sgl;
602
603 case FluidSystem::oilPhaseIdx:
604 throw std::runtime_error("Min saturation not implemented for oil phase.");
605
606 default:
607 throw std::runtime_error("Unknown phaseIdx .");
608 }
609 return -1.0;
610}
611
612template <class FluidSystem, class MaterialLawManager>
613typename FluidSystem::Scalar
614maxSaturations(const MaterialLawManager& materialLawManager,
615 const int phase, const int cell)
616{
617 const auto& scaledDrainageInfo =
618 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
619
620 // Find minimum and maximum saturations.
621 switch(phase) {
622 case FluidSystem::waterPhaseIdx:
623 return scaledDrainageInfo.Swu;
624
625 case FluidSystem::gasPhaseIdx:
626 return scaledDrainageInfo.Sgu;
627
628 case FluidSystem::oilPhaseIdx:
629 throw std::runtime_error("Max saturation not implemented for oil phase.");
630
631 default:
632 throw std::runtime_error("Unknown phaseIdx .");
633 }
634 return -1.0;
635}
636
637template <class FluidSystem, class MaterialLawManager>
638typename FluidSystem::Scalar
639satFromPc(const MaterialLawManager& materialLawManager,
640 const int phase,
641 const int cell,
642 const typename FluidSystem::Scalar targetPc,
643 const bool increasing)
644{
645 using Scalar = typename FluidSystem::Scalar;
646 // Find minimum and maximum saturations.
647 Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
648 Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
649
650 // Create the equation f(s) = pc(s) - targetPc
651 const PcEq<FluidSystem, MaterialLawManager> f(materialLawManager, phase, cell, targetPc);
652 Scalar f0 = f(s0);
653 Scalar f1 = f(s1);
654 if (!std::isfinite(f0 + f1))
655 throw std::logic_error(fmt::format("The capillary pressure values {} and {} are not finite", f0, f1));
656
657 if (f0 <= 0.0)
658 return s0;
659 else if (f1 >= 0.0)
660 return s1;
661
662 const Scalar tol = 1e-10;
663 // should at least converge in 2 times bisection but some safety here:
664 const int maxIter = -2*static_cast<int>(std::log2(tol)) + 10;
665 int usedIterations = -1;
666 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
667 return root;
668}
669
670template<class FluidSystem, class MaterialLawManager>
671typename FluidSystem::Scalar
672satFromSumOfPcs(const MaterialLawManager& materialLawManager,
673 const int phase1,
674 const int phase2,
675 const int cell,
676 const typename FluidSystem::Scalar targetPc)
677{
678 using Scalar = typename FluidSystem::Scalar;
679 // Find minimum and maximum saturations.
680 Scalar s0 = minSaturations<FluidSystem>(materialLawManager, phase1, cell);
681 Scalar s1 = maxSaturations<FluidSystem>(materialLawManager, phase1, cell);
682
683 // Create the equation f(s) = pc1(s) + pc2(1-s) - targetPc
684 const PcEqSum<FluidSystem, MaterialLawManager> f(materialLawManager, phase1, phase2, cell, targetPc);
685 Scalar f0 = f(s0);
686 Scalar f1 = f(s1);
687 if (f0 <= 0.0)
688 return s0;
689 else if (f1 >= 0.0)
690 return s1;
691
692 assert(f0 > 0.0 && f1 < 0.0);
693 const Scalar tol = 1e-10;
694 // should at least converge in 2 times bisection but some safety here:
695 const int maxIter = -2*static_cast<int>(std::log2(tol)) + 10;
696 int usedIterations = -1;
697 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
698 return root;
699}
700
701template<class FluidSystem, class MaterialLawManager>
702typename FluidSystem::Scalar
703satFromDepth(const MaterialLawManager& materialLawManager,
704 const typename FluidSystem::Scalar cellDepth,
705 const typename FluidSystem::Scalar contactDepth,
706 const int phase,
707 const int cell,
708 const bool increasing)
709{
710 using Scalar = typename FluidSystem::Scalar;
711 const Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
712 const Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
713
714 if (cellDepth <= contactDepth) {
715 return s0;
716 }
717 else {
718 return s1;
719 }
720}
721
722template<class FluidSystem, class MaterialLawManager>
723bool isConstPc(const MaterialLawManager& materialLawManager,
724 const int phase,
725 const int cell)
726{
727 using Scalar = typename FluidSystem::Scalar;
728 // Create the equation f(s) = pc(s);
729 const PcEq<FluidSystem, MaterialLawManager> f(materialLawManager, phase, cell, 0);
730 const Scalar f0 = f(minSaturations<FluidSystem>(materialLawManager, phase, cell));
731 const Scalar f1 = f(maxSaturations<FluidSystem>(materialLawManager, phase, cell));
732 return std::abs(f0 - f1) < std::numeric_limits<Scalar>::epsilon();
733}
734
735}
736}
Auxiliary routines that to solve the ODEs that emerge from the hydrostatic equilibrium problem.
Scalar datum() const
Definition: EquilibrationHelpers_impl.hpp:426
int pvtIdx() const
Definition: EquilibrationHelpers_impl.hpp:503
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:408
Scalar pcgoGoc() const
Definition: EquilibrationHelpers_impl.hpp:456
int equilibrationAccuracy() const
Definition: EquilibrationHelpers_impl.hpp:462
const CalcEvaporation & evaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:476
Scalar zgoc() const
Definition: EquilibrationHelpers_impl.hpp:450
Scalar pressure() const
Definition: EquilibrationHelpers_impl.hpp:432
const TabulatedFunction & saltVdTable() const
Definition: EquilibrationHelpers_impl.hpp:490
Scalar pcowWoc() const
Definition: EquilibrationHelpers_impl.hpp:444
const TabulatedFunction & tempVdTable() const
Definition: EquilibrationHelpers_impl.hpp:497
Scalar zwoc() const
Definition: EquilibrationHelpers_impl.hpp:438
const CalcDissolution & dissolutionCalculator() const
Definition: EquilibrationHelpers_impl.hpp:469
const CalcWaterEvaporation & waterEvaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:483
Definition: EquilibrationHelpers.hpp:216
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:113
PBVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pbub)
Definition: EquilibrationHelpers_impl.hpp:102
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:218
Definition: EquilibrationHelpers.hpp:269
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:150
PDVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pdew)
Definition: EquilibrationHelpers_impl.hpp:139
RsConst(const Scalar rs_constant, const Scalar pb_constant)
Definition: EquilibrationHelpers_impl.hpp:359
Scalar satRs(const Scalar press, const Scalar temp) const
Definition: EquilibrationHelpers_impl.hpp:376
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:385
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:614
Definition: EquilibrationHelpers.hpp:101
Definition: EquilibrationHelpers.hpp:439
RsSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:264
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:273
Definition: EquilibrationHelpers.hpp:162
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:75
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:164
RsVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rs)
Definition: EquilibrationHelpers_impl.hpp:64
Definition: EquilibrationHelpers.hpp:500
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:305
RvSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:296
Definition: EquilibrationHelpers.hpp:322
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:187
RvVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rv)
Definition: EquilibrationHelpers_impl.hpp:176
Definition: EquilibrationHelpers.hpp:560
RvwSatAtContact(const int pvtRegionIdx, const Scalar pContact, const Scalar T_contact)
Definition: EquilibrationHelpers_impl.hpp:328
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:337
Definition: EquilibrationHelpers.hpp:376
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:230
RvwVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rvw)
Definition: EquilibrationHelpers_impl.hpp:219
BlackOilFluidSystem< double > FluidSystemSimple
Definition: EquilibrationHelpers_impl.hpp:45
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:703
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:639
SimpleModularFluidState< double, 3, 3, FluidSystemSimple, false, false, false, false, true, false, false, false > SatOnlyFluidState
Definition: EquilibrationHelpers_impl.hpp:59
FluidSystem::Scalar minSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:589
FluidSystem::Scalar maxSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:614
bool isConstPc(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Return true if capillary pressure function is constant.
Definition: EquilibrationHelpers_impl.hpp:723
FluidSystem::Scalar satFromSumOfPcs(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const typename FluidSystem::Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:672
Definition: blackoilbioeffectsmodules.hh:45
Definition: EquilibrationHelpers.hpp:801
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:524
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:802
PcEq(const MaterialLawManager &materialLawManager, const int phase, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:510
Definition: EquilibrationHelpers.hpp:842
PcEqSum(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:550
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:566
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:843