opm-common
Loading...
Searching...
No Matches
BrineH2Pvt.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/*
4This file is part of the Open Porous Media project (OPM).
5
6OPM is free software: you can redistribute it and/or modify
7it under the terms of the GNU General Public License as published by
8the Free Software Foundation, either version 2 of the License, or
9(at your option) any later version.
10
11OPM is distributed in the hope that it will be useful,
12but WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19Consult the COPYING file in the top-level source directory of this
20module for the precise wording of the license and the list of
21copyright holders.
22*/
27#ifndef OPM_BRINE_H2_PVT_HPP
28#define OPM_BRINE_H2_PVT_HPP
29
31
39
40#include <cstddef>
41#include <vector>
42
43namespace Opm {
44
45#if HAVE_ECL_INPUT
46class EclipseState;
47class Schedule;
48#endif
49
53template <class Scalar>
54class BrineH2Pvt
55{
56 static const bool extrapolate = true;
57public:
58 using H2O = SimpleHuDuanH2O<Scalar>;
60 using H2 = ::Opm::H2<Scalar>;
61
62 // The binary coefficients for brine and H2 used by this fluid system
63 using BinaryCoeffBrineH2 = BinaryCoeff::Brine_H2<Scalar, H2O, H2>;
64
65 explicit BrineH2Pvt() = default;
66
67 explicit BrineH2Pvt(const std::vector<Scalar>& salinity,
68 Scalar T_ref = 288.71, //(273.15 + 15.56)
69 Scalar P_ref = 101325);
70
71#if HAVE_ECL_INPUT
76 void initFromState(const EclipseState& eclState, const Schedule&);
77#endif
78
79 void setNumRegions(std::size_t numRegions);
80
81 void setVapPars(const Scalar, const Scalar)
82 {
83 }
84
88 void setReferenceDensities(unsigned regionIdx,
89 Scalar rhoRefBrine,
90 Scalar rhoRefH2,
91 Scalar /*rhoRefWater*/);
92
96 void initEnd()
97 {
98 }
99
106 void setEnableDissolvedGas(bool yesno)
107 { enableDissolution_ = yesno; }
108
116 { enableSaltConcentration_ = yesno; }
117
121 unsigned numRegions() const
122 { return brineReferenceDensity_.size(); }
123
127 Scalar hVap(unsigned /*regionIdx*/) const
128 { return 0.0; }
129
130 template <class Evaluation>
131 Evaluation internalEnergy(unsigned regionIdx,
132 const Evaluation& temperature,
133 const Evaluation& pressure,
134 const Evaluation& Rs,
135 const Evaluation& saltConcentration) const
136 {
137 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
138 pressure, saltConcentration);
139 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
140 return liquidEnthalpyBrineH2_(temperature,
141 pressure,
142 salinity,
143 xlH2)
144 - pressure / density_(regionIdx, temperature, pressure, Rs, salinity);
145 }
146
150 template <class Evaluation>
151 Evaluation internalEnergy(unsigned regionIdx,
152 const Evaluation& temperature,
153 const Evaluation& pressure,
154 const Evaluation& Rs) const
155 {
156 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
157 return liquidEnthalpyBrineH2_(temperature,
158 pressure,
159 Evaluation(salinity_[regionIdx]),
160 xlH2)
161 - pressure / density_(regionIdx, temperature, pressure,
162 Rs, Evaluation(salinity_[regionIdx]));
163 }
164
168 template <class Evaluation>
169 Evaluation viscosity(unsigned regionIdx,
170 const Evaluation& temperature,
171 const Evaluation& pressure,
172 const Evaluation& /*Rs*/) const
173 {
174 //TODO: The viscosity does not yet depend on the composition
175 return saturatedViscosity(regionIdx, temperature, pressure);
176 }
177
181 template <class Evaluation>
182 Evaluation saturatedViscosity(unsigned regionIdx,
183 const Evaluation& temperature,
184 const Evaluation& pressure,
185 const Evaluation& saltConcentration) const
186 {
187 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
188 pressure, saltConcentration);
189 return Brine::liquidViscosity(temperature, pressure, salinity);
190 }
191
195 template <class Evaluation>
196 Evaluation viscosity(unsigned regionIdx,
197 const Evaluation& temperature,
198 const Evaluation& pressure,
199 const Evaluation& /*Rsw*/,
200 const Evaluation& saltConcentration) const
201 {
202 //TODO: The viscosity does not yet depend on the composition
203 return saturatedViscosity(regionIdx, temperature, pressure, saltConcentration);
204 }
205
209 template <class Evaluation>
210 Evaluation saturatedViscosity(unsigned regionIdx,
211 const Evaluation& temperature,
212 const Evaluation& pressure) const
213 {
214 return Brine::liquidViscosity(temperature, pressure, Evaluation(salinity_[regionIdx]));
215 }
216
220 template <class Evaluation>
221 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
222 const Evaluation& temperature,
223 const Evaluation& pressure,
224 const Evaluation& saltconcentration) const
225 {
226 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
227 pressure, saltconcentration);
228 Evaluation rsSat = rsSat_(regionIdx, temperature, pressure, salinity);
229 return (1.0 - convertRsToXoG_(rsSat,regionIdx))
230 * density_(regionIdx, temperature, pressure, rsSat, salinity)
231 / brineReferenceDensity_[regionIdx];
232 }
233
237 template <class Evaluation>
238 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
239 const Evaluation& temperature,
240 const Evaluation& pressure,
241 const Evaluation& Rs,
242 const Evaluation& saltConcentration) const
243 {
244 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
245 pressure, saltConcentration);
246 return (1.0 - convertRsToXoG_(Rs,regionIdx))
247 * density_(regionIdx, temperature, pressure, Rs, salinity)
248 / brineReferenceDensity_[regionIdx];
249 }
250
254 template <class Evaluation>
255 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
256 const Evaluation& temperature,
257 const Evaluation& pressure,
258 const Evaluation& Rs) const
259 {
260 return (1.0 - convertRsToXoG_(Rs, regionIdx))
261 * density_(regionIdx, temperature, pressure, Rs, Evaluation(salinity_[regionIdx]))
262 / brineReferenceDensity_[regionIdx];
263 }
264
268 template <class FluidState, class LhsEval = typename FluidState::Scalar>
269 std::pair<LhsEval, LhsEval>
270 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
271 {
272 // Deal with the possibility that we are in a two-phase H2STORE with OIL and GAS as phases.
273 const bool waterIsActive = fluidState.phaseIsActive(FluidState::waterPhaseIdx);
274 const int myPhaseIdx = waterIsActive ? FluidState::waterPhaseIdx : FluidState::oilPhaseIdx;
275 const LhsEval& Rsw = waterIsActive ? decay<LhsEval>(fluidState.Rsw()) : decay<LhsEval>(fluidState.Rs());
276
277 const LhsEval& T = decay<LhsEval>(fluidState.temperature(myPhaseIdx));
278 const LhsEval& p = decay<LhsEval>(fluidState.pressure(myPhaseIdx));
279 const LhsEval& saltConcentration
280 = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
281 // TODO: The viscosity does not yet depend on the composition
282 return { this->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration) ,
283 this->saturatedViscosity(regionIdx, T, p, saltConcentration) };
284 }
285
290 template <class Evaluation>
291 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
292 const Evaluation& temperature,
293 const Evaluation& pressure) const
294 {
295 Evaluation rsSat = rsSat_(regionIdx, temperature,
296 pressure, Evaluation(salinity_[regionIdx]));
297 return (1.0 - convertRsToXoG_(rsSat, regionIdx))
298 * density_(regionIdx, temperature, pressure, rsSat,
299 Evaluation(salinity_[regionIdx]))
300 / brineReferenceDensity_[regionIdx];
301 }
302
309 template <class Evaluation>
310 Evaluation saturationPressure(unsigned /*regionIdx*/,
311 const Evaluation& /*temperature*/,
312 const Evaluation& /*Rs*/) const
313 {
314 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
315 "has not been implemented yet!");
316 }
317
324 template <class Evaluation>
325 Evaluation saturationPressure(unsigned /*regionIdx*/,
326 const Evaluation& /*temperature*/,
327 const Evaluation& /*Rs*/,
328 const Evaluation& /*saltConcentration*/) const
329 {
330 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
331 "has not been implemented yet!");
332 }
333
337 template <class Evaluation>
338 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
339 const Evaluation& temperature,
340 const Evaluation& pressure,
341 const Evaluation& /*oilSaturation*/,
342 const Evaluation& /*maxOilSaturation*/) const
343 {
344 //TODO support VAPPARS
345 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
346 }
347
351 template <class Evaluation>
352 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
353 const Evaluation& temperature,
354 const Evaluation& pressure,
355 const Evaluation& saltConcentration) const
356 {
357 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
358 pressure, saltConcentration);
359 return rsSat_(regionIdx, temperature, pressure, salinity);
360 }
361
365 template <class Evaluation>
366 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
367 const Evaluation& temperature,
368 const Evaluation& pressure) const
369 {
370 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
371 }
372
373 Scalar oilReferenceDensity(unsigned regionIdx) const
374 { return brineReferenceDensity_[regionIdx]; }
375
376 Scalar waterReferenceDensity(unsigned regionIdx) const
377 { return brineReferenceDensity_[regionIdx]; }
378
379 Scalar gasReferenceDensity(unsigned regionIdx) const
380 { return h2ReferenceDensity_[regionIdx]; }
381
382 Scalar salinity(unsigned regionIdx) const
383 { return salinity_[regionIdx]; }
384
388 template <class Evaluation>
389 Evaluation diffusionCoefficient(const Evaluation& temperature,
390 const Evaluation& pressure,
391 unsigned /*compIdx*/) const
392 {
393 // Diffusion coefficient of H2 in pure water according to
394 // Ferrell & Himmelbau, AIChE Journal, 13(4), 1967 (Eq. 23)
395 // Some intermediate calculations and definitions
396
397 // molar volume at normal boiling point (20.271 K and 1 atm) in cm2/mol
398 const Scalar vm = 28.45;
399 const Scalar sigma = 2.96 * 1e-8; // Lennard-Jones 6-12 potential in cm (1 Å = 1e-8 cm)
400 const Scalar avogadro = 6.022e23; // Avogrado's number in mol^-1
401 const Scalar alpha = sigma / pow((vm / avogadro), 1 / 3); // Eq. (19)
402 const Scalar lambda = 1.729; // quantum parameter [-]
403 const Evaluation mu_pure = H2O::liquidViscosity(temperature, pressure, extrapolate) * 1e3; // [cP]
404 const Evaluation mu_brine = Brine::liquidViscosity(temperature, pressure,
405 Evaluation(salinity_[0])) * 1e3;
406
407 // Diffusion coeff in pure water in cm2/s
408 const Evaluation D_pure = ((4.8e-7 * temperature) / pow(mu_pure, alpha)) *
409 pow((1 + pow(lambda, 2)) / vm, 0.6);
410
411 // Diffusion coefficient in brine using Ratcliff and Holdcroft,
412 // J. G. Trans. Inst. Chem. Eng, 1963. OBS: Value
413 // of n is noted as the recommended single value according to
414 // Akita, Ind. Eng. Chem. Fundam., 1981.
415 const Evaluation log_D_brine = log10(D_pure) - 0.637 * log10(mu_brine / mu_pure);
416
417 return pow(Evaluation(10), log_D_brine) * 1e-4; // convert from cm2/s to m2/s
418 }
419
420private:
428 template <class LhsEval>
429 LhsEval density_(unsigned regionIdx,
430 const LhsEval& temperature,
431 const LhsEval& pressure,
432 const LhsEval& Rs,
433 const LhsEval& salinity) const
434 {
435 // convert Rs to mole fraction (via mass fraction)
436 LhsEval xlH2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx), salinity);
437
438 // calculate the density of solution
439 LhsEval result = liquidDensity_(temperature,
440 pressure,
441 xlH2,
442 salinity);
443
445 return result;
446 }
447
456 template <class LhsEval>
457 LhsEval liquidDensity_(const LhsEval& T,
458 const LhsEval& pl,
459 const LhsEval& xlH2,
460 const LhsEval& salinity) const
461 {
462 // check input variables
463 Valgrind::CheckDefined(T);
464 Valgrind::CheckDefined(pl);
465 Valgrind::CheckDefined(xlH2);
466
467 // check if pressure and temperature is valid
468 if (!extrapolate && T < 273.15) {
469 const std::string msg =
470 "Liquid density for Brine and H2 is only "
471 "defined above 273.15K (is " +
472 std::to_string(getValue(T)) + "K)";
473 throw NumericalProblem(msg);
474 }
475 if (!extrapolate && pl >= 2.5e8) {
476 const std::string msg =
477 "Liquid density for Brine and H2 is only "
478 "defined below 250MPa (is " +
479 std::to_string(getValue(pl)) + "Pa)";
480 throw NumericalProblem(msg);
481 }
482
483 // calculate individual contribution to density
484 const LhsEval& rho_pure = H2O::liquidDensity(T, pl, extrapolate);
485 const LhsEval& rho_brine = Brine::liquidDensity(T, pl, salinity, rho_pure);
486 const LhsEval& rho_lH2 = liquidDensityWaterH2_(T, pl, xlH2);
487 const LhsEval& contribH2 = rho_lH2 - rho_pure;
488
489 return rho_brine + contribH2;
490 }
491
501 template <class LhsEval>
502 LhsEval liquidDensityWaterH2_(const LhsEval& temperature,
503 const LhsEval& pl,
504 const LhsEval& xlH2) const
505 {
506 // molar masses
507 Scalar M_H2 = H2::molarMass();
508 Scalar M_H2O = H2O::molarMass();
509
510 // density of pure water
511 const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl, extrapolate);
512
513 // (apparent) molar volume of H2, Eq. (14) in Li et al. (2018)
514 const LhsEval& A1 = 51.1904 - 0.208062 * temperature +
515 3.4427e-4 * temperature * temperature;
516 const LhsEval& A2 = -0.022;
517 // pressure in [MPa] and Vphi in [m3/mol] (from [cm3/mol])
518 const LhsEval& V_phi = (A1 + A2 * (pl / 1e6)) / 1e6;
519
520 // density of solution, Eq. (19) in Garcia (2001)
521 const LhsEval xlH2O = 1.0 - xlH2;
522 const LhsEval& M_T = M_H2O * xlH2O + M_H2 * xlH2;
523 const LhsEval& rho_aq = 1 / (xlH2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
524
525 return rho_aq;
526 }
527
535 template <class LhsEval>
536 LhsEval convertRsToXoG_(const LhsEval& Rs, unsigned regionIdx) const
537 {
538 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
539 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
540
541 const LhsEval& rho_oG = Rs*rho_gRef;
542 return rho_oG/(rho_oRef + rho_oG);
543 }
544
551 template <class LhsEval>
552 LhsEval convertXoGToxoG_(const LhsEval& XoG, const LhsEval& salinity) const
553 {
554 Scalar M_H2 = H2::molarMass();
555 LhsEval M_Brine = Brine::molarMass(salinity);
556 return XoG*M_Brine / (M_H2*(1 - XoG) + XoG*M_Brine);
557 }
558
565 template <class LhsEval>
566 LhsEval convertxoGToXoG(const LhsEval& xoG, const LhsEval& salinity) const
567 {
568 Scalar M_H2 = H2::molarMass();
569 LhsEval M_Brine = Brine::molarMass(salinity);
570
571 return xoG*M_H2 / (xoG*(M_H2 - M_Brine) + M_Brine);
572 }
573
581 template <class LhsEval>
582 LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) const
583 {
584 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
585 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
586
587 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
588 }
589
597 template <class LhsEval>
598 LhsEval rsSat_(unsigned regionIdx,
599 const LhsEval& temperature,
600 const LhsEval& pressure,
601 const LhsEval& salinity) const
602 {
603 // Return Rs=0.0 if dissolution is disabled
604 if (!enableDissolution_)
605 return 0.0;
606
607 // calulate the equilibrium composition for the given temperature and pressure
608 LhsEval xlH2 = BinaryCoeffBrineH2::calculateMoleFractions(temperature, pressure,
609 salinity, extrapolate);
610
611 // normalize the phase compositions
612 xlH2 = max(0.0, min(1.0, xlH2));
613
614 return convertXoGToRs(convertxoGToXoG(xlH2, salinity), regionIdx);
615 }
616
617 template <class LhsEval>
618 static LhsEval liquidEnthalpyBrineH2_(const LhsEval& T,
619 const LhsEval& p,
620 const LhsEval& salinity,
621 const LhsEval& X_H2_w)
622 {
623 /* X_H2_w : mass fraction of H2 in brine */
624 /*NOTE: The heat of dissolution of H2 in brine is not included at the moment!*/
625
626 /*Numerical coefficents from PALLISER*/
627 static constexpr Scalar f[] = {
628 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
629 };
630
631 /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
632 static constexpr Scalar a[4][3] = {
633 { 9633.6, -4080.0, +286.49 },
634 { +166.58, +68.577, -4.6856 },
635 { -0.90963, -0.36524, +0.249667E-1 },
636 { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
637 };
638
639 LhsEval theta, h_NaCl;
640 LhsEval h_ls1, d_h;
641 LhsEval delta_h;
642 LhsEval hg, hw;
643
644 // Temperature in Celsius
645 theta = T - 273.15;
646
647 // Regularization
648 Scalar scalarTheta = scalarValue(theta);
649 Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
650
651 LhsEval S = salinity;
652 if (S > S_lSAT) {
653 S = S_lSAT;
654 }
655
656 hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
657
658 /*DAUBERT and DANNER*/
659 /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
660 +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
661
662 LhsEval m = 1E3/58.44 * S/(1-S);
663 d_h = 0;
664
665 for (int i = 0; i <=3 ; ++i) {
666 for (int j = 0; j <= 2; ++j) {
667 d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * pow(m, j);
668 }
669 }
670 /* heat of dissolution for halite according to Michaelides 1971 */
671 delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
672
673 /* Enthalpy of brine without H2 */
674 h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
675
676 /* enthalpy contribution of H2 gas (kJ/kg) */
677 hg = H2::gasEnthalpy(T, p, extrapolate) / 1E3;
678
679 /* Enthalpy of brine with dissolved H2 */
680 return (h_ls1 - X_H2_w*hw + hg*X_H2_w)*1E3; /*J/kg*/
681 }
682
683 template <class LhsEval>
684 const LhsEval salinityFromConcentration(unsigned regionIdx,
685 const LhsEval&T,
686 const LhsEval& P,
687 const LhsEval& saltConcentration) const
688 {
689 if (enableSaltConcentration_) {
690 return saltConcentration/H2O::liquidDensity(T, P, true);
691 }
692
693 return salinity(regionIdx);
694 }
695
696 std::vector<Scalar> brineReferenceDensity_{};
697 std::vector<Scalar> h2ReferenceDensity_{};
698 std::vector<Scalar> salinity_{};
699 bool enableDissolution_ = true;
700 bool enableSaltConcentration_ = false;
701}; // end class BrineH2Pvt
702
703} // end namespace Opm
704
705#endif
A class for the brine fluid properties.
Binary coefficients for brine and H2.
Provides the OPM specific exception classes.
Properties of pure molecular hydrogen .
A simple version of pure water with density from Hu et al.
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Some templates to wrap the valgrind client request macros.
OPM_HOST_DEVICE bool CheckDefined(const T &value)
Make valgrind complain if any of the memory occupied by an object is undefined.
Definition Valgrind.hpp:76
Binary coefficients for brine and H2.
Definition Brine_H2.hpp:41
static Evaluation calculateMoleFractions(const Evaluation &temperature, const Evaluation &pg, const Evaluation &salinity, bool extrapolate=false)
Returns the mol (!) fraction of H2 in the liquid phase for a given temperature, pressure,...
Definition Brine_H2.hpp:57
A class for the brine fluid properties.
Definition BrineDynamic.hpp:49
static OPM_HOST_DEVICE Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &, const Evaluation &salinity)
The dynamic liquid viscosity of the pure component.
Definition BrineDynamic.hpp:385
static OPM_HOST_DEVICE Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure, const Evaluation &salinity, bool extrapolate=false)
The density of the liquid component at a given pressure in and temperature in .
Definition BrineDynamic.hpp:283
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of oil saturated gas at given pressure.
Definition BrineH2Pvt.hpp:210
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the formation volume factor [-] of the fluid phase.
Definition BrineH2Pvt.hpp:255
void setEnableSaltConcentration(bool yesno)
Specify whether the PVT model should consider salt concentration from the fluidstate or a fixed salin...
Definition BrineH2Pvt.hpp:115
void initEnd()
Finish initializing the oil phase PVT properties.
Definition BrineH2Pvt.hpp:96
Evaluation saturationPressure(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the saturation pressure of the brine phase [Pa] depending on its mass fraction of the gas com...
Definition BrineH2Pvt.hpp:325
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs, const Evaluation &saltConcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition BrineH2Pvt.hpp:238
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition BrineH2Pvt.hpp:151
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltConcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition BrineH2Pvt.hpp:182
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &, const Evaluation &) const
Returns the gas dissolution factor [m^3/m^3] of the liquid phase.
Definition BrineH2Pvt.hpp:338
void setEnableDissolvedGas(bool yesno)
Specify whether the PVT model should consider that the H2 component can dissolve in the brine phase.
Definition BrineH2Pvt.hpp:106
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition BrineH2Pvt.hpp:169
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition BrineH2Pvt.hpp:270
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition BrineH2Pvt.hpp:121
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of brine saturated with H2 at a given pressure.
Definition BrineH2Pvt.hpp:291
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned) const
Diffusion coefficient of H2 in water.
Definition BrineH2Pvt.hpp:389
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns thegas dissoluiton factor [m^3/m^3] of the liquid phase.
Definition BrineH2Pvt.hpp:366
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefBrine, Scalar rhoRefH2, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition BrineH2Pvt.cpp:128
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition BrineH2Pvt.hpp:221
Evaluation saturationPressure(unsigned, const Evaluation &, const Evaluation &) const
Returns the saturation pressure of the brine phase [Pa] depending on its mass fraction of the gas com...
Definition BrineH2Pvt.hpp:310
Scalar hVap(unsigned) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition BrineH2Pvt.hpp:127
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &, const Evaluation &saltConcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition BrineH2Pvt.hpp:196
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltConcentration) const
Returns the gas dissoluiton factor [m^3/m^3] of the liquid phase.
Definition BrineH2Pvt.hpp:352
static Scalar molarMass()
Definition Component.hpp:93
Definition EclipseState.hpp:62
Properties of pure molecular hydrogen .
Definition H2.hpp:90
static constexpr Scalar molarMass()
The molar mass in of molecular hydrogen.
Definition H2.hpp:109
static const Evaluation gasEnthalpy(Evaluation temperature, Evaluation pressure, bool extrapolate=false)
Specific enthalpy of pure hydrogen gas.
Definition H2.hpp:267
Definition Schedule.hpp:101
A simple version of pure water with density from Hu et al.
Definition SimpleHuDuanH2O.hpp:66
static OPM_HOST_DEVICE Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate)
The density of pure water at a given pressure and temperature .
Definition SimpleHuDuanH2O.hpp:316
static OPM_HOST_DEVICE Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate)
The dynamic viscosity of pure water.
Definition SimpleHuDuanH2O.hpp:361
static OPM_HOST_DEVICE Scalar molarMass()
The molar mass in of water.
Definition SimpleHuDuanH2O.hpp:103
static OPM_HOST_DEVICE Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of liquid water .
Definition SimpleHuDuanH2O.hpp:202
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Scalar Brine< Scalar, H2O >::salinity
Default value for the salinity of the brine (dimensionless).
Definition Brine.hpp:391