opm-common
Loading...
Searching...
No Matches
BlackOilFluidSystem_macrotemplate.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 2 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*/
27
33#include "opm/material/fluidsystems/blackoilpvt/NullOilPvt.hpp"
34
35#include <opm/common/ErrorMacros.hpp>
36#include <opm/common/TimingMacros.hpp>
37#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
38#include <opm/common/utility/gpuDecorators.hpp>
40
43
49#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
50
51#include <array>
52#include <cstddef>
53#include <memory>
54#include <stdexcept>
55#include <string>
56#include <string_view>
57#include <type_traits>
58#include <vector>
59
60namespace Opm {
61
62class EclipseState;
63class Schedule;
64
71template <class Scalar, class IndexTraits = BlackOilDefaultFluidSystemIndices,
72 template<typename> typename Storage = VectorWithDefaultAllocator>
73class FLUIDSYSTEM_CLASSNAME : public BaseFluidSystem<Scalar, FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, Storage> >
74{
75 using ThisType = FLUIDSYSTEM_CLASSNAME;
76
77public:
78
79 // The logic here is the following: We test if we are instantiating on the CPU
80 // std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar> == true
81 // and if so, we use the multiplexer classes, if not, we use the concrete types.
82 using GasPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
85 using OilPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
87 NullOilPvt<Scalar>>; // Currently do not support on GPU
88 using WaterPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
91 using IndexTraitsType = IndexTraits;
92
93 #ifdef COMPILING_STATIC_FLUID_SYSTEM
95 template <class EvaluationT>
96 struct ParameterCache : public NullParameterCache<EvaluationT>
97 {
98 using Evaluation = EvaluationT;
99
100 public:
101 explicit ParameterCache(unsigned regionIdx = 0)
102 : regionIdx_(regionIdx)
103 {
104 }
105
113 template <class OtherCache>
114 void assignPersistentData(const OtherCache& other)
115 {
116 regionIdx_ = other.regionIndex();
117 }
118
126 unsigned regionIndex() const
127 { return regionIdx_; }
128
136 void setRegionIndex(unsigned val)
137 { regionIdx_ = val; }
138
139 private:
140 unsigned regionIdx_;
141 };
142
143 #else
144
145 // We want to use the exact same ParameterCache class for both the static and nonstatic versions
146 template<class EvaluationT>
147 using ParameterCache =
148 typename FLUIDSYSTEM_CLASSNAME_STATIC<Scalar,
149 IndexTraits,
150 Storage>::template ParameterCache<EvaluationT>;
151 #endif
152
153 #ifndef COMPILING_STATIC_FLUID_SYSTEM
159 template<template<typename> typename StorageT>
163 , reservoirTemperature_(other.reservoirTemperature_)
164 , gasPvt_(other.gasPvt_)
165 , oilPvt_(other.oilPvt_)
166 , waterPvt_(other.waterPvt_)
167 , enableDissolvedGas_(other.enableDissolvedGas_)
168 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
169 , enableVaporizedOil_(other.enableVaporizedOil_)
170 , enableVaporizedWater_(other.enableVaporizedWater_)
171 , enableDiffusion_(other.enableDiffusion_)
172 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
173 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
174 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
175 , phaseUsageInfo_(other.phaseUsageInfo_)
176 , isInitialized_(other.isInitialized_)
177 , useSaturatedTables_(other.useSaturatedTables_)
178 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
179 {
180 }
181
182 FLUIDSYSTEM_CLASSNAME(Scalar _surfacePressure_,
183 Scalar _surfaceTemperature_,
184 Scalar _reservoirTemperature_,
185 const GasPvt& _gasPvt_,
186 const OilPvt& _oilPvt_,
187 const WaterPvt& _waterPvt_,
188 bool _enableDissolvedGas_,
189 bool _enableDissolvedGasInWater_,
190 bool _enableVaporizedOil_,
191 bool _enableVaporizedWater_,
192 bool _enableConstantRs_,
193 bool _enableDiffusion_,
194 Storage<std::array<Scalar, 3>>&& _referenceDensity_,
195 Storage<std::array<Scalar, 3>>&& _molarMass_,
196 Storage<std::array<Scalar, 3 * 3>>&& _diffusionCoefficients_,
197 const PhaseUsageInfo<IndexTraits>& _phaseUsageInfo_,
198 bool _isInitialized_,
199 bool _useSaturatedTables_,
200 bool _enthalpy_eq_energy_)
201 : surfacePressure(_surfacePressure_)
202 , surfaceTemperature(_surfaceTemperature_)
203 , reservoirTemperature_(_reservoirTemperature_)
204 , gasPvt_(_gasPvt_)
205 , oilPvt_(_oilPvt_)
206 , waterPvt_(_waterPvt_)
207 , enableDissolvedGas_(_enableDissolvedGas_)
208 , enableDissolvedGasInWater_(_enableDissolvedGasInWater_)
209 , enableVaporizedOil_(_enableVaporizedOil_)
210 , enableVaporizedWater_(_enableVaporizedWater_)
211 , enableConstantRs_(_enableConstantRs_)
212 , enableDiffusion_(_enableDiffusion_)
213 , referenceDensity_(std::move(_referenceDensity_))
214 , molarMass_(std::move(_molarMass_))
215 , diffusionCoefficients_(std::move(_diffusionCoefficients_))
216 , phaseUsageInfo_(_phaseUsageInfo_)
217 , isInitialized_(_isInitialized_)
218 , useSaturatedTables_(_useSaturatedTables_)
219 , enthalpy_eq_energy_(_enthalpy_eq_energy_)
220 {
221 }
222
223 #if HAVE_CUDA
224 template <class ScalarT, class IndexTraitsT>
225 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>
226 gpuistl::copy_to_gpu(const FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT>& oldFluidSystem);
227
228 template <class ScalarT, class IndexTraitsT>
229 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuView>
230 gpuistl::make_view(FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>& oldFluidSystem);
231
232 #endif // HAVE_CUDA
233 #endif // COMPILING_STATIC_FLUID_SYSTEM
234
235 #ifdef COMPILING_STATIC_FLUID_SYSTEM
243 template<template<typename> typename StorageT = VectorWithDefaultAllocator>
244 static FLUIDSYSTEM_CLASSNAME_NONSTATIC<Scalar, IndexTraits, StorageT>& getNonStaticInstance()
245 {
246 static FLUIDSYSTEM_CLASSNAME_NONSTATIC<Scalar, IndexTraits, StorageT> instance{FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, Storage>()};
247 return instance;
248
249 }
250 #endif
251
252 /****************************************
253 * Initialization
254 ****************************************/
258 STATIC_OR_NOTHING void initFromState(const EclipseState& eclState, const Schedule& schedule);
259
268 STATIC_OR_NOTHING void initBegin(std::size_t numPvtRegions);
269
276 STATIC_OR_DEVICE void setEnableDissolvedGas(bool yesno)
277 { enableDissolvedGas_ = yesno; }
278
285 STATIC_OR_DEVICE void setEnableVaporizedOil(bool yesno)
286 { enableVaporizedOil_ = yesno; }
287
294 STATIC_OR_DEVICE void setEnableVaporizedWater(bool yesno)
295 { enableVaporizedWater_ = yesno; }
296
303 STATIC_OR_DEVICE void setEnableDissolvedGasInWater(bool yesno)
304 { enableDissolvedGasInWater_ = yesno; }
305
312 STATIC_OR_DEVICE void setEnableConstantRs(bool yesno)
313 { enableConstantRs_ = yesno; }
314
320 STATIC_OR_DEVICE void setEnableDiffusion(bool yesno)
321 { enableDiffusion_ = yesno; }
322
328 STATIC_OR_DEVICE void setUseSaturatedTables(bool yesno)
329 { useSaturatedTables_ = yesno; }
330
334 STATIC_OR_DEVICE void setGasPvt(std::shared_ptr<GasPvt> pvtObj)
335 { gasPvt_ = *pvtObj; }
336
340 STATIC_OR_DEVICE void setOilPvt(std::shared_ptr<OilPvt> pvtObj)
341 { oilPvt_ = *pvtObj; }
342
346 STATIC_OR_DEVICE void setWaterPvt(std::shared_ptr<WaterPvt> pvtObj)
347 { waterPvt_ = *pvtObj; }
348
349 STATIC_OR_DEVICE void setVapPars(const Scalar par1, const Scalar par2)
350 {
351 if (gasPvt_.isActive()) {
352 gasPvt_.setVapPars(par1, par2);
353 }
354 if (oilPvt_.isActive()) {
355 oilPvt_.setVapPars(par1, par2);
356 }
357 if (waterPvt_.isActive()) {
358 waterPvt_.setVapPars(par1, par2);
359 }
360 }
361
370 STATIC_OR_DEVICE void setReferenceDensities(Scalar rhoOil,
371 Scalar rhoWater,
372 Scalar rhoGas,
373 unsigned regionIdx);
374
378 STATIC_OR_DEVICE void initEnd();
379
380 STATIC_OR_DEVICE bool isInitialized() NOTHING_OR_CONST
381 { return isInitialized_; }
382
383 /****************************************
384 * Generic phase properties
385 ****************************************/
386
388 static constexpr unsigned numPhases = IndexTraits::numPhases;
389
391 static constexpr unsigned waterPhaseIdx = IndexTraits::waterPhaseIdx;
393 static constexpr unsigned oilPhaseIdx = IndexTraits::oilPhaseIdx;
395 static constexpr unsigned gasPhaseIdx = IndexTraits::gasPhaseIdx;
396
398 STATIC_OR_NOTHING Scalar surfacePressure;
399
401 STATIC_OR_NOTHING Scalar surfaceTemperature;
402
404 STATIC_OR_NOTHING std::string_view phaseName(unsigned phaseIdx) NOTHING_OR_CONST;
405
407 STATIC_OR_DEVICE bool isLiquid(unsigned phaseIdx) NOTHING_OR_CONST
408 {
409 assert(phaseIdx < numPhases);
410 return phaseIdx != gasPhaseIdx;
411 }
412
413 /****************************************
414 * Generic component related properties
415 ****************************************/
416
418 static constexpr unsigned numComponents = IndexTraits::numComponents;
419
421 static constexpr int oilCompIdx = IndexTraits::oilCompIdx;
423 static constexpr int waterCompIdx = IndexTraits::waterCompIdx;
425 static constexpr int gasCompIdx = IndexTraits::gasCompIdx;
426
427public:
428
430 STATIC_OR_DEVICE const PhaseUsageInfo<IndexTraits>& phaseUsage() NOTHING_OR_CONST
431 { return phaseUsageInfo_; }
432
434 STATIC_OR_DEVICE unsigned numActivePhases() NOTHING_OR_CONST
435 { return phaseUsageInfo_.numActivePhases(); }
436
438 STATIC_OR_DEVICE bool phaseIsActive(unsigned phaseIdx) NOTHING_OR_CONST
439 {
440 return phaseUsageInfo_.phaseIsActive(phaseIdx);
441 }
442
444 STATIC_OR_DEVICE unsigned solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST;
445
447 STATIC_OR_DEVICE unsigned soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST;
448
450 STATIC_OR_NOTHING std::string_view componentName(unsigned compIdx) NOTHING_OR_CONST;
451
453 STATIC_OR_DEVICE Scalar molarMass(unsigned compIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
454 { return molarMass_[regionIdx][compIdx]; }
455
457 STATIC_OR_DEVICE bool isIdealMixture(unsigned /*phaseIdx*/)
458 {
459 // fugacity coefficients are only pressure dependent -> we
460 // have an ideal mixture
461 return true;
462 }
463
465 STATIC_OR_DEVICE bool isCompressible(unsigned /*phaseIdx*/) NOTHING_OR_CONST
466 { return true; /* all phases are compressible */ }
467
469 STATIC_OR_DEVICE bool isIdealGas(unsigned /*phaseIdx*/) NOTHING_OR_CONST
470 { return false; }
471
472
473 /****************************************
474 * Black-oil specific properties
475 ****************************************/
481 STATIC_OR_DEVICE std::size_t numRegions() NOTHING_OR_CONST
482 { return molarMass_.size(); }
483
490 STATIC_OR_DEVICE bool enableDissolvedGas() NOTHING_OR_CONST
491 { return enableDissolvedGas_; }
492
493
500 STATIC_OR_DEVICE bool enableDissolvedGasInWater() NOTHING_OR_CONST
501 { return enableDissolvedGasInWater_; }
502
509 STATIC_OR_DEVICE bool enableVaporizedOil() NOTHING_OR_CONST
510 { return enableVaporizedOil_; }
511
518 STATIC_OR_DEVICE bool enableVaporizedWater() NOTHING_OR_CONST
519 { return enableVaporizedWater_; }
520
524 STATIC_OR_DEVICE bool enableConstantRs() NOTHING_OR_CONST
525 { return enableConstantRs_; }
526
532 STATIC_OR_DEVICE bool enableDiffusion() NOTHING_OR_CONST
533 { return enableDiffusion_; }
534
540 STATIC_OR_DEVICE bool useSaturatedTables() NOTHING_OR_CONST
541 { return useSaturatedTables_; }
542
548 STATIC_OR_DEVICE Scalar referenceDensity(unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
549 { return referenceDensity_[regionIdx][phaseIdx]; }
550
551 /****************************************
552 * thermodynamic quantities (generic version)
553 ****************************************/
555 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
556 STATIC_OR_DEVICE LhsEval density(const FluidState& fluidState,
557 const ParameterCache<ParamCacheEval>& paramCache,
558 unsigned phaseIdx) NOTHING_OR_CONST
559 { return density<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
560
562 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
563 STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState& fluidState,
564 const ParameterCache<ParamCacheEval>& paramCache,
565 unsigned phaseIdx,
566 unsigned compIdx) NOTHING_OR_CONST
567 {
569 phaseIdx,
570 compIdx,
571 paramCache.regionIndex());
572 }
573
575 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
576 STATIC_OR_DEVICE LhsEval viscosity(const FluidState& fluidState,
577 const ParameterCache<ParamCacheEval>& paramCache,
578 unsigned phaseIdx) NOTHING_OR_CONST
579 { return viscosity<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
580
582 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
583 STATIC_OR_DEVICE LhsEval enthalpy(const FluidState& fluidState,
584 const ParameterCache<ParamCacheEval>& paramCache,
585 unsigned phaseIdx)
586 { return enthalpy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
587
588 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
589 STATIC_OR_DEVICE LhsEval internalEnergy(const FluidState& fluidState,
590 const ParameterCache<ParamCacheEval>& paramCache,
591 unsigned phaseIdx) NOTHING_OR_CONST
592 { return internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
593
594 /****************************************
595 * thermodynamic quantities (black-oil specific version: Note that the PVT region
596 * index is explicitly passed instead of a parameter cache object)
597 ****************************************/
599 template <class FluidState, class LhsEval = typename FluidState::ValueType>
600 STATIC_OR_DEVICE LhsEval density(const FluidState& fluidState,
601 unsigned phaseIdx,
602 unsigned regionIdx) NOTHING_OR_CONST
603 {
604 assert(phaseIdx <= numPhases);
605 assert(regionIdx <= numRegions());
606
607 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
608 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
609 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
610
611 switch (phaseIdx) {
612 case oilPhaseIdx: {
613 if (enableConstantRs()) {
614 // dead oil but positive constant Rs
615 const LhsEval& Rs = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
616 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
617
618 return
619 bo*referenceDensity(oilPhaseIdx, regionIdx)
620 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
621 }
622
623 if (enableDissolvedGas()) {
624 // miscible oil
625 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
626 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
627
628 return
629 bo*referenceDensity(oilPhaseIdx, regionIdx)
630 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
631 }
632
633 // immiscible oil
634 const LhsEval Rs(0.0);
635 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
636
637 return referenceDensity(phaseIdx, regionIdx)*bo;
638 }
639
640 case gasPhaseIdx: {
642 // gas containing vaporized oil and vaporized water
643 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
644 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
645 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
646
647 return
648 bg*referenceDensity(gasPhaseIdx, regionIdx)
649 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
650 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
651 }
652 if (enableVaporizedOil()) {
653 // miscible gas
654 const LhsEval Rvw(0.0);
655 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
656 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
657
658 return
659 bg*referenceDensity(gasPhaseIdx, regionIdx)
660 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
661 }
662 if (enableVaporizedWater()) {
663 // gas containing vaporized water
664 const LhsEval Rv(0.0);
665 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
666 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
667
668 return
669 bg*referenceDensity(gasPhaseIdx, regionIdx)
670 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
671 }
672
673 // immiscible gas
674 const LhsEval Rv(0.0);
675 const LhsEval Rvw(0.0);
676 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
677 return bg*referenceDensity(phaseIdx, regionIdx);
678 }
679
680 case waterPhaseIdx:
682 // gas miscible in water
683 const LhsEval& Rsw =BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
684 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
685 return
686 bw*referenceDensity(waterPhaseIdx, regionIdx)
687 + Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
688 }
689 const LhsEval Rsw(0.0);
690 return
692 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
693 }
694
695 throw std::logic_error("Unhandled phase index " + std::to_string(phaseIdx));
696 }
697
707 template <class FluidState, class LhsEval = typename FluidState::ValueType>
708 STATIC_OR_DEVICE LhsEval saturatedDensity(const FluidState& fluidState,
709 unsigned phaseIdx,
710 unsigned regionIdx) NOTHING_OR_CONST
711 {
712 assert(phaseIdx <= numPhases);
713 assert(regionIdx <= numRegions());
714
715 const auto& p = fluidState.pressure(phaseIdx);
716 const auto& T = fluidState.temperature(phaseIdx);
717
718 switch (phaseIdx) {
719 case oilPhaseIdx: {
720 if (enableConstantRs()) {
721 // dead oil but positive constant Rs
722 const LhsEval& Rs = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
723 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
724
725 return
726 bo*referenceDensity(oilPhaseIdx, regionIdx)
727 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
728 }
729
730 if (enableDissolvedGas()) {
731 // miscible oil
732 const LhsEval& Rs = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, oilPhaseIdx, regionIdx);
733 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
734
735 return
736 bo*referenceDensity(oilPhaseIdx, regionIdx)
737 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
738 }
739
740 // immiscible oil
741 const LhsEval Rs(0.0);
742 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
743 return referenceDensity(phaseIdx, regionIdx)*bo;
744 }
745
746 case gasPhaseIdx: {
748 // gas containing vaporized oil and vaporized water
749 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
750 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
751 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
752
753 return
754 bg*referenceDensity(gasPhaseIdx, regionIdx)
755 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
756 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx) ;
757 }
758
759 if (enableVaporizedOil()) {
760 // miscible gas
761 const LhsEval Rvw(0.0);
762 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
763 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
764
765 return
766 bg*referenceDensity(gasPhaseIdx, regionIdx)
767 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
768 }
769
770 if (enableVaporizedWater()) {
771 // gas containing vaporized water
772 const LhsEval Rv(0.0);
773 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
774 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
775
776 return
777 bg*referenceDensity(gasPhaseIdx, regionIdx)
778 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
779 }
780
781 // immiscible gas
782 const LhsEval Rv(0.0);
783 const LhsEval Rvw(0.0);
784 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
785
786 return referenceDensity(phaseIdx, regionIdx)*bg;
787
788 }
789
790 case waterPhaseIdx:
791 {
793 // miscible in water
794 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
795 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, waterPhaseIdx, regionIdx);
796 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
797 return
798 bw*referenceDensity(waterPhaseIdx, regionIdx)
799 + Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
800 }
801 return
804 }
805 }
806
807 throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
808 }
809
818 template <class FluidState, class LhsEval = typename FluidState::ValueType>
819 STATIC_OR_DEVICE LhsEval inverseFormationVolumeFactor(const FluidState& fluidState,
820 unsigned phaseIdx,
821 unsigned regionIdx) NOTHING_OR_CONST
822 {
823 OPM_TIMEBLOCK_LOCAL(inverseFormationVolumeFactor, Subsystem::PvtProps);
824 assert(phaseIdx <= numPhases);
825 assert(regionIdx <= numRegions());
826
827 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
828 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
829
830 switch (phaseIdx) {
831 case oilPhaseIdx: {
832 if (enableDissolvedGas()) {
833 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
834 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
835 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
836 {
837 return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
838 } else {
839 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
840 }
841 }
842
843 const LhsEval Rs(0.0);
844 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
845 }
846 case gasPhaseIdx: {
848 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
849 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
850 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
851 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
852 && fluidState.saturation(oilPhaseIdx) > 0.0
853 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
854 {
855 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
856 } else {
857 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
858 }
859 }
860
861 if (enableVaporizedOil()) {
862 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
863 if (useSaturatedTables() && fluidState.saturation(oilPhaseIdx) > 0.0
864 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
865 {
866 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
867 } else {
868 const LhsEval Rvw(0.0);
869 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
870 }
871 }
872
873 if (enableVaporizedWater()) {
874 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
875 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
876 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
877 {
878 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
879 } else {
880 const LhsEval Rv(0.0);
881 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
882 }
883 }
884
885 const LhsEval Rv(0.0);
886 const LhsEval Rvw(0.0);
887 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
888 }
889 case waterPhaseIdx:
890 {
891 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
893 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
894 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
895 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
896 {
897 return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
898 } else {
899 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
900 }
901 }
902 const LhsEval Rsw(0.0);
903 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
904 }
905 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
906 }
907 }
908
909 template <class FluidState, class LhsEval = typename FluidState::ValueType>
910 STATIC_OR_DEVICE std::pair<LhsEval, LhsEval>
911 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState,
912 unsigned phaseIdx,
913 unsigned regionIdx)
914 {
915 switch (phaseIdx) {
916 case oilPhaseIdx:
917 return oilPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
918 case gasPhaseIdx:
919 return gasPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
920 case waterPhaseIdx:
921 return waterPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
922 default:
923 throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
924 }
925 }
926
936 template <class FluidState, class LhsEval = typename FluidState::ValueType>
937 STATIC_OR_DEVICE LhsEval saturatedInverseFormationVolumeFactor(const FluidState& fluidState,
938 unsigned phaseIdx,
939 unsigned regionIdx) NOTHING_OR_CONST
940 {
941 OPM_TIMEBLOCK_LOCAL(saturatedInverseFormationVolumeFactor, Subsystem::PvtProps);
942 assert(phaseIdx <= numPhases);
943 assert(regionIdx <= numRegions());
944
945 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
946 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
947 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
948
949 switch (phaseIdx) {
950 case oilPhaseIdx: return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
951 case gasPhaseIdx: return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
952 case waterPhaseIdx: return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
953 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
954 }
955 }
956
958 template <class FluidState, class LhsEval = typename FluidState::ValueType>
959 STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState& fluidState,
960 unsigned phaseIdx,
961 unsigned compIdx,
962 unsigned regionIdx) NOTHING_OR_CONST
963 {
964 assert(phaseIdx <= numPhases);
965 assert(compIdx <= numComponents);
966 assert(regionIdx <= numRegions());
967
968 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
969 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
970
971 // for the fugacity coefficient of the oil component in the oil phase, we use
972 // some pseudo-realistic value for the vapor pressure to ease physical
973 // interpretation of the results
974 const LhsEval phi_oO = 20e3/p;
975
976 // for the gas component in the gas phase, assume it to be an ideal gas
977 constexpr const Scalar phi_gG = 1.0;
978
979 // for the fugacity coefficient of the water component in the water phase, we use
980 // the same approach as for the oil component in the oil phase
981 const LhsEval phi_wW = 30e3/p;
982
983 switch (phaseIdx) {
984 case gasPhaseIdx: // fugacity coefficients for all components in the gas phase
985 switch (compIdx) {
986 case gasCompIdx:
987 return phi_gG;
988
989 // for the oil component, we calculate the Rv value for saturated gas and Rs
990 // for saturated oil, and compute the fugacity coefficients at the
991 // equilibrium. for this, we assume that in equilibrium the fugacities of the
992 // oil component is the same in both phases.
993 case oilCompIdx: {
994 if (!enableVaporizedOil())
995 // if there's no vaporized oil, the gas phase is assumed to be
996 // immiscible with the oil component
997 return phi_gG*1e6;
998
999 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1000 const auto& X_gOSat = convertRvToXgO(R_vSat, regionIdx);
1001 const auto& x_gOSat = convertXgOToxgO(X_gOSat, regionIdx);
1002
1003 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1004 const auto& X_oGSat = convertRsToXoG(R_sSat, regionIdx);
1005 const auto& x_oGSat = convertXoGToxoG(X_oGSat, regionIdx);
1006 const auto& x_oOSat = 1.0 - x_oGSat;
1007
1008 const auto& p_o = decay<LhsEval>(fluidState.pressure(oilPhaseIdx));
1009 const auto& p_g = decay<LhsEval>(fluidState.pressure(gasPhaseIdx));
1010
1011 return phi_oO*p_o*x_oOSat / (p_g*x_gOSat);
1012 }
1013
1014 case waterCompIdx:
1015 // the water component is assumed to be never miscible with the gas phase
1016 return phi_gG*1e6;
1017
1018 default:
1019 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
1020 }
1021
1022 case oilPhaseIdx: // fugacity coefficients for all components in the oil phase
1023 switch (compIdx) {
1024 case oilCompIdx:
1025 return phi_oO;
1026
1027 // for the oil and water components, we proceed analogous to the gas and
1028 // water components in the gas phase
1029 case gasCompIdx: {
1030 if (!enableDissolvedGas())
1031 // if there's no dissolved gas, the oil phase is assumed to be
1032 // immiscible with the gas component
1033 return phi_oO*1e6;
1034
1035 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1036 const auto& X_gOSat = convertRvToXgO(R_vSat, regionIdx);
1037 const auto& x_gOSat = convertXgOToxgO(X_gOSat, regionIdx);
1038 const auto& x_gGSat = 1.0 - x_gOSat;
1039
1040 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1041 const auto& X_oGSat = convertRsToXoG(R_sSat, regionIdx);
1042 const auto& x_oGSat = convertXoGToxoG(X_oGSat, regionIdx);
1043
1044 const auto& p_o = decay<LhsEval>(fluidState.pressure(oilPhaseIdx));
1045 const auto& p_g = decay<LhsEval>(fluidState.pressure(gasPhaseIdx));
1046
1047 return phi_gG*p_g*x_gGSat / (p_o*x_oGSat);
1048 }
1049
1050 case waterCompIdx:
1051 return phi_oO*1e6;
1052
1053 default:
1054 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
1055 }
1056
1057 case waterPhaseIdx: // fugacity coefficients for all components in the water phase
1058 // the water phase fugacity coefficients are pretty simple: because the water
1059 // phase is assumed to consist entirely from the water component, we just
1060 // need to make sure that the fugacity coefficients for the other components
1061 // are a few orders of magnitude larger than the one of the water
1062 // component. (i.e., the affinity of the gas and oil components to the water
1063 // phase is lower by a few orders of magnitude)
1064 switch (compIdx) {
1065 case waterCompIdx: return phi_wW;
1066 case oilCompIdx: return 1.1e6*phi_wW;
1067 case gasCompIdx: return 1e6*phi_wW;
1068 default:
1069 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
1070 }
1071
1072 default:
1073 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
1074 }
1075
1076 throw std::logic_error("Unhandled phase or component index");
1077 }
1078
1080 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1081 STATIC_OR_DEVICE LhsEval viscosity(const FluidState& fluidState,
1082 unsigned phaseIdx,
1083 unsigned regionIdx) NOTHING_OR_CONST
1084 {
1085 OPM_TIMEBLOCK_LOCAL(viscosity, Subsystem::PvtProps);
1086 assert(phaseIdx <= numPhases);
1087 assert(regionIdx <= numRegions());
1088
1089 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1090 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1091
1092 switch (phaseIdx) {
1093 case oilPhaseIdx: {
1094 if (enableDissolvedGas()) {
1095 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1096 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
1097 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
1098 {
1099 return oilPvt_.saturatedViscosity(regionIdx, T, p);
1100 } else {
1101 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1102 }
1103 }
1104
1105 const LhsEval Rs(0.0);
1106 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1107 }
1108
1109 case gasPhaseIdx: {
1111 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1112 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1113 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
1114 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
1115 && fluidState.saturation(oilPhaseIdx) > 0.0
1116 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1117 {
1118 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1119 } else {
1120 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1121 }
1122 }
1123 if (enableVaporizedOil()) {
1124 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1125 if (useSaturatedTables() && fluidState.saturation(oilPhaseIdx) > 0.0
1126 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1127 {
1128 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1129 } else {
1130 const LhsEval Rvw(0.0);
1131 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1132 }
1133 }
1134 if (enableVaporizedWater()) {
1135 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1136 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
1137 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1138 {
1139 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1140 } else {
1141 const LhsEval Rv(0.0);
1142 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1143 }
1144 }
1145
1146 const LhsEval Rv(0.0);
1147 const LhsEval Rvw(0.0);
1148 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1149 }
1150
1151 case waterPhaseIdx:
1152 {
1153 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1155 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1156 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
1157 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
1158 {
1159 return waterPvt_.saturatedViscosity(regionIdx, T, p, saltConcentration);
1160 } else {
1161 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1162 }
1163 }
1164 const LhsEval Rsw(0.0);
1165 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1166 }
1167 }
1168
1169 OPM_THROW(std::logic_error, "Unhandled phase index "+std::to_string(phaseIdx));
1170 }
1171
1172 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1173 STATIC_OR_DEVICE LhsEval internalEnergy(const FluidState& fluidState,
1174 const unsigned phaseIdx,
1175 const unsigned regionIdx) NOTHING_OR_CONST
1176 {
1177 const auto p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1178 const auto T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1179
1180 switch (phaseIdx) {
1181 case oilPhaseIdx:
1182 if (!oilPvt_.mixingEnergy()) {
1183 return oilPvt_.internalEnergy
1184 (regionIdx, T, p,
1185 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1186 }
1187 break;
1188
1189 case waterPhaseIdx:
1190 if (!waterPvt_.mixingEnergy()) {
1191 return waterPvt_.internalEnergy
1192 (regionIdx, T, p,
1193 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1194 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1195 }
1196 break;
1197
1198 case gasPhaseIdx:
1199 if (!gasPvt_.mixingEnergy()) {
1200 return gasPvt_.internalEnergy
1201 (regionIdx, T, p,
1202 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1203 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1204 }
1205 break;
1206
1207 default:
1208 throw std::logic_error {
1209 "Phase index " + std::to_string(phaseIdx) + " does not support internal energy"
1210 };
1211 }
1212
1213 return internalMixingTotalEnergy<FluidState,LhsEval>(fluidState, phaseIdx, regionIdx)
1214 / density<FluidState,LhsEval>(fluidState, phaseIdx, regionIdx);
1215 }
1216
1217
1218 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1219 STATIC_OR_DEVICE LhsEval internalMixingTotalEnergy(const FluidState& fluidState,
1220 unsigned phaseIdx,
1221 unsigned regionIdx) NOTHING_OR_CONST
1222 {
1223 assert(phaseIdx <= numPhases);
1224 assert(regionIdx <= numRegions());
1225 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1226 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1227 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1228 // to avoid putting all thermal into the interface of the multiplexer
1229 switch (phaseIdx) {
1230 case oilPhaseIdx: {
1231 auto oilEnergy = oilPvt_.internalEnergy(regionIdx, T, p,
1232 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1233 assert(oilPvt_.mixingEnergy());
1234 //mixing energy adsed
1235 if (enableDissolvedGas()) {
1236 // miscible oil
1237 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1238 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1239 const auto& gasEnergy =
1240 gasPvt_.internalEnergy(regionIdx, T, p,
1241 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1242 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1243 const auto hVapG = gasPvt_.hVap(regionIdx);// pressure correction ? assume equal to energy change
1244 return
1245 oilEnergy*bo*referenceDensity(oilPhaseIdx, regionIdx)
1246 + (gasEnergy-hVapG)*Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
1247 }
1248
1249 // immiscible oil
1250 const LhsEval Rs(0.0);
1251 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1252
1253 return oilEnergy*referenceDensity(phaseIdx, regionIdx)*bo;
1254 }
1255
1256 case gasPhaseIdx: {
1257 const auto& gasEnergy =
1258 gasPvt_.internalEnergy(regionIdx, T, p,
1259 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1260 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1261 assert(gasPvt_.mixingEnergy());
1263 const auto& oilEnergy =
1264 oilPvt_.internalEnergy(regionIdx, T, p,
1265 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1266 const auto waterEnergy =
1267 waterPvt_.internalEnergy(regionIdx, T, p,
1268 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1269 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1270 // gas containing vaporized oil and vaporized water
1271 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1272 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1273 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1274 const auto hVapO = oilPvt_.hVap(regionIdx);
1275 const auto hVapW = waterPvt_.hVap(regionIdx);
1276 return
1277 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1278 + (oilEnergy+hVapO)*Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
1279 + (waterEnergy+hVapW)*Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
1280 }
1281 if (enableVaporizedOil()) {
1282 const auto& oilEnergy =
1283 oilPvt_.internalEnergy(regionIdx, T, p,
1284 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1285 // miscible gas
1286 const LhsEval Rvw(0.0);
1287 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1288 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1289 const auto hVapO = oilPvt_.hVap(regionIdx);
1290 return
1291 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1292 + (oilEnergy+hVapO)*Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
1293 }
1294 if (enableVaporizedWater()) {
1295 // gas containing vaporized water
1296 const LhsEval Rv(0.0);
1297 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1298 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1299 const auto waterEnergy =
1300 waterPvt_.internalEnergy(regionIdx, T, p,
1301 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1302 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1303 const auto hVapW = waterPvt_.hVap(regionIdx);
1304 return
1305 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1306 + (waterEnergy+hVapW)*Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
1307 }
1308
1309 // immiscible gas
1310 const LhsEval Rv(0.0);
1311 const LhsEval Rvw(0.0);
1312 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1313 return gasEnergy*bg*referenceDensity(phaseIdx, regionIdx);
1314 }
1315
1316 case waterPhaseIdx:
1317 const auto waterEnergy =
1318 waterPvt_.internalEnergy(regionIdx, T, p,
1319 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1320 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1321 assert(waterPvt_.mixingEnergy());
1323 const auto& gasEnergy =
1324 gasPvt_.internalEnergy(regionIdx, T, p,
1325 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1326 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1327 // gas miscible in water
1328 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, waterPhaseIdx, regionIdx);
1329 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1330 return
1331 waterEnergy*bw*referenceDensity(waterPhaseIdx, regionIdx)
1332 + gasEnergy*Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
1333 }
1334 const LhsEval Rsw(0.0);
1335 return
1336 waterEnergy*referenceDensity(waterPhaseIdx, regionIdx)
1337 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1338 }
1339 throw std::logic_error("Unhandled phase index " + std::to_string(phaseIdx));
1340 }
1341
1342
1343
1345 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1346 STATIC_OR_DEVICE LhsEval enthalpy(const FluidState& fluidState,
1347 unsigned phaseIdx,
1348 unsigned regionIdx) NOTHING_OR_CONST
1349 {
1350 // should preferably not be used values should be taken from intensive quantities fluid state.
1351 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1352 auto energy = internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1353 if(!enthalpy_eq_energy_){
1354 // used for simplified models
1355 energy += p/density<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1356 }
1357 return energy;
1358 }
1359
1366 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1367 STATIC_OR_DEVICE LhsEval saturatedVaporizationFactor(const FluidState& fluidState,
1368 unsigned phaseIdx,
1369 unsigned regionIdx) NOTHING_OR_CONST
1370 {
1371 assert(phaseIdx <= numPhases);
1372 assert(regionIdx <= numRegions());
1373
1374 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1375 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1376 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
1377
1378 switch (phaseIdx) {
1379 case oilPhaseIdx: return 0.0;
1380 case gasPhaseIdx: return gasPvt_.saturatedWaterVaporizationFactor(regionIdx, T, p, saltConcentration);
1381 case waterPhaseIdx: return 0.0;
1382 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1383 }
1384 }
1385
1392 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1393 STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState& fluidState,
1394 unsigned phaseIdx,
1395 unsigned regionIdx,
1396 const LhsEval& maxOilSaturation) NOTHING_OR_CONST
1397 {
1398 OPM_TIMEBLOCK_LOCAL(saturatedDissolutionFactor, Subsystem::PvtProps);
1399 assert(phaseIdx <= numPhases);
1400 assert(regionIdx <= numRegions());
1401
1402 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1403 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1404 const auto& So = (phaseIdx == waterPhaseIdx) ? 0 : decay<LhsEval>(fluidState.saturation(oilPhaseIdx));
1405
1406 switch (phaseIdx) {
1407 case oilPhaseIdx: return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p, So, maxOilSaturation);
1408 case gasPhaseIdx: return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p, So, maxOilSaturation);
1409 case waterPhaseIdx: return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1410 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1411 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1412 }
1413 }
1414
1423 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1424 STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState& fluidState,
1425 unsigned phaseIdx,
1426 unsigned regionIdx) NOTHING_OR_CONST
1427 {
1428 OPM_TIMEBLOCK_LOCAL(saturatedDissolutionFactor, Subsystem::PvtProps);
1429 assert(phaseIdx <= numPhases);
1430 assert(regionIdx <= numRegions());
1431
1432 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1433 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1434
1435 switch (phaseIdx) {
1436 case oilPhaseIdx: return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1437 case gasPhaseIdx: return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1438 case waterPhaseIdx: return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1439 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1440 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1441 }
1442 }
1443
1447 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1448 STATIC_OR_DEVICE LhsEval bubblePointPressure(const FluidState& fluidState,
1449 unsigned regionIdx) NOTHING_OR_CONST
1450 {
1451 return saturationPressure(fluidState, oilPhaseIdx, regionIdx);
1452 }
1453
1454
1458 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1459 STATIC_OR_DEVICE LhsEval dewPointPressure(const FluidState& fluidState,
1460 unsigned regionIdx) NOTHING_OR_CONST
1461 {
1462 return saturationPressure(fluidState, gasPhaseIdx, regionIdx);
1463 }
1464
1475 template <class FluidState, class LhsEval = typename FluidState::ValueType>
1476 STATIC_OR_DEVICE LhsEval saturationPressure(const FluidState& fluidState,
1477 unsigned phaseIdx,
1478 unsigned regionIdx) NOTHING_OR_CONST
1479 {
1480 assert(phaseIdx <= numPhases);
1481 assert(regionIdx <= numRegions());
1482
1483 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1484
1485 switch (phaseIdx) {
1486 case oilPhaseIdx: return oilPvt_.saturationPressure(regionIdx, T, BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1487 case gasPhaseIdx: return gasPvt_.saturationPressure(regionIdx, T, BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1488 case waterPhaseIdx: return waterPvt_.saturationPressure(regionIdx, T,
1489 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1490 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1491 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1492 }
1493 }
1494
1495 /****************************************
1496 * Auxiliary and convenience methods for the black-oil model
1497 ****************************************/
1502 template <class LhsEval>
1503 STATIC_OR_DEVICE LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) NOTHING_OR_CONST
1504 {
1505 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1506 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1507
1508 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
1509 }
1510
1515 template <class LhsEval>
1516 STATIC_OR_DEVICE LhsEval convertXwGToRsw(const LhsEval& XwG, unsigned regionIdx) NOTHING_OR_CONST
1517 {
1518 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1519 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1520
1521 return XwG/(1.0 - XwG)*(rho_wRef/rho_gRef);
1522 }
1523
1528 template <class LhsEval>
1529 STATIC_OR_DEVICE LhsEval convertXgOToRv(const LhsEval& XgO, unsigned regionIdx) NOTHING_OR_CONST
1530 {
1531 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1532 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1533
1534 return XgO/(1.0 - XgO)*(rho_gRef/rho_oRef);
1535 }
1536
1541 template <class LhsEval>
1542 STATIC_OR_DEVICE LhsEval convertXgWToRvw(const LhsEval& XgW, unsigned regionIdx) NOTHING_OR_CONST
1543 {
1544 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1545 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1546
1547 return XgW/(1.0 - XgW)*(rho_gRef/rho_wRef);
1548 }
1549
1550
1555 template <class LhsEval>
1556 STATIC_OR_DEVICE LhsEval convertRsToXoG(const LhsEval& Rs, unsigned regionIdx) NOTHING_OR_CONST
1557 {
1558 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1559 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1560
1561 const LhsEval& rho_oG = Rs*rho_gRef;
1562 return rho_oG/(rho_oRef + rho_oG);
1563 }
1564
1569 template <class LhsEval>
1570 STATIC_OR_DEVICE LhsEval convertRswToXwG(const LhsEval& Rsw, unsigned regionIdx) NOTHING_OR_CONST
1571 {
1572 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1573 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1574
1575 const LhsEval& rho_wG = Rsw*rho_gRef;
1576 return rho_wG/(rho_wRef + rho_wG);
1577 }
1578
1583 template <class LhsEval>
1584 STATIC_OR_DEVICE LhsEval convertRvToXgO(const LhsEval& Rv, unsigned regionIdx) NOTHING_OR_CONST
1585 {
1586 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1587 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1588
1589 const LhsEval& rho_gO = Rv*rho_oRef;
1590 return rho_gO/(rho_gRef + rho_gO);
1591 }
1592
1597 template <class LhsEval>
1598 STATIC_OR_DEVICE LhsEval convertRvwToXgW(const LhsEval& Rvw, unsigned regionIdx) NOTHING_OR_CONST
1599 {
1600 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1601 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1602
1603 const LhsEval& rho_gW = Rvw*rho_wRef;
1604 return rho_gW/(rho_gRef + rho_gW);
1605 }
1606
1610 template <class LhsEval>
1611 STATIC_OR_DEVICE LhsEval convertXgWToxgW(const LhsEval& XgW, unsigned regionIdx) NOTHING_OR_CONST
1612 {
1613 Scalar MW = molarMass_[regionIdx][waterCompIdx];
1614 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1615
1616 return XgW*MG / (MW*(1 - XgW) + XgW*MG);
1617 }
1618
1622 template <class LhsEval>
1623 STATIC_OR_DEVICE LhsEval convertXwGToxwG(const LhsEval& XwG, unsigned regionIdx) NOTHING_OR_CONST
1624 {
1625 Scalar MW = molarMass_[regionIdx][waterCompIdx];
1626 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1627
1628 return XwG*MW / (MG*(1 - XwG) + XwG*MW);
1629 }
1630
1634 template <class LhsEval>
1635 STATIC_OR_DEVICE LhsEval convertXoGToxoG(const LhsEval& XoG, unsigned regionIdx) NOTHING_OR_CONST
1636 {
1637 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1638 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1639
1640 return XoG*MO / (MG*(1 - XoG) + XoG*MO);
1641 }
1642
1646 template <class LhsEval>
1647 STATIC_OR_DEVICE LhsEval convertxoGToXoG(const LhsEval& xoG, unsigned regionIdx) NOTHING_OR_CONST
1648 {
1649 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1650 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1651
1652 return xoG*MG / (xoG*(MG - MO) + MO);
1653 }
1654
1658 template <class LhsEval>
1659 STATIC_OR_DEVICE LhsEval convertXgOToxgO(const LhsEval& XgO, unsigned regionIdx) NOTHING_OR_CONST
1660 {
1661 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1662 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1663
1664 return XgO*MG / (MO*(1 - XgO) + XgO*MG);
1665 }
1666
1670 template <class LhsEval>
1671 STATIC_OR_DEVICE LhsEval convertxgOToXgO(const LhsEval& xgO, unsigned regionIdx) NOTHING_OR_CONST
1672 {
1673 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1674 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1675
1676 return xgO*MO / (xgO*(MO - MG) + MG);
1677 }
1678
1686 STATIC_OR_DEVICE const GasPvt& gasPvt() NOTHING_OR_CONST
1687 { return gasPvt_; }
1688
1696 STATIC_OR_DEVICE const OilPvt& oilPvt() NOTHING_OR_CONST
1697 { return oilPvt_; }
1698
1706 STATIC_OR_DEVICE const WaterPvt& waterPvt() NOTHING_OR_CONST
1707 { return waterPvt_; }
1708
1714 STATIC_OR_DEVICE Scalar reservoirTemperature(unsigned = 0) NOTHING_OR_CONST
1715 { return reservoirTemperature_; }
1716
1722 STATIC_OR_DEVICE void setReservoirTemperature(Scalar value) NOTHING_OR_CONST
1723 { reservoirTemperature_ = value; }
1724
1725 STATIC_OR_DEVICE short activeToCanonicalPhaseIdx(unsigned activePhaseIdx) NOTHING_OR_CONST;
1726
1727 STATIC_OR_DEVICE short canonicalToActivePhaseIdx(unsigned phaseIdx) NOTHING_OR_CONST;
1728
1729 STATIC_OR_DEVICE short activeToCanonicalCompIdx(unsigned activeCompIdx) NOTHING_OR_CONST;
1730
1731 STATIC_OR_DEVICE short canonicalToActiveCompIdx(unsigned compIdx) NOTHING_OR_CONST;
1732
1733 STATIC_OR_DEVICE short activePhaseToActiveCompIdx(unsigned activePhaseIdx) NOTHING_OR_CONST;
1734
1735 STATIC_OR_DEVICE short activeCompToActivePhaseIdx(unsigned activeCompIdx) NOTHING_OR_CONST;
1736
1738 STATIC_OR_DEVICE Scalar diffusionCoefficient(unsigned compIdx, unsigned phaseIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
1739 { return diffusionCoefficients_[regionIdx][numPhases*compIdx + phaseIdx]; }
1740
1742 STATIC_OR_DEVICE void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
1743 { diffusionCoefficients_[regionIdx][numPhases*compIdx + phaseIdx] = coefficient ; }
1744
1748 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
1749 STATIC_OR_DEVICE LhsEval diffusionCoefficient(const FluidState& fluidState,
1750 const ParameterCache<ParamCacheEval>& paramCache,
1751 unsigned phaseIdx,
1752 unsigned compIdx) NOTHING_OR_CONST
1753 {
1754 // diffusion is disabled by the user
1755 if(!enableDiffusion())
1756 return 0.0;
1757
1758 // diffusion coefficients are set, and we use them
1759 if(!diffusionCoefficients_.empty()) {
1760 return diffusionCoefficient(compIdx, phaseIdx, paramCache.regionIndex());
1761 }
1762
1763 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1764 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1765 const unsigned regionIdx = paramCache.regionIndex();
1766
1767 switch (phaseIdx) {
1768 case oilPhaseIdx: return oilPvt().diffusionCoefficient(T, p, compIdx);
1769 case gasPhaseIdx: return gasPvt().diffusionCoefficient(T, p, compIdx);
1770 case waterPhaseIdx: return waterPvt().diffusionCoefficient(T, p, compIdx, regionIdx);
1771 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1772 }
1773 }
1774 STATIC_OR_DEVICE void setEnergyEqualEnthalpy(bool enthalpy_eq_energy){
1775 enthalpy_eq_energy_ = enthalpy_eq_energy;
1776 }
1777
1778 STATIC_OR_DEVICE bool enthalpyEqualEnergy() NOTHING_OR_CONST{
1779 return enthalpy_eq_energy_;
1780 }
1781
1782private:
1783 STATIC_OR_NOTHING void resizeArrays_(std::size_t numRegions);
1784
1785 STATIC_OR_NOTHING Scalar reservoirTemperature_;
1786
1787 STATIC_OR_NOTHING GasPvt gasPvt_;
1788 STATIC_OR_NOTHING OilPvt oilPvt_;
1789 STATIC_OR_NOTHING WaterPvt waterPvt_;
1790
1791 STATIC_OR_NOTHING bool enableDissolvedGas_;
1792 STATIC_OR_NOTHING bool enableDissolvedGasInWater_;
1793 STATIC_OR_NOTHING bool enableVaporizedOil_;
1794 STATIC_OR_NOTHING bool enableVaporizedWater_;
1795 STATIC_OR_NOTHING bool enableConstantRs_;
1796 STATIC_OR_NOTHING bool enableDiffusion_;
1797
1798 // HACK for GCC 4.4: the array size has to be specified using the literal value '3'
1799 // here, because GCC 4.4 seems to be unable to determine the number of phases from
1800 // the BlackOil fluid system in the attribute declaration below...
1801 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numPhases=*/3> > referenceDensity_;
1802 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numComponents=*/3> > molarMass_;
1803 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numComponents=*/3 * /*numPhases=*/3> > diffusionCoefficients_;
1804
1805 STATIC_OR_NOTHING PhaseUsageInfo<IndexTraits> phaseUsageInfo_;
1806
1807 STATIC_OR_NOTHING bool isInitialized_;
1808 STATIC_OR_NOTHING bool useSaturatedTables_;
1809 STATIC_OR_NOTHING bool enthalpy_eq_energy_;
1810
1811 #ifndef COMPILING_STATIC_FLUID_SYSTEM
1812 template<template<typename> typename StorageT>
1813 explicit FLUIDSYSTEM_CLASSNAME(const FLUIDSYSTEM_CLASSNAME_STATIC<Scalar, IndexTraits, StorageT>& other)
1816 , reservoirTemperature_(other.reservoirTemperature_)
1817 , gasPvt_(other.gasPvt_)
1818 , oilPvt_(other.oilPvt_)
1819 , waterPvt_(other.waterPvt_)
1820 , enableDissolvedGas_(other.enableDissolvedGas_)
1821 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
1822 , enableVaporizedOil_(other.enableVaporizedOil_)
1823 , enableVaporizedWater_(other.enableVaporizedWater_)
1824 , enableConstantRs_(other.enableConstantRs_)
1825 , enableDiffusion_(other.enableDiffusion_)
1826 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
1827 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
1828 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
1829 , phaseUsageInfo_(other.phaseUsageInfo_)
1830 , isInitialized_(other.isInitialized_)
1831 , useSaturatedTables_(other.useSaturatedTables_)
1832 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
1833 {
1834 OPM_ERROR_IF(!other.isInitialized(), "The fluid system must be initialized before it can be copied.");
1835 }
1836
1837 template<class ScalarT, class IndexTraitsT, template<typename> typename StorageT>
1838 friend class FLUIDSYSTEM_CLASSNAME_STATIC;
1839 #else
1840 template<class ScalarT, class IndexTraitsT, template<typename> typename StorageT>
1841 friend class FLUIDSYSTEM_CLASSNAME_NONSTATIC;
1842 #endif
1843};
1844
1845template <class Scalar, class IndexTraits, template<typename> typename Storage>
1847initBegin(std::size_t numPvtRegions)
1848{
1849 isInitialized_ = false;
1850 useSaturatedTables_ = true;
1851
1852 enableDissolvedGas_ = true;
1853 enableDissolvedGasInWater_ = false;
1854 enableVaporizedOil_ = false;
1855 enableVaporizedWater_ = false;
1856 enableConstantRs_ = false;
1857 enableDiffusion_ = false;
1858
1859 surfaceTemperature = 273.15 + 15.56; // [K]
1860 surfacePressure = 1.01325e5; // [Pa]
1862
1863 phaseUsageInfo_.initFromPhases(Phases{true, true, true});
1864
1865 resizeArrays_(numPvtRegions);
1866}
1867
1868template <class Scalar, class IndexTraits, template<typename> typename Storage>
1871 Scalar rhoWater,
1872 Scalar rhoGas,
1873 unsigned regionIdx)
1874{
1875 referenceDensity_[regionIdx][oilPhaseIdx] = rhoOil;
1876 referenceDensity_[regionIdx][waterPhaseIdx] = rhoWater;
1877 referenceDensity_[regionIdx][gasPhaseIdx] = rhoGas;
1878}
1879
1880template <class Scalar, class IndexTraits, template<typename> typename Storage>
1882{
1883 // calculate the final 2D functions which are used for interpolation.
1884 const std::size_t num_regions = molarMass_.size();
1885 for (unsigned regionIdx = 0; regionIdx < num_regions; ++regionIdx) {
1886 // calculate molar masses
1887
1888 // water is simple: 18 g/mol
1889 molarMass_[regionIdx][waterCompIdx] = 18e-3;
1890
1892 // for gas, we take the density at standard conditions and assume it to be ideal
1895 Scalar rho_g = referenceDensity_[/*regionIdx=*/0][gasPhaseIdx];
1896 molarMass_[regionIdx][gasCompIdx] = Constants<Scalar>::R*T*rho_g / p;
1897 }
1898 else
1899 // hydrogen gas. we just set this do avoid NaNs later
1900 molarMass_[regionIdx][gasCompIdx] = 2e-3;
1901
1902 // finally, for oil phase, we take the molar mass from the spe9 paper
1903 molarMass_[regionIdx][oilCompIdx] = 175e-3; // kg/mol
1904 }
1905
1906 isInitialized_ = true;
1907}
1908
1909template <class Scalar, class IndexTraits, template<typename> typename Storage>
1911phaseName(unsigned phaseIdx) NOTHING_OR_CONST
1912{
1913 switch (phaseIdx) {
1914 case waterPhaseIdx:
1915 return "water";
1916 case oilPhaseIdx:
1917 return "oil";
1918 case gasPhaseIdx:
1919 return "gas";
1920
1921 default:
1922 OPM_THROW(std::logic_error, "Phase index " + std::to_string(phaseIdx) + " is unknown");
1923 }
1924}
1925
1926template <class Scalar, class IndexTraits, template<typename> typename Storage>
1928solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
1929{
1930 switch (phaseIdx) {
1931 case waterPhaseIdx:
1932 return waterCompIdx;
1933 case oilPhaseIdx:
1934 return oilCompIdx;
1935 case gasPhaseIdx:
1936 return gasCompIdx;
1937
1938 default:
1939 OPM_THROW(std::logic_error, "Phase index " + std::to_string(phaseIdx) + " is unknown");
1940 }
1941}
1942
1943template <class Scalar, class IndexTraits, template<typename> typename Storage>
1945soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
1946{
1947 switch (phaseIdx) {
1948 case waterPhaseIdx:
1950 return gasCompIdx;
1951 OPM_THROW(std::logic_error, "The water phase does not have any solutes in the black oil model!");
1952 case oilPhaseIdx:
1953 return gasCompIdx;
1954 case gasPhaseIdx:
1955 if (enableVaporizedWater()) {
1956 return waterCompIdx;
1957 }
1958 return oilCompIdx;
1959
1960 default:
1961 OPM_THROW(std::logic_error, "Phase index " + std::to_string(phaseIdx) + " is unknown");
1962 }
1963}
1964
1965template <class Scalar, class IndexTraits, template<typename> typename Storage>
1967componentName(unsigned compIdx) NOTHING_OR_CONST
1968{
1969 switch (compIdx) {
1970 case waterCompIdx:
1971 return "Water";
1972 case oilCompIdx:
1973 return "Oil";
1974 case gasCompIdx:
1975 return "Gas";
1976
1977 default:
1978 OPM_THROW(std::logic_error, "Component index " + std::to_string(compIdx) + " is unknown");
1979 }
1980}
1981
1982template <class Scalar, class IndexTraits, template<typename> typename Storage>
1984activeToCanonicalPhaseIdx(unsigned activePhaseIdx) NOTHING_OR_CONST
1985{
1986 return phaseUsageInfo_.activeToCanonicalPhaseIdx(activePhaseIdx);
1987}
1988
1989template <class Scalar, class IndexTraits, template<typename> typename Storage>
1990NOTHING_OR_DEVICE short FLUIDSYSTEM_CLASSNAME<Scalar,IndexTraits, Storage>::
1991canonicalToActivePhaseIdx(unsigned phaseIdx) NOTHING_OR_CONST
1992{
1993 return phaseUsageInfo_.canonicalToActivePhaseIdx(phaseIdx);
1994}
1995
1996template <class Scalar, class IndexTraits, template<typename> typename Storage>
1998activeToCanonicalCompIdx(unsigned activeCompIdx) NOTHING_OR_CONST
1999{
2000 return phaseUsageInfo_.activeToCanonicalCompIdx(activeCompIdx);
2001}
2002
2003template <class Scalar, class IndexTraits, template<typename> typename Storage>
2005canonicalToActiveCompIdx(unsigned compIdx) NOTHING_OR_CONST
2006{
2007 return phaseUsageInfo_.canonicalToActiveCompIdx(compIdx);
2008}
2009
2010template <class Scalar, class IndexTraits, template<typename> typename Storage>
2012activePhaseToActiveCompIdx(unsigned activePhaseIdx) NOTHING_OR_CONST
2013{
2014 return phaseUsageInfo_.activePhaseToActiveCompIdx(activePhaseIdx);
2015}
2016
2017template <class Scalar, class IndexTraits, template<typename> typename Storage>
2019activeCompToActivePhaseIdx(unsigned activeCompIdx) NOTHING_OR_CONST
2020{
2021 return phaseUsageInfo_.activeCompToActivePhaseIdx(activeCompIdx);
2022}
2023
2024template <class Scalar, class IndexTraits, template<typename> typename Storage>
2026resizeArrays_(std::size_t numRegions)
2027{
2028 molarMass_.resize(numRegions);
2029 referenceDensity_.resize(numRegions);
2030}
2031
2032#ifdef COMPILING_STATIC_FLUID_SYSTEM
2034
2035#define DECLARE_INSTANCE(T) \
2036template<> PhaseUsageInfo<BlackOilDefaultFluidSystemIndices> BOFS<T>::phaseUsageInfo_;\
2037template<> T BOFS<T>::surfaceTemperature; \
2038template<> T BOFS<T>::surfacePressure; \
2039template<> T BOFS<T>::reservoirTemperature_; \
2040template<> bool BOFS<T>::enableDissolvedGas_; \
2041template<> bool BOFS<T>::enableDissolvedGasInWater_; \
2042template<> bool BOFS<T>::enableVaporizedOil_; \
2043template<> bool BOFS<T>::enableVaporizedWater_; \
2044template<> bool BOFS<T>::enableConstantRs_; \
2045template<> bool BOFS<T>::enableDiffusion_; \
2046template<> BOFS<T>::OilPvt BOFS<T>::oilPvt_; \
2047template<> BOFS<T>::GasPvt BOFS<T>::gasPvt_; \
2048template<> BOFS<T>::WaterPvt BOFS<T>::waterPvt_; \
2049template<> std::vector<std::array<T, 3>> BOFS<T>::referenceDensity_; \
2050template<> std::vector<std::array<T, 3>> BOFS<T>::molarMass_; \
2051template<> std::vector<std::array<T, 9>> BOFS<T>::diffusionCoefficients_; \
2052template<> bool BOFS<T>::isInitialized_; \
2053template<> bool BOFS<T>::useSaturatedTables_; \
2054template<> bool BOFS<T>::enthalpy_eq_energy_;
2055
2056DECLARE_INSTANCE(float)
2057DECLARE_INSTANCE(double)
2058
2059#undef DECLARE_INSTANCE
2060#endif
2061
2062#ifndef COMPILING_STATIC_FLUID_SYSTEM
2063#if HAVE_CUDA
2064namespace gpuistl
2065{
2066
2067template <class Scalar, class IndexTraits>
2068FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>
2069copy_to_gpu(const FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits>& oldFluidSystem) {
2070
2071 using GpuBuffer3Array = GpuBuffer<std::array<Scalar, 3>>;
2072 using GpuBuffer9Array = GpuBuffer<std::array<Scalar, 9>>;
2073
2074 OPM_ERROR_IF(oldFluidSystem.gasPvt_.approach() != GasPvtApproach::Co2Gas,
2075 fmt::format(fmt::runtime("Incompatible gas PVT approach. Given {}, expected {} (GasPvtApproach::Co2Gas)."),
2076 int(oldFluidSystem.gasPvt_.approach()), int(GasPvtApproach::Co2Gas)));
2077
2078 OPM_ERROR_IF(oldFluidSystem.waterPvt_.approach() != WaterPvtApproach::BrineCo2,
2079 fmt::format(fmt::runtime("Incompatible water PVT approach. Given {}, expected {} (WaterPvtApproach::BrineCo2)."),
2080 int(oldFluidSystem.waterPvt_.approach()), int(WaterPvtApproach::BrineCo2)));
2081
2082 OPM_ERROR_IF(oldFluidSystem.oilPvt_.isActive(),
2083 "We currently do not support an active oil phase for the FluidSystem on GPU.");
2084
2085 auto newGasPvt = copy_to_gpu(oldFluidSystem.gasPvt_.template getRealPvt<GasPvtApproach::Co2Gas>());
2086
2087 auto newOilPvt = NullOilPvt<Scalar>();
2088
2089 auto newWaterPvt = copy_to_gpu(oldFluidSystem.waterPvt_.template getRealPvt<WaterPvtApproach::BrineCo2>());
2090
2091 auto newReferenceDensity = GpuBuffer3Array(oldFluidSystem.referenceDensity_);
2092 auto newMolarMass = GpuBuffer3Array(oldFluidSystem.molarMass_);
2093 auto newDiffusionCoefficients = GpuBuffer9Array(oldFluidSystem.diffusionCoefficients_);
2094
2095 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>(
2096 oldFluidSystem.surfacePressure,
2097 oldFluidSystem.surfaceTemperature,
2098 oldFluidSystem.reservoirTemperature_,
2099 newGasPvt,
2100 newOilPvt,
2101 newWaterPvt,
2102 oldFluidSystem.enableDissolvedGas_,
2103 oldFluidSystem.enableDissolvedGasInWater_,
2104 oldFluidSystem.enableVaporizedOil_,
2105 oldFluidSystem.enableVaporizedWater_,
2106 oldFluidSystem.enableConstantRs_,
2107 oldFluidSystem.enableDiffusion_,
2108 std::move(newReferenceDensity),
2109 std::move(newMolarMass),
2110 std::move(newDiffusionCoefficients),
2111 oldFluidSystem.phaseUsageInfo_,
2112 oldFluidSystem.isInitialized_,
2113 oldFluidSystem.useSaturatedTables_,
2114 oldFluidSystem.enthalpy_eq_energy_
2115 );
2116}
2117
2118template <class Scalar, class IndexTraits>
2119FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>
2120make_view(FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>& oldFluidSystem)
2121{
2122 using Array3 = std::array<Scalar, 3>;
2123 using Array9 = std::array<Scalar, 9>;
2124
2125 auto newGasPvt = make_view(oldFluidSystem.gasPvt_);
2126 auto newOilPvt = NullOilPvt<Scalar>();
2127 auto newWaterPvt = make_view(oldFluidSystem.waterPvt_);
2128
2129 auto newReferenceDensity = make_view<Array3>(oldFluidSystem.referenceDensity_);
2130 auto newMolarMass = make_view<Array3>(oldFluidSystem.molarMass_);
2131 auto newDiffusionCoefficients = make_view<Array9>(oldFluidSystem.diffusionCoefficients_);
2132
2133 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>(
2134 oldFluidSystem.surfacePressure,
2135 oldFluidSystem.surfaceTemperature,
2136 oldFluidSystem.reservoirTemperature_,
2137 newGasPvt,
2138 newOilPvt,
2139 newWaterPvt,
2140 oldFluidSystem.enableDissolvedGas_,
2141 oldFluidSystem.enableDissolvedGasInWater_,
2142 oldFluidSystem.enableVaporizedOil_,
2143 oldFluidSystem.enableVaporizedWater_,
2144 oldFluidSystem.enableConstantRs_,
2145 oldFluidSystem.enableDiffusion_,
2146 std::move(newReferenceDensity),
2147 std::move(newMolarMass),
2148 std::move(newDiffusionCoefficients),
2149 oldFluidSystem.phaseUsageInfo_,
2150 oldFluidSystem.isInitialized_,
2151 oldFluidSystem.useSaturatedTables_,
2152 oldFluidSystem.enthalpy_eq_energy_
2153 );
2154}
2155}
2156#endif // HAVE_CUDA
2157#endif // COMPILING_STATIC_FLUID_SYSTEM
2158} // namespace Opm
The base class for all fluid systems.
The class which specifies the default phase and component indices for the black-oil fluid system.
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
This macro generates a class HasMember_${MEMBER_NAME} which can be used for template specialization.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
A parameter cache which does nothing.
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, GPUContainerType > copy_to_gpu(const PiecewiseLinearTwoPhaseMaterialParams< TraitsT > &params)
Move a PiecewiseLinearTwoPhaseMaterialParams-object to the GPU.
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:285
PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ViewType > make_view(PiecewiseLinearTwoPhaseMaterialParams< TraitsT, ContainerType > &params)
this function is intented to make a GPU friendly view of the PiecewiseLinearTwoPhaseMaterialParams
Definition PiecewiseLinearTwoPhaseMaterialParams.hpp:312
Some templates to wrap the valgrind client request macros.
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
The base class for all fluid systems.
Definition BaseFluidSystem.hpp:43
The class which specifies the default phase and component indices for the black-oil fluid system.
Definition BlackOilDefaultFluidSystemIndices.hpp:39
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
Definition BrineCo2Pvt.hpp:80
This class represents the Pressure-Volume-Temperature relations of the gas phase for CO2.
Definition Co2GasPvt.hpp:75
static constexpr Scalar R
The ideal gas constant [J/(mol K)].
Definition Constants.hpp:47
Definition EclipseState.hpp:66
A fluid system which uses the black-oil model assumptions to calculate termodynamically meaningful qu...
Definition BlackOilFluidSystem_macrotemplate.hpp:74
STATIC_OR_DEVICE Scalar diffusionCoefficient(unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem_macrotemplate.hpp:1738
STATIC_OR_DEVICE unsigned numActivePhases() NOTHING_OR_CONST
Returns the number of active fluid phases (i.e., usually three).
Definition BlackOilFluidSystem_macrotemplate.hpp:434
STATIC_OR_DEVICE bool enableConstantRs() NOTHING_OR_CONST
Returns whether constant Rs tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:524
STATIC_OR_DEVICE LhsEval convertRvwToXgW(const LhsEval &Rvw, unsigned regionIdx) NOTHING_OR_CONST
Convert an water vaporization factor to the corresponding mass fraction of the water component in the...
Definition BlackOilFluidSystem_macrotemplate.hpp:1598
STATIC_OR_DEVICE bool useSaturatedTables() NOTHING_OR_CONST
Returns whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:540
STATIC_OR_NOTHING void initBegin(std::size_t numPvtRegions)
Begin the initialization of the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1847
STATIC_OR_DEVICE LhsEval saturatedDensity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Compute the density of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:708
STATIC_OR_DEVICE LhsEval convertRvToXgO(const LhsEval &Rv, unsigned regionIdx) NOTHING_OR_CONST
Convert an oil vaporization factor to the corresponding mass fraction of the oil component in the gas...
Definition BlackOilFluidSystem_macrotemplate.hpp:1584
STATIC_OR_DEVICE void initEnd()
Finish initializing the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1881
STATIC_OR_DEVICE LhsEval dewPointPressure(const FluidState &fluidState, unsigned regionIdx) NOTHING_OR_CONST
Returns the dew point pressure $P_d$ using the current Rv.
Definition BlackOilFluidSystem_macrotemplate.hpp:1459
STATIC_OR_DEVICE LhsEval convertXgOToRv(const LhsEval &XgO, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the oil component in the gas phase to the corresponding oil vaporization...
Definition BlackOilFluidSystem_macrotemplate.hpp:1529
STATIC_OR_DEVICE LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BlackOilFluidSystem_macrotemplate.hpp:583
STATIC_OR_DEVICE const OilPvt & oilPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the oil phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1696
STATIC_OR_DEVICE bool isIdealMixture(unsigned)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition BlackOilFluidSystem_macrotemplate.hpp:457
STATIC_OR_DEVICE void setEnableVaporizedOil(bool yesno)
Specify whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:285
STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx, const LhsEval &maxOilSaturation) NOTHING_OR_CONST
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1393
STATIC_OR_DEVICE bool enableVaporizedOil() NOTHING_OR_CONST
Returns whether the fluid system should consider that the oil component can dissolve in the gas phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:509
STATIC_OR_DEVICE Scalar reservoirTemperature(unsigned=0) NOTHING_OR_CONST
Set the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1714
STATIC_OR_DEVICE LhsEval saturatedVaporizationFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the water vaporization factor of saturated phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1367
STATIC_OR_DEVICE LhsEval convertXoGToRs(const LhsEval &XoG, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the gas component in the oil phase to the corresponding gas dissolution ...
Definition BlackOilFluidSystem_macrotemplate.hpp:1503
STATIC_OR_DEVICE LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx) NOTHING_OR_CONST
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem_macrotemplate.hpp:576
STATIC_OR_DEVICE LhsEval convertXgWToxgW(const LhsEval &XgW, unsigned regionIdx) NOTHING_OR_CONST
Convert a water mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1611
STATIC_OR_DEVICE void setWaterPvt(std::shared_ptr< WaterPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the water phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:346
STATIC_OR_DEVICE LhsEval convertXgWToRvw(const LhsEval &XgW, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the water component in the gas phase to the corresponding water vaporiza...
Definition BlackOilFluidSystem_macrotemplate.hpp:1542
static constexpr unsigned gasPhaseIdx
Index of the gas phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:395
static constexpr int oilCompIdx
Index of the oil component.
Definition BlackOilFluidSystem_macrotemplate.hpp:421
FLUIDSYSTEM_CLASSNAME(const FLUIDSYSTEM_CLASSNAME< Scalar, IndexTraits, StorageT > &other)
Default copy constructor.
Definition BlackOilFluidSystem_macrotemplate.hpp:160
STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState &fluidState, unsigned phaseIdx, unsigned compIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:959
STATIC_OR_DEVICE void setEnableVaporizedWater(bool yesno)
Specify whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:294
STATIC_OR_DEVICE LhsEval convertXoGToxoG(const LhsEval &XoG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mass fraction in the oil phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1635
STATIC_OR_DEVICE LhsEval density(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:600
STATIC_OR_DEVICE void setUseSaturatedTables(bool yesno)
Specify whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:328
STATIC_OR_DEVICE bool phaseIsActive(unsigned phaseIdx) NOTHING_OR_CONST
Returns whether a fluid phase is active.
Definition BlackOilFluidSystem_macrotemplate.hpp:438
STATIC_OR_DEVICE LhsEval convertxgOToXgO(const LhsEval &xgO, unsigned regionIdx) NOTHING_OR_CONST
Convert a oil mole fraction in the gas phase the corresponding mass fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1671
STATIC_OR_DEVICE LhsEval saturationPressure(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the saturation pressure of a given phase [Pa] depending on its composition.
Definition BlackOilFluidSystem_macrotemplate.hpp:1476
STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the dissolution factor of a saturated fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1424
STATIC_OR_DEVICE LhsEval inverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the formation volume factor of an "undersaturated" fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:819
STATIC_OR_DEVICE void setEnableConstantRs(bool yesno)
Specify whether the fluid system should use constant Rs tables.
Definition BlackOilFluidSystem_macrotemplate.hpp:312
STATIC_OR_DEVICE bool enableDiffusion() NOTHING_OR_CONST
Returns whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:532
static constexpr unsigned numComponents
Number of chemical species in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:418
STATIC_OR_DEVICE bool isCompressible(unsigned) NOTHING_OR_CONST
Returns true if and only if a fluid phase is assumed to be compressible.
Definition BlackOilFluidSystem_macrotemplate.hpp:465
STATIC_OR_DEVICE LhsEval convertXgOToxgO(const LhsEval &XgO, unsigned regionIdx) NOTHING_OR_CONST
Convert a oil mass fraction in the gas phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1659
STATIC_OR_DEVICE LhsEval convertxoGToXoG(const LhsEval &xoG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mole fraction in the oil phase the corresponding mass fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1647
STATIC_OR_DEVICE Scalar referenceDensity(unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the density of a fluid phase at surface pressure [kg/m^3].
Definition BlackOilFluidSystem_macrotemplate.hpp:548
STATIC_OR_DEVICE LhsEval convertXwGToxwG(const LhsEval &XwG, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas mass fraction in the water phase the corresponding mole fraction.
Definition BlackOilFluidSystem_macrotemplate.hpp:1623
STATIC_OR_DEVICE void setEnableDiffusion(bool yesno)
Specify whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:320
STATIC_OR_DEVICE void setOilPvt(std::shared_ptr< OilPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the oil phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:340
STATIC_OR_DEVICE std::size_t numRegions() NOTHING_OR_CONST
Returns the number of PVT regions which are considered.
Definition BlackOilFluidSystem_macrotemplate.hpp:481
STATIC_OR_NOTHING std::string_view phaseName(unsigned phaseIdx) NOTHING_OR_CONST
Return the human readable name of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:1911
STATIC_OR_DEVICE const PhaseUsageInfo< IndexTraits > & phaseUsage() NOTHING_OR_CONST
Returns a const reference to PhaseUsageInfo.
Definition BlackOilFluidSystem_macrotemplate.hpp:430
STATIC_OR_DEVICE bool enableDissolvedGasInWater() NOTHING_OR_CONST
Returns whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:500
STATIC_OR_DEVICE LhsEval bubblePointPressure(const FluidState &fluidState, unsigned regionIdx) NOTHING_OR_CONST
Returns the bubble point pressure $P_b$ using the current Rs.
Definition BlackOilFluidSystem_macrotemplate.hpp:1448
STATIC_OR_DEVICE LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx) NOTHING_OR_CONST
Calculate the density [kg/m^3] of a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:556
STATIC_OR_DEVICE const GasPvt & gasPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the gas phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1686
static constexpr unsigned numPhases
Number of fluid phases in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:388
STATIC_OR_DEVICE LhsEval viscosity(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BlackOilFluidSystem_macrotemplate.hpp:1081
STATIC_OR_DEVICE unsigned solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "primary" component of a phase (solvent)
Definition BlackOilFluidSystem_macrotemplate.hpp:1928
static constexpr unsigned waterPhaseIdx
Index of the water phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:391
static constexpr int gasCompIdx
Index of the gas component.
Definition BlackOilFluidSystem_macrotemplate.hpp:425
STATIC_OR_DEVICE bool enableDissolvedGas() NOTHING_OR_CONST
Returns whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:490
STATIC_OR_DEVICE void setReservoirTemperature(Scalar value) NOTHING_OR_CONST
Return the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1722
STATIC_OR_DEVICE unsigned soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "secondary" component of a phase (solute)
Definition BlackOilFluidSystem_macrotemplate.hpp:1945
static constexpr unsigned oilPhaseIdx
Index of the oil phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:393
STATIC_OR_DEVICE LhsEval enthalpy(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition BlackOilFluidSystem_macrotemplate.hpp:1346
STATIC_OR_NOTHING std::string_view componentName(unsigned compIdx) NOTHING_OR_CONST
Return the human readable name of a component.
Definition BlackOilFluidSystem_macrotemplate.hpp:1967
STATIC_OR_DEVICE const WaterPvt & waterPvt() NOTHING_OR_CONST
Return a reference to the low-level object which calculates the water phase quantities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1706
STATIC_OR_NOTHING Scalar surfaceTemperature
The temperature at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:401
STATIC_OR_DEVICE bool isIdealGas(unsigned) NOTHING_OR_CONST
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition BlackOilFluidSystem_macrotemplate.hpp:469
STATIC_OR_DEVICE LhsEval convertRsToXoG(const LhsEval &Rs, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the o...
Definition BlackOilFluidSystem_macrotemplate.hpp:1556
STATIC_OR_DEVICE void setGasPvt(std::shared_ptr< GasPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the gas phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:334
STATIC_OR_DEVICE void setEnableDissolvedGasInWater(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the water pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:303
STATIC_OR_DEVICE void setReferenceDensities(Scalar rhoOil, Scalar rhoWater, Scalar rhoGas, unsigned regionIdx)
Initialize the values of the reference densities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1870
STATIC_OR_DEVICE bool enableVaporizedWater() NOTHING_OR_CONST
Returns whether the fluid system should consider that the water component can dissolve in the gas pha...
Definition BlackOilFluidSystem_macrotemplate.hpp:518
STATIC_OR_DEVICE LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx, unsigned compIdx) NOTHING_OR_CONST
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition BlackOilFluidSystem_macrotemplate.hpp:1749
STATIC_OR_NOTHING void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the fluid system using an ECL deck object.
STATIC_OR_DEVICE LhsEval convertRswToXwG(const LhsEval &Rsw, unsigned regionIdx) NOTHING_OR_CONST
Convert a gas dissolution factor to the the corresponding mass fraction of the gas component in the w...
Definition BlackOilFluidSystem_macrotemplate.hpp:1570
STATIC_OR_DEVICE Scalar molarMass(unsigned compIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Return the molar mass of a component in [kg/mol].
Definition BlackOilFluidSystem_macrotemplate.hpp:453
STATIC_OR_NOTHING Scalar surfacePressure
The pressure at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:398
STATIC_OR_DEVICE void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Definition BlackOilFluidSystem_macrotemplate.hpp:1742
STATIC_OR_DEVICE void setEnableDissolvedGas(bool yesno)
Specify whether the fluid system should consider that the gas component can dissolve in the oil phase...
Definition BlackOilFluidSystem_macrotemplate.hpp:276
STATIC_OR_DEVICE bool isLiquid(unsigned phaseIdx) NOTHING_OR_CONST
Return whether a phase is liquid.
Definition BlackOilFluidSystem_macrotemplate.hpp:407
STATIC_OR_DEVICE LhsEval saturatedInverseFormationVolumeFactor(const FluidState &fluidState, unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
Returns the formation volume factor of a "saturated" fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:937
STATIC_OR_DEVICE LhsEval convertXwGToRsw(const LhsEval &XwG, unsigned regionIdx) NOTHING_OR_CONST
Convert the mass fraction of the gas component in the water phase to the corresponding gas dissolutio...
Definition BlackOilFluidSystem_macrotemplate.hpp:1516
STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx, unsigned compIdx) NOTHING_OR_CONST
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:563
static constexpr int waterCompIdx
Index of the water component.
Definition BlackOilFluidSystem_macrotemplate.hpp:423
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
Definition GasPvtMultiplexer.hpp:109
Null object for oil PVT calculations.
Definition NullOilPvt.hpp:35
A parameter cache which does nothing.
Definition NullParameterCache.hpp:40
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Definition OilPvtMultiplexer.hpp:110
Definition PhaseUsageInfo.hpp:42
Definition Runspec.hpp:46
Definition Schedule.hpp:101
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition WaterPvtMultiplexer.hpp:88
Convience header to include the gpuistl headers if HAVE_CUDA is defined.
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30