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/TimingMacros.hpp>
36#include <opm/common/utility/VectorWithDefaultAllocator.hpp>
37#include <opm/common/utility/gpuDecorators.hpp>
39
42
48#include <opm/material/fluidsystems/PhaseUsageInfo.hpp>
49
50#include <array>
51#include <cstddef>
52#include <memory>
53#include <stdexcept>
54#include <string>
55#include <string_view>
56#include <type_traits>
57#include <vector>
58
59namespace Opm {
60
61#if HAVE_ECL_INPUT
62class EclipseState;
63class Schedule;
64#endif
65
66
67
68
75template <class Scalar, class IndexTraits = BlackOilDefaultFluidSystemIndices,
76 template<typename> typename Storage = VectorWithDefaultAllocator>
77class FLUIDSYSTEM_CLASSNAME : public BaseFluidSystem<Scalar, FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, Storage> >
78{
79 using ThisType = FLUIDSYSTEM_CLASSNAME;
80
81public:
82
83 // The logic here is the following: We test if we are instantiating on the CPU
84 // std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar> == true
85 // and if so, we use the multiplexer classes, if not, we use the concrete types.
86 using GasPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
89 using OilPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
91 NullOilPvt<Scalar>>; // Currently do not support on GPU
92 using WaterPvt = std::conditional_t<std::is_same_v<Storage<Scalar>, VectorWithDefaultAllocator<Scalar>>,
95 using IndexTraitsType = IndexTraits;
96
97 #ifdef COMPILING_STATIC_FLUID_SYSTEM
99 template <class EvaluationT>
100 struct ParameterCache : public NullParameterCache<EvaluationT>
101 {
102 using Evaluation = EvaluationT;
103
104 public:
105 explicit ParameterCache(unsigned regionIdx = 0)
106 : regionIdx_(regionIdx)
107 {
108 }
109
117 template <class OtherCache>
118 void assignPersistentData(const OtherCache& other)
119 {
120 regionIdx_ = other.regionIndex();
121 }
122
130 unsigned regionIndex() const
131 { return regionIdx_; }
132
140 void setRegionIndex(unsigned val)
141 { regionIdx_ = val; }
142
143 private:
144 unsigned regionIdx_;
145 };
146
147 #else
148
149 // We want to use the exact same ParameterCache class for both the static and nonstatic versions
150 template<class EvaluationT>
151 using ParameterCache =
152 typename FLUIDSYSTEM_CLASSNAME_STATIC<Scalar,
153 IndexTraits,
154 Storage>::template ParameterCache<EvaluationT>;
155 #endif
156
157 #ifndef COMPILING_STATIC_FLUID_SYSTEM
163 template<template<typename> typename StorageT>
167 , reservoirTemperature_(other.reservoirTemperature_)
168 , gasPvt_(other.gasPvt_)
169 , oilPvt_(other.oilPvt_)
170 , waterPvt_(other.waterPvt_)
171 , enableDissolvedGas_(other.enableDissolvedGas_)
172 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
173 , enableVaporizedOil_(other.enableVaporizedOil_)
174 , enableVaporizedWater_(other.enableVaporizedWater_)
175 , enableDiffusion_(other.enableDiffusion_)
176 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
177 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
178 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
179 , phaseUsageInfo_(other.phaseUsageInfo_)
180 , isInitialized_(other.isInitialized_)
181 , useSaturatedTables_(other.useSaturatedTables_)
182 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
183 {
184 }
185
186 FLUIDSYSTEM_CLASSNAME(Scalar _surfacePressure_,
187 Scalar _surfaceTemperature_,
188 Scalar _reservoirTemperature_,
189 const GasPvt& _gasPvt_,
190 const OilPvt& _oilPvt_,
191 const WaterPvt& _waterPvt_,
192 bool _enableDissolvedGas_,
193 bool _enableDissolvedGasInWater_,
194 bool _enableVaporizedOil_,
195 bool _enableVaporizedWater_,
196 bool _enableDiffusion_,
197 Storage<std::array<Scalar, 3>>&& _referenceDensity_,
198 Storage<std::array<Scalar, 3>>&& _molarMass_,
199 Storage<std::array<Scalar, 3 * 3>>&& _diffusionCoefficients_,
200 const PhaseUsageInfo<IndexTraits>& _phaseUsageInfo_,
201 bool _isInitialized_,
202 bool _useSaturatedTables_,
203 bool _enthalpy_eq_energy_)
204 : surfacePressure(_surfacePressure_)
205 , surfaceTemperature(_surfaceTemperature_)
206 , reservoirTemperature_(_reservoirTemperature_)
207 , gasPvt_(_gasPvt_)
208 , oilPvt_(_oilPvt_)
209 , waterPvt_(_waterPvt_)
210 , enableDissolvedGas_(_enableDissolvedGas_)
211 , enableDissolvedGasInWater_(_enableDissolvedGasInWater_)
212 , enableVaporizedOil_(_enableVaporizedOil_)
213 , enableVaporizedWater_(_enableVaporizedWater_)
214 , enableDiffusion_(_enableDiffusion_)
215 , referenceDensity_(std::move(_referenceDensity_))
216 , molarMass_(std::move(_molarMass_))
217 , diffusionCoefficients_(std::move(_diffusionCoefficients_))
218 , phaseUsageInfo_(_phaseUsageInfo_)
219 , isInitialized_(_isInitialized_)
220 , useSaturatedTables_(_useSaturatedTables_)
221 , enthalpy_eq_energy_(_enthalpy_eq_energy_)
222 {
223 }
224
225 #if HAVE_CUDA
226 template <class ScalarT, class IndexTraitsT>
227 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>
228 gpuistl::copy_to_gpu(const FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT>& oldFluidSystem);
229
230 template <class ScalarT, class IndexTraitsT>
231 friend FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuView>
232 gpuistl::make_view(FLUIDSYSTEM_CLASSNAME<ScalarT, IndexTraitsT, gpuistl::GpuBuffer>& oldFluidSystem);
233
234 #endif // HAVE_CUDA
235 #endif // COMPILING_STATIC_FLUID_SYSTEM
236
237 #ifdef COMPILING_STATIC_FLUID_SYSTEM
245 template<template<typename> typename StorageT = VectorWithDefaultAllocator>
246 static FLUIDSYSTEM_CLASSNAME_NONSTATIC<Scalar, IndexTraits, StorageT>& getNonStaticInstance()
247 {
248 static FLUIDSYSTEM_CLASSNAME_NONSTATIC<Scalar, IndexTraits, StorageT> instance{FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, Storage>()};
249 return instance;
250
251 }
252 #endif
253
254 /****************************************
255 * Initialization
256 ****************************************/
257#if HAVE_ECL_INPUT
261 STATIC_OR_NOTHING void initFromState(const EclipseState& eclState, const Schedule& schedule);
262#endif // HAVE_ECL_INPUT
263
272 STATIC_OR_NOTHING void initBegin(std::size_t numPvtRegions);
273
280 STATIC_OR_DEVICE void setEnableDissolvedGas(bool yesno)
281 { enableDissolvedGas_ = yesno; }
282
289 STATIC_OR_DEVICE void setEnableVaporizedOil(bool yesno)
290 { enableVaporizedOil_ = yesno; }
291
298 STATIC_OR_DEVICE void setEnableVaporizedWater(bool yesno)
299 { enableVaporizedWater_ = yesno; }
300
307 STATIC_OR_DEVICE void setEnableDissolvedGasInWater(bool yesno)
308 { enableDissolvedGasInWater_ = yesno; }
309
314 STATIC_OR_DEVICE void setEnableDiffusion(bool yesno)
315 { enableDiffusion_ = yesno; }
316
322 STATIC_OR_DEVICE void setUseSaturatedTables(bool yesno)
323 { useSaturatedTables_ = yesno; }
324
328 STATIC_OR_DEVICE void setGasPvt(std::shared_ptr<GasPvt> pvtObj)
329 { gasPvt_ = *pvtObj; }
330
334 STATIC_OR_DEVICE void setOilPvt(std::shared_ptr<OilPvt> pvtObj)
335 { oilPvt_ = *pvtObj; }
336
340 STATIC_OR_DEVICE void setWaterPvt(std::shared_ptr<WaterPvt> pvtObj)
341 { waterPvt_ = *pvtObj; }
342
343 STATIC_OR_DEVICE void setVapPars(const Scalar par1, const Scalar par2)
344 {
345 if (gasPvt_.isActive()) {
346 gasPvt_.setVapPars(par1, par2);
347 }
348 if (oilPvt_.isActive()) {
349 oilPvt_.setVapPars(par1, par2);
350 }
351 if (waterPvt_.isActive()) {
352 waterPvt_.setVapPars(par1, par2);
353 }
354 }
355
364 STATIC_OR_DEVICE void setReferenceDensities(Scalar rhoOil,
365 Scalar rhoWater,
366 Scalar rhoGas,
367 unsigned regionIdx);
368
372 STATIC_OR_DEVICE void initEnd();
373
374 STATIC_OR_DEVICE bool isInitialized() NOTHING_OR_CONST
375 { return isInitialized_; }
376
377 /****************************************
378 * Generic phase properties
379 ****************************************/
380
382 static constexpr unsigned numPhases = IndexTraits::numPhases;
383
385 static constexpr unsigned waterPhaseIdx = IndexTraits::waterPhaseIdx;
387 static constexpr unsigned oilPhaseIdx = IndexTraits::oilPhaseIdx;
389 static constexpr unsigned gasPhaseIdx = IndexTraits::gasPhaseIdx;
390
392 STATIC_OR_NOTHING Scalar surfacePressure;
393
395 STATIC_OR_NOTHING Scalar surfaceTemperature;
396
398 STATIC_OR_NOTHING std::string_view phaseName(unsigned phaseIdx) NOTHING_OR_CONST;
399
401 STATIC_OR_DEVICE bool isLiquid(unsigned phaseIdx) NOTHING_OR_CONST
402 {
403 assert(phaseIdx < numPhases);
404 return phaseIdx != gasPhaseIdx;
405 }
406
407 /****************************************
408 * Generic component related properties
409 ****************************************/
410
412 static constexpr unsigned numComponents = IndexTraits::numComponents;
413
415 static constexpr int oilCompIdx = IndexTraits::oilCompIdx;
417 static constexpr int waterCompIdx = IndexTraits::waterCompIdx;
419 static constexpr int gasCompIdx = IndexTraits::gasCompIdx;
420
421public:
422
424 STATIC_OR_DEVICE const PhaseUsageInfo<IndexTraits>& phaseUsage() NOTHING_OR_CONST
425 { return phaseUsageInfo_; }
426
428 STATIC_OR_DEVICE unsigned numActivePhases() NOTHING_OR_CONST
429 { return phaseUsageInfo_.numActivePhases(); }
430
432 STATIC_OR_DEVICE bool phaseIsActive(unsigned phaseIdx) NOTHING_OR_CONST
433 {
434 return phaseUsageInfo_.phaseIsActive(phaseIdx);
435 }
436
438 STATIC_OR_DEVICE unsigned solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST;
439
441 STATIC_OR_DEVICE unsigned soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST;
442
444 STATIC_OR_NOTHING std::string_view componentName(unsigned compIdx) NOTHING_OR_CONST;
445
447 STATIC_OR_DEVICE Scalar molarMass(unsigned compIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
448 { return molarMass_[regionIdx][compIdx]; }
449
451 STATIC_OR_DEVICE bool isIdealMixture(unsigned /*phaseIdx*/)
452 {
453 // fugacity coefficients are only pressure dependent -> we
454 // have an ideal mixture
455 return true;
456 }
457
459 STATIC_OR_DEVICE bool isCompressible(unsigned /*phaseIdx*/) NOTHING_OR_CONST
460 { return true; /* all phases are compressible */ }
461
463 STATIC_OR_DEVICE bool isIdealGas(unsigned /*phaseIdx*/) NOTHING_OR_CONST
464 { return false; }
465
466
467 /****************************************
468 * Black-oil specific properties
469 ****************************************/
475 STATIC_OR_DEVICE std::size_t numRegions() NOTHING_OR_CONST
476 { return molarMass_.size(); }
477
484 STATIC_OR_DEVICE bool enableDissolvedGas() NOTHING_OR_CONST
485 { return enableDissolvedGas_; }
486
487
494 STATIC_OR_DEVICE bool enableDissolvedGasInWater() NOTHING_OR_CONST
495 { return enableDissolvedGasInWater_; }
496
503 STATIC_OR_DEVICE bool enableVaporizedOil() NOTHING_OR_CONST
504 { return enableVaporizedOil_; }
505
512 STATIC_OR_DEVICE bool enableVaporizedWater() NOTHING_OR_CONST
513 { return enableVaporizedWater_; }
514
520 STATIC_OR_DEVICE bool enableDiffusion() NOTHING_OR_CONST
521 { return enableDiffusion_; }
522
528 STATIC_OR_DEVICE bool useSaturatedTables() NOTHING_OR_CONST
529 { return useSaturatedTables_; }
530
536 STATIC_OR_DEVICE Scalar referenceDensity(unsigned phaseIdx, unsigned regionIdx) NOTHING_OR_CONST
537 { return referenceDensity_[regionIdx][phaseIdx]; }
538
539 /****************************************
540 * thermodynamic quantities (generic version)
541 ****************************************/
543 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
544 STATIC_OR_DEVICE LhsEval density(const FluidState& fluidState,
545 const ParameterCache<ParamCacheEval>& paramCache,
546 unsigned phaseIdx) NOTHING_OR_CONST
547 { return density<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
548
550 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
551 STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState& fluidState,
552 const ParameterCache<ParamCacheEval>& paramCache,
553 unsigned phaseIdx,
554 unsigned compIdx) NOTHING_OR_CONST
555 {
557 phaseIdx,
558 compIdx,
559 paramCache.regionIndex());
560 }
561
563 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
564 STATIC_OR_DEVICE LhsEval viscosity(const FluidState& fluidState,
565 const ParameterCache<ParamCacheEval>& paramCache,
566 unsigned phaseIdx) NOTHING_OR_CONST
567 { return viscosity<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
568
570 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
571 STATIC_OR_DEVICE LhsEval enthalpy(const FluidState& fluidState,
572 const ParameterCache<ParamCacheEval>& paramCache,
573 unsigned phaseIdx)
574 { return enthalpy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
575
576 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
577 STATIC_OR_DEVICE LhsEval internalEnergy(const FluidState& fluidState,
578 const ParameterCache<ParamCacheEval>& paramCache,
579 unsigned phaseIdx) NOTHING_OR_CONST
580 { return internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, paramCache.regionIndex()); }
581
582 /****************************************
583 * thermodynamic quantities (black-oil specific version: Note that the PVT region
584 * index is explicitly passed instead of a parameter cache object)
585 ****************************************/
587 template <class FluidState, class LhsEval = typename FluidState::Scalar>
588 STATIC_OR_DEVICE LhsEval density(const FluidState& fluidState,
589 unsigned phaseIdx,
590 unsigned regionIdx) NOTHING_OR_CONST
591 {
592 assert(phaseIdx <= numPhases);
593 assert(regionIdx <= numRegions());
594
595 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
596 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
597 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
598
599 switch (phaseIdx) {
600 case oilPhaseIdx: {
601 if (enableDissolvedGas()) {
602 // miscible oil
603 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
604 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
605
606 return
607 bo*referenceDensity(oilPhaseIdx, regionIdx)
608 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
609 }
610
611 // immiscible oil
612 const LhsEval Rs(0.0);
613 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
614
615 return referenceDensity(phaseIdx, regionIdx)*bo;
616 }
617
618 case gasPhaseIdx: {
620 // gas containing vaporized oil and vaporized water
621 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
622 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
623 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
624
625 return
626 bg*referenceDensity(gasPhaseIdx, regionIdx)
627 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
628 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
629 }
630 if (enableVaporizedOil()) {
631 // miscible gas
632 const LhsEval Rvw(0.0);
633 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
634 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
635
636 return
637 bg*referenceDensity(gasPhaseIdx, regionIdx)
638 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
639 }
640 if (enableVaporizedWater()) {
641 // gas containing vaporized water
642 const LhsEval Rv(0.0);
643 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
644 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
645
646 return
647 bg*referenceDensity(gasPhaseIdx, regionIdx)
648 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
649 }
650
651 // immiscible gas
652 const LhsEval Rv(0.0);
653 const LhsEval Rvw(0.0);
654 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
655 return bg*referenceDensity(phaseIdx, regionIdx);
656 }
657
658 case waterPhaseIdx:
660 // gas miscible in water
661 const LhsEval& Rsw =BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
662 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
663 return
664 bw*referenceDensity(waterPhaseIdx, regionIdx)
665 + Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
666 }
667 const LhsEval Rsw(0.0);
668 return
670 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
671 }
672
673 throw std::logic_error("Unhandled phase index " + std::to_string(phaseIdx));
674 }
675
685 template <class FluidState, class LhsEval = typename FluidState::Scalar>
686 STATIC_OR_DEVICE LhsEval saturatedDensity(const FluidState& fluidState,
687 unsigned phaseIdx,
688 unsigned regionIdx) NOTHING_OR_CONST
689 {
690 assert(phaseIdx <= numPhases);
691 assert(regionIdx <= numRegions());
692
693 const auto& p = fluidState.pressure(phaseIdx);
694 const auto& T = fluidState.temperature(phaseIdx);
695
696 switch (phaseIdx) {
697 case oilPhaseIdx: {
698 if (enableDissolvedGas()) {
699 // miscible oil
700 const LhsEval& Rs = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, oilPhaseIdx, regionIdx);
701 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
702
703 return
704 bo*referenceDensity(oilPhaseIdx, regionIdx)
705 + Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
706 }
707
708 // immiscible oil
709 const LhsEval Rs(0.0);
710 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
711 return referenceDensity(phaseIdx, regionIdx)*bo;
712 }
713
714 case gasPhaseIdx: {
716 // gas containing vaporized oil and vaporized water
717 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
718 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
719 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
720
721 return
722 bg*referenceDensity(gasPhaseIdx, regionIdx)
723 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
724 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx) ;
725 }
726
727 if (enableVaporizedOil()) {
728 // miscible gas
729 const LhsEval Rvw(0.0);
730 const LhsEval& Rv = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
731 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
732
733 return
734 bg*referenceDensity(gasPhaseIdx, regionIdx)
735 + Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
736 }
737
738 if (enableVaporizedWater()) {
739 // gas containing vaporized water
740 const LhsEval Rv(0.0);
741 const LhsEval& Rvw = saturatedVaporizationFactor<FluidState, LhsEval>(fluidState, gasPhaseIdx, regionIdx);
742 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
743
744 return
745 bg*referenceDensity(gasPhaseIdx, regionIdx)
746 + Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
747 }
748
749 // immiscible gas
750 const LhsEval Rv(0.0);
751 const LhsEval Rvw(0.0);
752 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
753
754 return referenceDensity(phaseIdx, regionIdx)*bg;
755
756 }
757
758 case waterPhaseIdx:
759 {
761 // miscible in water
762 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
763 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, waterPhaseIdx, regionIdx);
764 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
765 return
766 bw*referenceDensity(waterPhaseIdx, regionIdx)
767 + Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
768 }
769 return
772 }
773 }
774
775 throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
776 }
777
786 template <class FluidState, class LhsEval = typename FluidState::Scalar>
787 STATIC_OR_DEVICE LhsEval inverseFormationVolumeFactor(const FluidState& fluidState,
788 unsigned phaseIdx,
789 unsigned regionIdx) NOTHING_OR_CONST
790 {
791 OPM_TIMEBLOCK_LOCAL(inverseFormationVolumeFactor, Subsystem::PvtProps);
792 assert(phaseIdx <= numPhases);
793 assert(regionIdx <= numRegions());
794
795 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
796 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
797
798 switch (phaseIdx) {
799 case oilPhaseIdx: {
800 if (enableDissolvedGas()) {
801 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
802 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
803 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
804 {
805 return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
806 } else {
807 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
808 }
809 }
810
811 const LhsEval Rs(0.0);
812 return oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
813 }
814 case gasPhaseIdx: {
816 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
817 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
818 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
819 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
820 && fluidState.saturation(oilPhaseIdx) > 0.0
821 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
822 {
823 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
824 } else {
825 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
826 }
827 }
828
829 if (enableVaporizedOil()) {
830 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
831 if (useSaturatedTables() && fluidState.saturation(oilPhaseIdx) > 0.0
832 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
833 {
834 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
835 } else {
836 const LhsEval Rvw(0.0);
837 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
838 }
839 }
840
841 if (enableVaporizedWater()) {
842 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
843 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
844 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
845 {
846 return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
847 } else {
848 const LhsEval Rv(0.0);
849 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
850 }
851 }
852
853 const LhsEval Rv(0.0);
854 const LhsEval Rvw(0.0);
855 return gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
856 }
857 case waterPhaseIdx:
858 {
859 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
861 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
862 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
863 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
864 {
865 return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
866 } else {
867 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
868 }
869 }
870 const LhsEval Rsw(0.0);
871 return waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
872 }
873 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
874 }
875 }
876
877 template <class FluidState, class LhsEval = typename FluidState::Scalar>
878 STATIC_OR_DEVICE std::pair<LhsEval, LhsEval>
879 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState,
880 unsigned phaseIdx,
881 unsigned regionIdx)
882 {
883 switch (phaseIdx) {
884 case oilPhaseIdx:
885 return oilPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
886 case gasPhaseIdx:
887 return gasPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
888 case waterPhaseIdx:
889 return waterPvt_.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
890 default:
891 throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
892 }
893 }
894
904 template <class FluidState, class LhsEval = typename FluidState::Scalar>
905 STATIC_OR_DEVICE LhsEval saturatedInverseFormationVolumeFactor(const FluidState& fluidState,
906 unsigned phaseIdx,
907 unsigned regionIdx) NOTHING_OR_CONST
908 {
909 OPM_TIMEBLOCK_LOCAL(saturatedInverseFormationVolumeFactor, Subsystem::PvtProps);
910 assert(phaseIdx <= numPhases);
911 assert(regionIdx <= numRegions());
912
913 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
914 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
915 const auto& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
916
917 switch (phaseIdx) {
918 case oilPhaseIdx: return oilPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
919 case gasPhaseIdx: return gasPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p);
920 case waterPhaseIdx: return waterPvt_.saturatedInverseFormationVolumeFactor(regionIdx, T, p, saltConcentration);
921 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
922 }
923 }
924
926 template <class FluidState, class LhsEval = typename FluidState::Scalar>
927 STATIC_OR_DEVICE LhsEval fugacityCoefficient(const FluidState& fluidState,
928 unsigned phaseIdx,
929 unsigned compIdx,
930 unsigned regionIdx) NOTHING_OR_CONST
931 {
932 assert(phaseIdx <= numPhases);
933 assert(compIdx <= numComponents);
934 assert(regionIdx <= numRegions());
935
936 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
937 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
938
939 // for the fugacity coefficient of the oil component in the oil phase, we use
940 // some pseudo-realistic value for the vapor pressure to ease physical
941 // interpretation of the results
942 const LhsEval phi_oO = 20e3/p;
943
944 // for the gas component in the gas phase, assume it to be an ideal gas
945 constexpr const Scalar phi_gG = 1.0;
946
947 // for the fugacity coefficient of the water component in the water phase, we use
948 // the same approach as for the oil component in the oil phase
949 const LhsEval phi_wW = 30e3/p;
950
951 switch (phaseIdx) {
952 case gasPhaseIdx: // fugacity coefficients for all components in the gas phase
953 switch (compIdx) {
954 case gasCompIdx:
955 return phi_gG;
956
957 // for the oil component, we calculate the Rv value for saturated gas and Rs
958 // for saturated oil, and compute the fugacity coefficients at the
959 // equilibrium. for this, we assume that in equilibrium the fugacities of the
960 // oil component is the same in both phases.
961 case oilCompIdx: {
962 if (!enableVaporizedOil())
963 // if there's no vaporized oil, the gas phase is assumed to be
964 // immiscible with the oil component
965 return phi_gG*1e6;
966
967 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
968 const auto& X_gOSat = convertRvToXgO(R_vSat, regionIdx);
969 const auto& x_gOSat = convertXgOToxgO(X_gOSat, regionIdx);
970
971 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
972 const auto& X_oGSat = convertRsToXoG(R_sSat, regionIdx);
973 const auto& x_oGSat = convertXoGToxoG(X_oGSat, regionIdx);
974 const auto& x_oOSat = 1.0 - x_oGSat;
975
976 const auto& p_o = decay<LhsEval>(fluidState.pressure(oilPhaseIdx));
977 const auto& p_g = decay<LhsEval>(fluidState.pressure(gasPhaseIdx));
978
979 return phi_oO*p_o*x_oOSat / (p_g*x_gOSat);
980 }
981
982 case waterCompIdx:
983 // the water component is assumed to be never miscible with the gas phase
984 return phi_gG*1e6;
985
986 default:
987 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
988 }
989
990 case oilPhaseIdx: // fugacity coefficients for all components in the oil phase
991 switch (compIdx) {
992 case oilCompIdx:
993 return phi_oO;
994
995 // for the oil and water components, we proceed analogous to the gas and
996 // water components in the gas phase
997 case gasCompIdx: {
998 if (!enableDissolvedGas())
999 // if there's no dissolved gas, the oil phase is assumed to be
1000 // immiscible with the gas component
1001 return phi_oO*1e6;
1002
1003 const auto& R_vSat = gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1004 const auto& X_gOSat = convertRvToXgO(R_vSat, regionIdx);
1005 const auto& x_gOSat = convertXgOToxgO(X_gOSat, regionIdx);
1006 const auto& x_gGSat = 1.0 - x_gOSat;
1007
1008 const auto& R_sSat = oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1009 const auto& X_oGSat = convertRsToXoG(R_sSat, regionIdx);
1010 const auto& x_oGSat = convertXoGToxoG(X_oGSat, regionIdx);
1011
1012 const auto& p_o = decay<LhsEval>(fluidState.pressure(oilPhaseIdx));
1013 const auto& p_g = decay<LhsEval>(fluidState.pressure(gasPhaseIdx));
1014
1015 return phi_gG*p_g*x_gGSat / (p_o*x_oGSat);
1016 }
1017
1018 case waterCompIdx:
1019 return phi_oO*1e6;
1020
1021 default:
1022 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
1023 }
1024
1025 case waterPhaseIdx: // fugacity coefficients for all components in the water phase
1026 // the water phase fugacity coefficients are pretty simple: because the water
1027 // phase is assumed to consist entirely from the water component, we just
1028 // need to make sure that the fugacity coefficients for the other components
1029 // are a few orders of magnitude larger than the one of the water
1030 // component. (i.e., the affinity of the gas and oil components to the water
1031 // phase is lower by a few orders of magnitude)
1032 switch (compIdx) {
1033 case waterCompIdx: return phi_wW;
1034 case oilCompIdx: return 1.1e6*phi_wW;
1035 case gasCompIdx: return 1e6*phi_wW;
1036 default:
1037 throw std::logic_error("Invalid component index "+std::to_string(compIdx));
1038 }
1039
1040 default:
1041 throw std::logic_error("Invalid phase index "+std::to_string(phaseIdx));
1042 }
1043
1044 throw std::logic_error("Unhandled phase or component index");
1045 }
1046
1048 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1049 STATIC_OR_DEVICE LhsEval viscosity(const FluidState& fluidState,
1050 unsigned phaseIdx,
1051 unsigned regionIdx) NOTHING_OR_CONST
1052 {
1053 OPM_TIMEBLOCK_LOCAL(viscosity, Subsystem::PvtProps);
1054 assert(phaseIdx <= numPhases);
1055 assert(regionIdx <= numRegions());
1056
1057 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1058 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1059
1060 switch (phaseIdx) {
1061 case oilPhaseIdx: {
1062 if (enableDissolvedGas()) {
1063 const auto& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1064 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
1065 && Rs >= (1.0 - 1e-10)*oilPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p)))
1066 {
1067 return oilPvt_.saturatedViscosity(regionIdx, T, p);
1068 } else {
1069 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1070 }
1071 }
1072
1073 const LhsEval Rs(0.0);
1074 return oilPvt_.viscosity(regionIdx, T, p, Rs);
1075 }
1076
1077 case gasPhaseIdx: {
1079 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1080 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1081 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
1082 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p))
1083 && fluidState.saturation(oilPhaseIdx) > 0.0
1084 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1085 {
1086 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1087 } else {
1088 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1089 }
1090 }
1091 if (enableVaporizedOil()) {
1092 const auto& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1093 if (useSaturatedTables() && fluidState.saturation(oilPhaseIdx) > 0.0
1094 && Rv >= (1.0 - 1e-10)*gasPvt_.saturatedOilVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1095 {
1096 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1097 } else {
1098 const LhsEval Rvw(0.0);
1099 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1100 }
1101 }
1102 if (enableVaporizedWater()) {
1103 const auto& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1104 if (useSaturatedTables() && fluidState.saturation(waterPhaseIdx) > 0.0
1105 && Rvw >= (1.0 - 1e-10)*gasPvt_.saturatedWaterVaporizationFactor(regionIdx, scalarValue(T), scalarValue(p)))
1106 {
1107 return gasPvt_.saturatedViscosity(regionIdx, T, p);
1108 } else {
1109 const LhsEval Rv(0.0);
1110 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1111 }
1112 }
1113
1114 const LhsEval Rv(0.0);
1115 const LhsEval Rvw(0.0);
1116 return gasPvt_.viscosity(regionIdx, T, p, Rv, Rvw);
1117 }
1118
1119 case waterPhaseIdx:
1120 {
1121 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1123 const auto& Rsw = BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1124 if (useSaturatedTables() && fluidState.saturation(gasPhaseIdx) > 0.0
1125 && Rsw >= (1.0 - 1e-10)*waterPvt_.saturatedGasDissolutionFactor(regionIdx, scalarValue(T), scalarValue(p), scalarValue(saltConcentration)))
1126 {
1127 return waterPvt_.saturatedViscosity(regionIdx, T, p, saltConcentration);
1128 } else {
1129 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1130 }
1131 }
1132 const LhsEval Rsw(0.0);
1133 return waterPvt_.viscosity(regionIdx, T, p, Rsw, saltConcentration);
1134 }
1135 }
1136
1137 throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1138 }
1139
1140 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1141 STATIC_OR_DEVICE LhsEval internalEnergy(const FluidState& fluidState,
1142 const unsigned phaseIdx,
1143 const unsigned regionIdx) NOTHING_OR_CONST
1144 {
1145 const auto p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1146 const auto T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1147
1148 switch (phaseIdx) {
1149 case oilPhaseIdx:
1150 if (!oilPvt_.mixingEnergy()) {
1151 return oilPvt_.internalEnergy
1152 (regionIdx, T, p,
1153 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1154 }
1155 break;
1156
1157 case waterPhaseIdx:
1158 if (!waterPvt_.mixingEnergy()) {
1159 return waterPvt_.internalEnergy
1160 (regionIdx, T, p,
1161 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1162 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1163 }
1164 break;
1165
1166 case gasPhaseIdx:
1167 if (!gasPvt_.mixingEnergy()) {
1168 return gasPvt_.internalEnergy
1169 (regionIdx, T, p,
1170 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1171 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1172 }
1173 break;
1174
1175 default:
1176 throw std::logic_error {
1177 "Phase index " + std::to_string(phaseIdx) + " does not support internal energy"
1178 };
1179 }
1180
1181 return internalMixingTotalEnergy<FluidState,LhsEval>(fluidState, phaseIdx, regionIdx)
1182 / density<FluidState,LhsEval>(fluidState, phaseIdx, regionIdx);
1183 }
1184
1185
1186 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1187 STATIC_OR_DEVICE LhsEval internalMixingTotalEnergy(const FluidState& fluidState,
1188 unsigned phaseIdx,
1189 unsigned regionIdx) NOTHING_OR_CONST
1190 {
1191 assert(phaseIdx <= numPhases);
1192 assert(regionIdx <= numRegions());
1193 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1194 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1195 const LhsEval& saltConcentration = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
1196 // to avoid putting all thermal into the interface of the multiplexer
1197 switch (phaseIdx) {
1198 case oilPhaseIdx: {
1199 auto oilEnergy = oilPvt_.internalEnergy(regionIdx, T, p,
1200 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1201 assert(oilPvt_.mixingEnergy());
1202 //mixing energy adsed
1203 if (enableDissolvedGas()) {
1204 // miscible oil
1205 const LhsEval& Rs = BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1206 const LhsEval& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1207 const auto& gasEnergy =
1208 gasPvt_.internalEnergy(regionIdx, T, p,
1209 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1210 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1211 const auto hVapG = gasPvt_.hVap(regionIdx);// pressure correction ? assume equal to energy change
1212 return
1213 oilEnergy*bo*referenceDensity(oilPhaseIdx, regionIdx)
1214 + (gasEnergy-hVapG)*Rs*bo*referenceDensity(gasPhaseIdx, regionIdx);
1215 }
1216
1217 // immiscible oil
1218 const LhsEval Rs(0.0);
1219 const auto& bo = oilPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rs);
1220
1221 return oilEnergy*referenceDensity(phaseIdx, regionIdx)*bo;
1222 }
1223
1224 case gasPhaseIdx: {
1225 const auto& gasEnergy =
1226 gasPvt_.internalEnergy(regionIdx, T, p,
1227 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1228 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1229 assert(gasPvt_.mixingEnergy());
1231 const auto& oilEnergy =
1232 oilPvt_.internalEnergy(regionIdx, T, p,
1233 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1234 const auto waterEnergy =
1235 waterPvt_.internalEnergy(regionIdx, T, p,
1236 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1237 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1238 // gas containing vaporized oil and vaporized water
1239 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1240 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1241 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1242 const auto hVapO = oilPvt_.hVap(regionIdx);
1243 const auto hVapW = waterPvt_.hVap(regionIdx);
1244 return
1245 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1246 + (oilEnergy+hVapO)*Rv*bg*referenceDensity(oilPhaseIdx, regionIdx)
1247 + (waterEnergy+hVapW)*Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
1248 }
1249 if (enableVaporizedOil()) {
1250 const auto& oilEnergy =
1251 oilPvt_.internalEnergy(regionIdx, T, p,
1252 BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1253 // miscible gas
1254 const LhsEval Rvw(0.0);
1255 const LhsEval& Rv = BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1256 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1257 const auto hVapO = oilPvt_.hVap(regionIdx);
1258 return
1259 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1260 + (oilEnergy+hVapO)*Rv*bg*referenceDensity(oilPhaseIdx, regionIdx);
1261 }
1262 if (enableVaporizedWater()) {
1263 // gas containing vaporized water
1264 const LhsEval Rv(0.0);
1265 const LhsEval& Rvw = BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx);
1266 const LhsEval& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1267 const auto waterEnergy =
1268 waterPvt_.internalEnergy(regionIdx, T, p,
1269 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1270 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1271 const auto hVapW = waterPvt_.hVap(regionIdx);
1272 return
1273 gasEnergy*bg*referenceDensity(gasPhaseIdx, regionIdx)
1274 + (waterEnergy+hVapW)*Rvw*bg*referenceDensity(waterPhaseIdx, regionIdx);
1275 }
1276
1277 // immiscible gas
1278 const LhsEval Rv(0.0);
1279 const LhsEval Rvw(0.0);
1280 const auto& bg = gasPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rv, Rvw);
1281 return gasEnergy*bg*referenceDensity(phaseIdx, regionIdx);
1282 }
1283
1284 case waterPhaseIdx:
1285 const auto waterEnergy =
1286 waterPvt_.internalEnergy(regionIdx, T, p,
1287 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1288 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1289 assert(waterPvt_.mixingEnergy());
1291 const auto& gasEnergy =
1292 gasPvt_.internalEnergy(regionIdx, T, p,
1293 BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1294 BlackOil::template getRvw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1295 // gas miscible in water
1296 const LhsEval& Rsw = saturatedDissolutionFactor<FluidState, LhsEval>(fluidState, waterPhaseIdx, regionIdx);
1297 const LhsEval& bw = waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1298 return
1299 waterEnergy*bw*referenceDensity(waterPhaseIdx, regionIdx)
1300 + gasEnergy*Rsw*bw*referenceDensity(gasPhaseIdx, regionIdx);
1301 }
1302 const LhsEval Rsw(0.0);
1303 return
1304 waterEnergy*referenceDensity(waterPhaseIdx, regionIdx)
1305 * waterPvt_.inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration);
1306 }
1307 throw std::logic_error("Unhandled phase index " + std::to_string(phaseIdx));
1308 }
1309
1310
1311
1313 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1314 STATIC_OR_DEVICE LhsEval enthalpy(const FluidState& fluidState,
1315 unsigned phaseIdx,
1316 unsigned regionIdx) NOTHING_OR_CONST
1317 {
1318 // should preferably not be used values should be taken from intensive quantities fluid state.
1319 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1320 auto energy = internalEnergy<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1321 if(!enthalpy_eq_energy_){
1322 // used for simplified models
1323 energy += p/density<FluidState, LhsEval>(fluidState, phaseIdx, regionIdx);
1324 }
1325 return energy;
1326 }
1327
1334 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1335 STATIC_OR_DEVICE LhsEval saturatedVaporizationFactor(const FluidState& fluidState,
1336 unsigned phaseIdx,
1337 unsigned regionIdx) NOTHING_OR_CONST
1338 {
1339 assert(phaseIdx <= numPhases);
1340 assert(regionIdx <= numRegions());
1341
1342 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1343 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1344 const auto& saltConcentration = decay<LhsEval>(fluidState.saltConcentration());
1345
1346 switch (phaseIdx) {
1347 case oilPhaseIdx: return 0.0;
1348 case gasPhaseIdx: return gasPvt_.saturatedWaterVaporizationFactor(regionIdx, T, p, saltConcentration);
1349 case waterPhaseIdx: return 0.0;
1350 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1351 }
1352 }
1353
1360 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1361 STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState& fluidState,
1362 unsigned phaseIdx,
1363 unsigned regionIdx,
1364 const LhsEval& maxOilSaturation) NOTHING_OR_CONST
1365 {
1366 OPM_TIMEBLOCK_LOCAL(saturatedDissolutionFactor, Subsystem::PvtProps);
1367 assert(phaseIdx <= numPhases);
1368 assert(regionIdx <= numRegions());
1369
1370 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1371 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1372 const auto& So = (phaseIdx == waterPhaseIdx) ? 0 : decay<LhsEval>(fluidState.saturation(oilPhaseIdx));
1373
1374 switch (phaseIdx) {
1375 case oilPhaseIdx: return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p, So, maxOilSaturation);
1376 case gasPhaseIdx: return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p, So, maxOilSaturation);
1377 case waterPhaseIdx: return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1378 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1379 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1380 }
1381 }
1382
1391 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1392 STATIC_OR_DEVICE LhsEval saturatedDissolutionFactor(const FluidState& fluidState,
1393 unsigned phaseIdx,
1394 unsigned regionIdx) NOTHING_OR_CONST
1395 {
1396 OPM_TIMEBLOCK_LOCAL(saturatedDissolutionFactor, Subsystem::PvtProps);
1397 assert(phaseIdx <= numPhases);
1398 assert(regionIdx <= numRegions());
1399
1400 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1401 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1402
1403 switch (phaseIdx) {
1404 case oilPhaseIdx: return oilPvt_.saturatedGasDissolutionFactor(regionIdx, T, p);
1405 case gasPhaseIdx: return gasPvt_.saturatedOilVaporizationFactor(regionIdx, T, p);
1406 case waterPhaseIdx: return waterPvt_.saturatedGasDissolutionFactor(regionIdx, T, p,
1407 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1408 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1409 }
1410 }
1411
1415 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1416 STATIC_OR_DEVICE LhsEval bubblePointPressure(const FluidState& fluidState,
1417 unsigned regionIdx) NOTHING_OR_CONST
1418 {
1419 return saturationPressure(fluidState, oilPhaseIdx, regionIdx);
1420 }
1421
1422
1426 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1427 STATIC_OR_DEVICE LhsEval dewPointPressure(const FluidState& fluidState,
1428 unsigned regionIdx) NOTHING_OR_CONST
1429 {
1430 return saturationPressure(fluidState, gasPhaseIdx, regionIdx);
1431 }
1432
1443 template <class FluidState, class LhsEval = typename FluidState::Scalar>
1444 STATIC_OR_DEVICE LhsEval saturationPressure(const FluidState& fluidState,
1445 unsigned phaseIdx,
1446 unsigned regionIdx) NOTHING_OR_CONST
1447 {
1448 assert(phaseIdx <= numPhases);
1449 assert(regionIdx <= numRegions());
1450
1451 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1452
1453 switch (phaseIdx) {
1454 case oilPhaseIdx: return oilPvt_.saturationPressure(regionIdx, T, BlackOil::template getRs_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1455 case gasPhaseIdx: return gasPvt_.saturationPressure(regionIdx, T, BlackOil::template getRv_<ThisType, FluidState, LhsEval>(fluidState, regionIdx));
1456 case waterPhaseIdx: return waterPvt_.saturationPressure(regionIdx, T,
1457 BlackOil::template getRsw_<ThisType, FluidState, LhsEval>(fluidState, regionIdx),
1458 BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx));
1459 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1460 }
1461 }
1462
1463 /****************************************
1464 * Auxiliary and convenience methods for the black-oil model
1465 ****************************************/
1470 template <class LhsEval>
1471 STATIC_OR_DEVICE LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) NOTHING_OR_CONST
1472 {
1473 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1474 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1475
1476 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
1477 }
1478
1483 template <class LhsEval>
1484 STATIC_OR_DEVICE LhsEval convertXwGToRsw(const LhsEval& XwG, unsigned regionIdx) NOTHING_OR_CONST
1485 {
1486 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1487 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1488
1489 return XwG/(1.0 - XwG)*(rho_wRef/rho_gRef);
1490 }
1491
1496 template <class LhsEval>
1497 STATIC_OR_DEVICE LhsEval convertXgOToRv(const LhsEval& XgO, unsigned regionIdx) NOTHING_OR_CONST
1498 {
1499 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1500 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1501
1502 return XgO/(1.0 - XgO)*(rho_gRef/rho_oRef);
1503 }
1504
1509 template <class LhsEval>
1510 STATIC_OR_DEVICE LhsEval convertXgWToRvw(const LhsEval& XgW, unsigned regionIdx) NOTHING_OR_CONST
1511 {
1512 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1513 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1514
1515 return XgW/(1.0 - XgW)*(rho_gRef/rho_wRef);
1516 }
1517
1518
1523 template <class LhsEval>
1524 STATIC_OR_DEVICE LhsEval convertRsToXoG(const LhsEval& Rs, unsigned regionIdx) NOTHING_OR_CONST
1525 {
1526 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1527 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1528
1529 const LhsEval& rho_oG = Rs*rho_gRef;
1530 return rho_oG/(rho_oRef + rho_oG);
1531 }
1532
1537 template <class LhsEval>
1538 STATIC_OR_DEVICE LhsEval convertRswToXwG(const LhsEval& Rsw, unsigned regionIdx) NOTHING_OR_CONST
1539 {
1540 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1541 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1542
1543 const LhsEval& rho_wG = Rsw*rho_gRef;
1544 return rho_wG/(rho_wRef + rho_wG);
1545 }
1546
1551 template <class LhsEval>
1552 STATIC_OR_DEVICE LhsEval convertRvToXgO(const LhsEval& Rv, unsigned regionIdx) NOTHING_OR_CONST
1553 {
1554 Scalar rho_oRef = referenceDensity_[regionIdx][oilPhaseIdx];
1555 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1556
1557 const LhsEval& rho_gO = Rv*rho_oRef;
1558 return rho_gO/(rho_gRef + rho_gO);
1559 }
1560
1565 template <class LhsEval>
1566 STATIC_OR_DEVICE LhsEval convertRvwToXgW(const LhsEval& Rvw, unsigned regionIdx) NOTHING_OR_CONST
1567 {
1568 Scalar rho_wRef = referenceDensity_[regionIdx][waterPhaseIdx];
1569 Scalar rho_gRef = referenceDensity_[regionIdx][gasPhaseIdx];
1570
1571 const LhsEval& rho_gW = Rvw*rho_wRef;
1572 return rho_gW/(rho_gRef + rho_gW);
1573 }
1574
1578 template <class LhsEval>
1579 STATIC_OR_DEVICE LhsEval convertXgWToxgW(const LhsEval& XgW, unsigned regionIdx) NOTHING_OR_CONST
1580 {
1581 Scalar MW = molarMass_[regionIdx][waterCompIdx];
1582 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1583
1584 return XgW*MG / (MW*(1 - XgW) + XgW*MG);
1585 }
1586
1590 template <class LhsEval>
1591 STATIC_OR_DEVICE LhsEval convertXwGToxwG(const LhsEval& XwG, unsigned regionIdx) NOTHING_OR_CONST
1592 {
1593 Scalar MW = molarMass_[regionIdx][waterCompIdx];
1594 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1595
1596 return XwG*MW / (MG*(1 - XwG) + XwG*MW);
1597 }
1598
1602 template <class LhsEval>
1603 STATIC_OR_DEVICE LhsEval convertXoGToxoG(const LhsEval& XoG, unsigned regionIdx) NOTHING_OR_CONST
1604 {
1605 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1606 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1607
1608 return XoG*MO / (MG*(1 - XoG) + XoG*MO);
1609 }
1610
1614 template <class LhsEval>
1615 STATIC_OR_DEVICE LhsEval convertxoGToXoG(const LhsEval& xoG, unsigned regionIdx) NOTHING_OR_CONST
1616 {
1617 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1618 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1619
1620 return xoG*MG / (xoG*(MG - MO) + MO);
1621 }
1622
1626 template <class LhsEval>
1627 STATIC_OR_DEVICE LhsEval convertXgOToxgO(const LhsEval& XgO, unsigned regionIdx) NOTHING_OR_CONST
1628 {
1629 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1630 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1631
1632 return XgO*MG / (MO*(1 - XgO) + XgO*MG);
1633 }
1634
1638 template <class LhsEval>
1639 STATIC_OR_DEVICE LhsEval convertxgOToXgO(const LhsEval& xgO, unsigned regionIdx) NOTHING_OR_CONST
1640 {
1641 Scalar MO = molarMass_[regionIdx][oilCompIdx];
1642 Scalar MG = molarMass_[regionIdx][gasCompIdx];
1643
1644 return xgO*MO / (xgO*(MO - MG) + MG);
1645 }
1646
1654 STATIC_OR_DEVICE const GasPvt& gasPvt() NOTHING_OR_CONST
1655 { return gasPvt_; }
1656
1664 STATIC_OR_DEVICE const OilPvt& oilPvt() NOTHING_OR_CONST
1665 { return oilPvt_; }
1666
1674 STATIC_OR_DEVICE const WaterPvt& waterPvt() NOTHING_OR_CONST
1675 { return waterPvt_; }
1676
1682 STATIC_OR_DEVICE Scalar reservoirTemperature(unsigned = 0) NOTHING_OR_CONST
1683 { return reservoirTemperature_; }
1684
1690 STATIC_OR_DEVICE void setReservoirTemperature(Scalar value) NOTHING_OR_CONST
1691 { reservoirTemperature_ = value; }
1692
1693 STATIC_OR_DEVICE short activeToCanonicalPhaseIdx(unsigned activePhaseIdx) NOTHING_OR_CONST;
1694
1695 STATIC_OR_DEVICE short canonicalToActivePhaseIdx(unsigned phaseIdx) NOTHING_OR_CONST;
1696
1697 STATIC_OR_DEVICE short activeToCanonicalCompIdx(unsigned activeCompIdx) NOTHING_OR_CONST;
1698
1699 STATIC_OR_DEVICE short canonicalToActiveCompIdx(unsigned compIdx) NOTHING_OR_CONST;
1700
1701 STATIC_OR_DEVICE short activePhaseToActiveCompIdx(unsigned activePhaseIdx) NOTHING_OR_CONST;
1702
1703 STATIC_OR_DEVICE short activeCompToActivePhaseIdx(unsigned activeCompIdx) NOTHING_OR_CONST;
1704
1706 STATIC_OR_DEVICE Scalar diffusionCoefficient(unsigned compIdx, unsigned phaseIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
1707 { return diffusionCoefficients_[regionIdx][numPhases*compIdx + phaseIdx]; }
1708
1710 STATIC_OR_DEVICE void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx = 0) NOTHING_OR_CONST
1711 { diffusionCoefficients_[regionIdx][numPhases*compIdx + phaseIdx] = coefficient ; }
1712
1716 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
1717 STATIC_OR_DEVICE LhsEval diffusionCoefficient(const FluidState& fluidState,
1718 const ParameterCache<ParamCacheEval>& paramCache,
1719 unsigned phaseIdx,
1720 unsigned compIdx) NOTHING_OR_CONST
1721 {
1722 // diffusion is disabled by the user
1723 if(!enableDiffusion())
1724 return 0.0;
1725
1726 // diffusion coefficients are set, and we use them
1727 if(!diffusionCoefficients_.empty()) {
1728 return diffusionCoefficient(compIdx, phaseIdx, paramCache.regionIndex());
1729 }
1730
1731 const auto& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
1732 const auto& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
1733
1734 switch (phaseIdx) {
1735 case oilPhaseIdx: return oilPvt().diffusionCoefficient(T, p, compIdx);
1736 case gasPhaseIdx: return gasPvt().diffusionCoefficient(T, p, compIdx);
1737 case waterPhaseIdx: return waterPvt().diffusionCoefficient(T, p, compIdx);
1738 default: throw std::logic_error("Unhandled phase index "+std::to_string(phaseIdx));
1739 }
1740 }
1741 STATIC_OR_DEVICE void setEnergyEqualEnthalpy(bool enthalpy_eq_energy){
1742 enthalpy_eq_energy_ = enthalpy_eq_energy;
1743 }
1744
1745 STATIC_OR_DEVICE bool enthalpyEqualEnergy(){
1746 return enthalpy_eq_energy_;
1747 }
1748
1749private:
1750 STATIC_OR_NOTHING void resizeArrays_(std::size_t numRegions);
1751
1752 STATIC_OR_NOTHING Scalar reservoirTemperature_;
1753
1754 STATIC_OR_NOTHING GasPvt gasPvt_;
1755 STATIC_OR_NOTHING OilPvt oilPvt_;
1756 STATIC_OR_NOTHING WaterPvt waterPvt_;
1757
1758 STATIC_OR_NOTHING bool enableDissolvedGas_;
1759 STATIC_OR_NOTHING bool enableDissolvedGasInWater_;
1760 STATIC_OR_NOTHING bool enableVaporizedOil_;
1761 STATIC_OR_NOTHING bool enableVaporizedWater_;
1762 STATIC_OR_NOTHING bool enableDiffusion_;
1763
1764 // HACK for GCC 4.4: the array size has to be specified using the literal value '3'
1765 // here, because GCC 4.4 seems to be unable to determine the number of phases from
1766 // the BlackOil fluid system in the attribute declaration below...
1767 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numPhases=*/3> > referenceDensity_;
1768 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numComponents=*/3> > molarMass_;
1769 STATIC_OR_NOTHING Storage<std::array<Scalar, /*numComponents=*/3 * /*numPhases=*/3> > diffusionCoefficients_;
1770
1771 STATIC_OR_NOTHING PhaseUsageInfo<IndexTraits> phaseUsageInfo_;
1772
1773 STATIC_OR_NOTHING bool isInitialized_;
1774 STATIC_OR_NOTHING bool useSaturatedTables_;
1775 STATIC_OR_NOTHING bool enthalpy_eq_energy_;
1776
1777 #ifndef COMPILING_STATIC_FLUID_SYSTEM
1778 template<template<typename> typename StorageT>
1779 explicit FLUIDSYSTEM_CLASSNAME(const FLUIDSYSTEM_CLASSNAME_STATIC<Scalar, IndexTraits, StorageT>& other)
1782 , reservoirTemperature_(other.reservoirTemperature_)
1783 , gasPvt_(other.gasPvt_)
1784 , oilPvt_(other.oilPvt_)
1785 , waterPvt_(other.waterPvt_)
1786 , enableDissolvedGas_(other.enableDissolvedGas_)
1787 , enableDissolvedGasInWater_(other.enableDissolvedGasInWater_)
1788 , enableVaporizedOil_(other.enableVaporizedOil_)
1789 , enableVaporizedWater_(other.enableVaporizedWater_)
1790 , enableDiffusion_(other.enableDiffusion_)
1791 , referenceDensity_(StorageT<typename decltype(referenceDensity_)::value_type>(other.referenceDensity_))
1792 , molarMass_(StorageT<typename decltype(molarMass_)::value_type>(other.molarMass_))
1793 , diffusionCoefficients_(StorageT<typename decltype(diffusionCoefficients_)::value_type>(other.diffusionCoefficients_))
1794 , phaseUsageInfo_(other.phaseUsageInfo_)
1795 , isInitialized_(other.isInitialized_)
1796 , useSaturatedTables_(other.useSaturatedTables_)
1797 , enthalpy_eq_energy_(other.enthalpy_eq_energy_)
1798 {
1799 OPM_ERROR_IF(!other.isInitialized(), "The fluid system must be initialized before it can be copied.");
1800 }
1801
1802 template<class ScalarT, class IndexTraitsT, template<typename> typename StorageT>
1803 friend class FLUIDSYSTEM_CLASSNAME_STATIC;
1804 #else
1805 template<class ScalarT, class IndexTraitsT, template<typename> typename StorageT>
1806 friend class FLUIDSYSTEM_CLASSNAME_NONSTATIC;
1807 #endif
1808};
1809
1810template <class Scalar, class IndexTraits, template<typename> typename Storage>
1812initBegin(std::size_t numPvtRegions)
1813{
1814 isInitialized_ = false;
1815 useSaturatedTables_ = true;
1816
1817 enableDissolvedGas_ = true;
1818 enableDissolvedGasInWater_ = false;
1819 enableVaporizedOil_ = false;
1820 enableVaporizedWater_ = false;
1821 enableDiffusion_ = false;
1822
1823 surfaceTemperature = 273.15 + 15.56; // [K]
1824 surfacePressure = 1.01325e5; // [Pa]
1826
1827 phaseUsageInfo_.initFromPhases(Phases{true, true, true});
1828
1829 resizeArrays_(numPvtRegions);
1830}
1831
1832template <class Scalar, class IndexTraits, template<typename> typename Storage>
1835 Scalar rhoWater,
1836 Scalar rhoGas,
1837 unsigned regionIdx)
1838{
1839 referenceDensity_[regionIdx][oilPhaseIdx] = rhoOil;
1840 referenceDensity_[regionIdx][waterPhaseIdx] = rhoWater;
1841 referenceDensity_[regionIdx][gasPhaseIdx] = rhoGas;
1842}
1843
1844template <class Scalar, class IndexTraits, template<typename> typename Storage>
1846{
1847 // calculate the final 2D functions which are used for interpolation.
1848 const std::size_t num_regions = molarMass_.size();
1849 for (unsigned regionIdx = 0; regionIdx < num_regions; ++regionIdx) {
1850 // calculate molar masses
1851
1852 // water is simple: 18 g/mol
1853 molarMass_[regionIdx][waterCompIdx] = 18e-3;
1854
1856 // for gas, we take the density at standard conditions and assume it to be ideal
1859 Scalar rho_g = referenceDensity_[/*regionIdx=*/0][gasPhaseIdx];
1860 molarMass_[regionIdx][gasCompIdx] = Constants<Scalar>::R*T*rho_g / p;
1861 }
1862 else
1863 // hydrogen gas. we just set this do avoid NaNs later
1864 molarMass_[regionIdx][gasCompIdx] = 2e-3;
1865
1866 // finally, for oil phase, we take the molar mass from the spe9 paper
1867 molarMass_[regionIdx][oilCompIdx] = 175e-3; // kg/mol
1868 }
1869
1870 isInitialized_ = true;
1871}
1872
1873template <class Scalar, class IndexTraits, template<typename> typename Storage>
1875phaseName(unsigned phaseIdx) NOTHING_OR_CONST
1876{
1877 switch (phaseIdx) {
1878 case waterPhaseIdx:
1879 return "water";
1880 case oilPhaseIdx:
1881 return "oil";
1882 case gasPhaseIdx:
1883 return "gas";
1884
1885 default:
1886 throw std::logic_error(std::string("Phase index ") + std::to_string(phaseIdx) + " is unknown");
1887 }
1888}
1889
1890template <class Scalar, class IndexTraits, template<typename> typename Storage>
1892solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
1893{
1894 switch (phaseIdx) {
1895 case waterPhaseIdx:
1896 return waterCompIdx;
1897 case oilPhaseIdx:
1898 return oilCompIdx;
1899 case gasPhaseIdx:
1900 return gasCompIdx;
1901
1902 default:
1903 throw std::logic_error(std::string("Phase index ") + std::to_string(phaseIdx) + " is unknown");
1904 }
1905}
1906
1907template <class Scalar, class IndexTraits, template<typename> typename Storage>
1909soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
1910{
1911 switch (phaseIdx) {
1912 case waterPhaseIdx:
1914 return gasCompIdx;
1915 throw std::logic_error("The water phase does not have any solutes in the black oil model!");
1916 case oilPhaseIdx:
1917 return gasCompIdx;
1918 case gasPhaseIdx:
1919 if (enableVaporizedWater()) {
1920 return waterCompIdx;
1921 }
1922 return oilCompIdx;
1923
1924 default:
1925 throw std::logic_error(std::string("Phase index ") + std::to_string(phaseIdx) + " is unknown");
1926 }
1927}
1928
1929template <class Scalar, class IndexTraits, template<typename> typename Storage>
1931componentName(unsigned compIdx) NOTHING_OR_CONST
1932{
1933 switch (compIdx) {
1934 case waterCompIdx:
1935 return "Water";
1936 case oilCompIdx:
1937 return "Oil";
1938 case gasCompIdx:
1939 return "Gas";
1940
1941 default:
1942 throw std::logic_error(std::string("Component index ") + std::to_string(compIdx) + " is unknown");
1943 }
1944}
1945
1946template <class Scalar, class IndexTraits, template<typename> typename Storage>
1948activeToCanonicalPhaseIdx(unsigned activePhaseIdx) NOTHING_OR_CONST
1949{
1950 return phaseUsageInfo_.activeToCanonicalPhaseIdx(activePhaseIdx);
1951}
1952
1953template <class Scalar, class IndexTraits, template<typename> typename Storage>
1954NOTHING_OR_DEVICE short FLUIDSYSTEM_CLASSNAME<Scalar,IndexTraits, Storage>::
1955canonicalToActivePhaseIdx(unsigned phaseIdx) NOTHING_OR_CONST
1956{
1957 return phaseUsageInfo_.canonicalToActivePhaseIdx(phaseIdx);
1958}
1959
1960template <class Scalar, class IndexTraits, template<typename> typename Storage>
1962activeToCanonicalCompIdx(unsigned activeCompIdx) NOTHING_OR_CONST
1963{
1964 return phaseUsageInfo_.activeToCanonicalCompIdx(activeCompIdx);
1965}
1966
1967template <class Scalar, class IndexTraits, template<typename> typename Storage>
1969canonicalToActiveCompIdx(unsigned compIdx) NOTHING_OR_CONST
1970{
1971 return phaseUsageInfo_.canonicalToActiveCompIdx(compIdx);
1972}
1973
1974template <class Scalar, class IndexTraits, template<typename> typename Storage>
1976activePhaseToActiveCompIdx(unsigned activePhaseIdx) NOTHING_OR_CONST
1977{
1978 return phaseUsageInfo_.activePhaseToActiveCompIdx(activePhaseIdx);
1979}
1980
1981template <class Scalar, class IndexTraits, template<typename> typename Storage>
1983activeCompToActivePhaseIdx(unsigned activeCompIdx) NOTHING_OR_CONST
1984{
1985 return phaseUsageInfo_.activeCompToActivePhaseIdx(activeCompIdx);
1986}
1987
1988template <class Scalar, class IndexTraits, template<typename> typename Storage>
1990resizeArrays_(std::size_t numRegions)
1991{
1992 molarMass_.resize(numRegions);
1993 referenceDensity_.resize(numRegions);
1994}
1995
1996#ifdef COMPILING_STATIC_FLUID_SYSTEM
1998
1999#define DECLARE_INSTANCE(T) \
2000template<> PhaseUsageInfo<BlackOilDefaultFluidSystemIndices> BOFS<T>::phaseUsageInfo_;\
2001template<> T BOFS<T>::surfaceTemperature; \
2002template<> T BOFS<T>::surfacePressure; \
2003template<> T BOFS<T>::reservoirTemperature_; \
2004template<> bool BOFS<T>::enableDissolvedGas_; \
2005template<> bool BOFS<T>::enableDissolvedGasInWater_; \
2006template<> bool BOFS<T>::enableVaporizedOil_; \
2007template<> bool BOFS<T>::enableVaporizedWater_; \
2008template<> bool BOFS<T>::enableDiffusion_; \
2009template<> BOFS<T>::OilPvt BOFS<T>::oilPvt_; \
2010template<> BOFS<T>::GasPvt BOFS<T>::gasPvt_; \
2011template<> BOFS<T>::WaterPvt BOFS<T>::waterPvt_; \
2012template<> std::vector<std::array<T, 3>> BOFS<T>::referenceDensity_; \
2013template<> std::vector<std::array<T, 3>> BOFS<T>::molarMass_; \
2014template<> std::vector<std::array<T, 9>> BOFS<T>::diffusionCoefficients_; \
2015template<> bool BOFS<T>::isInitialized_; \
2016template<> bool BOFS<T>::useSaturatedTables_; \
2017template<> bool BOFS<T>::enthalpy_eq_energy_;
2018
2019DECLARE_INSTANCE(float)
2020DECLARE_INSTANCE(double)
2021
2022#undef DECLARE_INSTANCE
2023#endif
2024
2025#ifndef COMPILING_STATIC_FLUID_SYSTEM
2026#if HAVE_CUDA
2027namespace gpuistl
2028{
2029
2030template <class Scalar, class IndexTraits>
2031FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>
2032copy_to_gpu(const FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits>& oldFluidSystem) {
2033
2034 using GpuBuffer3Array = GpuBuffer<std::array<Scalar, 3>>;
2035 using GpuBuffer9Array = GpuBuffer<std::array<Scalar, 9>>;
2036
2037 OPM_ERROR_IF(oldFluidSystem.gasPvt_.approach() != GasPvtApproach::Co2Gas,
2038 fmt::format("Incompatible gas PVT approach. Given {}, expected {} (GasPvtApproach::Co2Gas).",
2039 int(oldFluidSystem.gasPvt_.approach()), int(GasPvtApproach::Co2Gas)));
2040
2041 OPM_ERROR_IF(oldFluidSystem.waterPvt_.approach() != WaterPvtApproach::BrineCo2,
2042 fmt::format("Incompatible water PVT approach. Given {}, expected {} (WaterPvtApproach::BrineCo2).",
2043 int(oldFluidSystem.waterPvt_.approach()), int(WaterPvtApproach::BrineCo2)));
2044
2045 OPM_ERROR_IF(oldFluidSystem.oilPvt_.isActive(),
2046 "We currently do not support an active oil phase for the FluidSystem on GPU.");
2047
2048 auto newGasPvt = copy_to_gpu(oldFluidSystem.gasPvt_.template getRealPvt<GasPvtApproach::Co2Gas>());
2049
2050 auto newOilPvt = NullOilPvt<Scalar>();
2051
2052 auto newWaterPvt = copy_to_gpu(oldFluidSystem.waterPvt_.template getRealPvt<WaterPvtApproach::BrineCo2>());
2053
2054 auto newReferenceDensity = GpuBuffer3Array(oldFluidSystem.referenceDensity_);
2055 auto newMolarMass = GpuBuffer3Array(oldFluidSystem.molarMass_);
2056 auto newDiffusionCoefficients = GpuBuffer9Array(oldFluidSystem.diffusionCoefficients_);
2057
2058 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>(
2059 oldFluidSystem.surfacePressure,
2060 oldFluidSystem.surfaceTemperature,
2061 oldFluidSystem.reservoirTemperature_,
2062 newGasPvt,
2063 newOilPvt,
2064 newWaterPvt,
2065 oldFluidSystem.enableDissolvedGas_,
2066 oldFluidSystem.enableDissolvedGasInWater_,
2067 oldFluidSystem.enableVaporizedOil_,
2068 oldFluidSystem.enableVaporizedWater_,
2069 oldFluidSystem.enableDiffusion_,
2070 std::move(newReferenceDensity),
2071 std::move(newMolarMass),
2072 std::move(newDiffusionCoefficients),
2073 oldFluidSystem.phaseUsageInfo_,
2074 oldFluidSystem.isInitialized_,
2075 oldFluidSystem.useSaturatedTables_,
2076 oldFluidSystem.enthalpy_eq_energy_
2077 );
2078}
2079
2080template <class Scalar, class IndexTraits>
2081FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>
2082make_view(FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuBuffer>& oldFluidSystem)
2083{
2084 using Array3 = std::array<Scalar, 3>;
2085 using Array9 = std::array<Scalar, 9>;
2086
2087 auto newGasPvt = make_view(oldFluidSystem.gasPvt_);
2088 auto newOilPvt = NullOilPvt<Scalar>();
2089 auto newWaterPvt = make_view(oldFluidSystem.waterPvt_);
2090
2091 auto newReferenceDensity = make_view<Array3>(oldFluidSystem.referenceDensity_);
2092 auto newMolarMass = make_view<Array3>(oldFluidSystem.molarMass_);
2093 auto newDiffusionCoefficients = make_view<Array9>(oldFluidSystem.diffusionCoefficients_);
2094
2095 return FLUIDSYSTEM_CLASSNAME<Scalar, IndexTraits, GpuView>(
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.enableDiffusion_,
2107 std::move(newReferenceDensity),
2108 std::move(newMolarMass),
2109 std::move(newDiffusionCoefficients),
2110 oldFluidSystem.phaseUsageInfo_,
2111 oldFluidSystem.isInitialized_,
2112 oldFluidSystem.useSaturatedTables_,
2113 oldFluidSystem.enthalpy_eq_energy_
2114 );
2115}
2116}
2117#endif // HAVE_CUDA
2118#endif // COMPILING_STATIC_FLUID_SYSTEM
2119} // 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:45
Definition EclipseState.hpp:62
A fluid system which uses the black-oil model assumptions to calculate termodynamically meaningful qu...
Definition BlackOilFluidSystem_macrotemplate.hpp:78
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:1706
STATIC_OR_DEVICE unsigned numActivePhases() NOTHING_OR_CONST
Returns the number of active fluid phases (i.e., usually three).
Definition BlackOilFluidSystem_macrotemplate.hpp:428
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:1566
STATIC_OR_DEVICE bool useSaturatedTables() NOTHING_OR_CONST
Returns whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:528
STATIC_OR_NOTHING void initBegin(std::size_t numPvtRegions)
Begin the initialization of the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1812
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:686
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:1552
STATIC_OR_DEVICE void initEnd()
Finish initializing the black oil fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:1845
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:1427
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:1497
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:571
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:1664
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:451
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:289
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:1361
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:503
STATIC_OR_DEVICE Scalar reservoirTemperature(unsigned=0) NOTHING_OR_CONST
Set the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1682
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:1335
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:1471
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:564
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:1579
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:340
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:1510
static constexpr unsigned gasPhaseIdx
Index of the gas phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:389
static constexpr int oilCompIdx
Index of the oil component.
Definition BlackOilFluidSystem_macrotemplate.hpp:415
FLUIDSYSTEM_CLASSNAME(const FLUIDSYSTEM_CLASSNAME< Scalar, IndexTraits, StorageT > &other)
Default copy constructor.
Definition BlackOilFluidSystem_macrotemplate.hpp:164
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:927
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:298
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:1603
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:588
STATIC_OR_DEVICE void setUseSaturatedTables(bool yesno)
Specify whether the saturated tables should be used.
Definition BlackOilFluidSystem_macrotemplate.hpp:322
STATIC_OR_DEVICE bool phaseIsActive(unsigned phaseIdx) NOTHING_OR_CONST
Returns whether a fluid phase is active.
Definition BlackOilFluidSystem_macrotemplate.hpp:432
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:1639
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:1444
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:1392
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:787
STATIC_OR_DEVICE bool enableDiffusion() NOTHING_OR_CONST
Returns whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:520
static constexpr unsigned numComponents
Number of chemical species in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:412
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:459
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:1627
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:1615
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:536
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:1591
STATIC_OR_DEVICE void setEnableDiffusion(bool yesno)
Specify whether the fluid system should consider diffusion.
Definition BlackOilFluidSystem_macrotemplate.hpp:314
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:334
STATIC_OR_DEVICE std::size_t numRegions() NOTHING_OR_CONST
Returns the number of PVT regions which are considered.
Definition BlackOilFluidSystem_macrotemplate.hpp:475
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:1875
STATIC_OR_DEVICE const PhaseUsageInfo< IndexTraits > & phaseUsage() NOTHING_OR_CONST
Returns a const reference to PhaseUsageInfo.
Definition BlackOilFluidSystem_macrotemplate.hpp:424
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:494
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:1416
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:544
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:1654
static constexpr unsigned numPhases
Number of fluid phases in the fluid system.
Definition BlackOilFluidSystem_macrotemplate.hpp:382
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:1049
STATIC_OR_DEVICE unsigned solventComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "primary" component of a phase (solvent)
Definition BlackOilFluidSystem_macrotemplate.hpp:1892
static constexpr unsigned waterPhaseIdx
Index of the water phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:385
static constexpr int gasCompIdx
Index of the gas component.
Definition BlackOilFluidSystem_macrotemplate.hpp:419
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:484
STATIC_OR_DEVICE void setReservoirTemperature(Scalar value) NOTHING_OR_CONST
Return the temperature of the reservoir.
Definition BlackOilFluidSystem_macrotemplate.hpp:1690
STATIC_OR_DEVICE unsigned soluteComponentIndex(unsigned phaseIdx) NOTHING_OR_CONST
returns the index of "secondary" component of a phase (solute)
Definition BlackOilFluidSystem_macrotemplate.hpp:1909
static constexpr unsigned oilPhaseIdx
Index of the oil phase.
Definition BlackOilFluidSystem_macrotemplate.hpp:387
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:1314
STATIC_OR_NOTHING std::string_view componentName(unsigned compIdx) NOTHING_OR_CONST
Return the human readable name of a component.
Definition BlackOilFluidSystem_macrotemplate.hpp:1931
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:1674
STATIC_OR_NOTHING Scalar surfaceTemperature
The temperature at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:395
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:463
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:1524
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:328
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:307
STATIC_OR_DEVICE void setReferenceDensities(Scalar rhoOil, Scalar rhoWater, Scalar rhoGas, unsigned regionIdx)
Initialize the values of the reference densities.
Definition BlackOilFluidSystem_macrotemplate.hpp:1834
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:512
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:1717
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:1538
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:447
STATIC_OR_NOTHING Scalar surfacePressure
The pressure at the surface.
Definition BlackOilFluidSystem_macrotemplate.hpp:392
STATIC_OR_DEVICE void setDiffusionCoefficient(Scalar coefficient, unsigned compIdx, unsigned phaseIdx, unsigned regionIdx=0) NOTHING_OR_CONST
Definition BlackOilFluidSystem_macrotemplate.hpp:1710
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:280
STATIC_OR_DEVICE bool isLiquid(unsigned phaseIdx) NOTHING_OR_CONST
Return whether a phase is liquid.
Definition BlackOilFluidSystem_macrotemplate.hpp:401
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:905
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:1484
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:551
static constexpr int waterCompIdx
Index of the water component.
Definition BlackOilFluidSystem_macrotemplate.hpp:417
This class represents the Pressure-Volume-Temperature relations of the gas phase in the black-oil mod...
Definition GasPvtMultiplexer.hpp:111
Null object for oil PVT calculations.
Definition NullOilPvt.hpp:33
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:105
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:90
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