My Project
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
38
39#include <cstddef>
40#include <vector>
41
42namespace Opm {
43
44#if HAVE_ECL_INPUT
45class EclipseState;
46class Schedule;
47#endif
48
52template <class Scalar>
54{
55 static const bool extrapolate = true;
56public:
59 using H2 = ::Opm::H2<Scalar>;
60
61 // The binary coefficients for brine and H2 used by this fluid system
63
64 explicit BrineH2Pvt() = default;
65
66 explicit BrineH2Pvt(const std::vector<Scalar>& salinity,
67 Scalar T_ref = 288.71, //(273.15 + 15.56)
68 Scalar P_ref = 101325);
69
70#if HAVE_ECL_INPUT
75 void initFromState(const EclipseState& eclState, const Schedule&);
76#endif
77
78 void setNumRegions(std::size_t numRegions);
79
80 void setVapPars(const Scalar, const Scalar)
81 {
82 }
83
87 void setReferenceDensities(unsigned regionIdx,
88 Scalar rhoRefBrine,
89 Scalar rhoRefH2,
90 Scalar /*rhoRefWater*/);
91
95 void initEnd()
96 {
97 }
98
105 void setEnableDissolvedGas(bool yesno)
106 { enableDissolution_ = yesno; }
107
115 { enableSaltConcentration_ = yesno; }
116
120 unsigned numRegions() const
121 { return brineReferenceDensity_.size(); }
122
126 Scalar hVap(unsigned /*regionIdx*/) const
127 { return 0.0; }
128
129 template <class Evaluation>
130 Evaluation internalEnergy(unsigned regionIdx,
131 const Evaluation& temperature,
132 const Evaluation& pressure,
133 const Evaluation& Rs,
134 const Evaluation& saltConcentration) const
135 {
136 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
137 pressure, saltConcentration);
138 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
139 return liquidEnthalpyBrineH2_(temperature,
140 pressure,
141 salinity,
142 xlH2)
143 - pressure / density_(regionIdx, temperature, pressure, Rs, salinity);
144 }
145
149 template <class Evaluation>
150 Evaluation internalEnergy(unsigned regionIdx,
151 const Evaluation& temperature,
152 const Evaluation& pressure,
153 const Evaluation& Rs) const
154 {
155 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
156 return liquidEnthalpyBrineH2_(temperature,
157 pressure,
158 Evaluation(salinity_[regionIdx]),
159 xlH2)
160 - pressure / density_(regionIdx, temperature, pressure,
161 Rs, Evaluation(salinity_[regionIdx]));
162 }
163
167 template <class Evaluation>
168 Evaluation viscosity(unsigned regionIdx,
169 const Evaluation& temperature,
170 const Evaluation& pressure,
171 const Evaluation& /*Rs*/) const
172 {
173 //TODO: The viscosity does not yet depend on the composition
174 return saturatedViscosity(regionIdx, temperature, pressure);
175 }
176
180 template <class Evaluation>
181 Evaluation saturatedViscosity(unsigned regionIdx,
182 const Evaluation& temperature,
183 const Evaluation& pressure,
184 const Evaluation& saltConcentration) const
185 {
186 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
187 pressure, saltConcentration);
188 return Brine::liquidViscosity(temperature, pressure, salinity);
189 }
190
194 template <class Evaluation>
195 Evaluation viscosity(unsigned regionIdx,
196 const Evaluation& temperature,
197 const Evaluation& pressure,
198 const Evaluation& /*Rsw*/,
199 const Evaluation& saltConcentration) const
200 {
201 //TODO: The viscosity does not yet depend on the composition
202 return saturatedViscosity(regionIdx, temperature, pressure, saltConcentration);
203 }
204
208 template <class Evaluation>
209 Evaluation saturatedViscosity(unsigned regionIdx,
210 const Evaluation& temperature,
211 const Evaluation& pressure) const
212 {
213 return Brine::liquidViscosity(temperature, pressure, Evaluation(salinity_[regionIdx]));
214 }
215
219 template <class Evaluation>
220 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
221 const Evaluation& temperature,
222 const Evaluation& pressure,
223 const Evaluation& saltconcentration) const
224 {
225 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
226 pressure, saltconcentration);
227 Evaluation rsSat = rsSat_(regionIdx, temperature, pressure, salinity);
228 return (1.0 - convertRsToXoG_(rsSat,regionIdx))
229 * density_(regionIdx, temperature, pressure, rsSat, salinity)
230 / brineReferenceDensity_[regionIdx];
231 }
232
236 template <class Evaluation>
237 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
238 const Evaluation& temperature,
239 const Evaluation& pressure,
240 const Evaluation& Rs,
241 const Evaluation& saltConcentration) const
242 {
243 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
244 pressure, saltConcentration);
245 return (1.0 - convertRsToXoG_(Rs,regionIdx))
246 * density_(regionIdx, temperature, pressure, Rs, salinity)
247 / brineReferenceDensity_[regionIdx];
248 }
249
253 template <class Evaluation>
254 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
255 const Evaluation& temperature,
256 const Evaluation& pressure,
257 const Evaluation& Rs) const
258 {
259 return (1.0 - convertRsToXoG_(Rs, regionIdx))
260 * density_(regionIdx, temperature, pressure, Rs, Evaluation(salinity_[regionIdx]))
261 / brineReferenceDensity_[regionIdx];
262 }
263
268 template <class Evaluation>
269 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
270 const Evaluation& temperature,
271 const Evaluation& pressure) const
272 {
273 Evaluation rsSat = rsSat_(regionIdx, temperature,
274 pressure, Evaluation(salinity_[regionIdx]));
275 return (1.0 - convertRsToXoG_(rsSat, regionIdx))
276 * density_(regionIdx, temperature, pressure, rsSat,
277 Evaluation(salinity_[regionIdx]))
278 / brineReferenceDensity_[regionIdx];
279 }
280
287 template <class Evaluation>
288 Evaluation saturationPressure(unsigned /*regionIdx*/,
289 const Evaluation& /*temperature*/,
290 const Evaluation& /*Rs*/) const
291 {
292 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
293 "has not been implemented yet!");
294 }
295
302 template <class Evaluation>
303 Evaluation saturationPressure(unsigned /*regionIdx*/,
304 const Evaluation& /*temperature*/,
305 const Evaluation& /*Rs*/,
306 const Evaluation& /*saltConcentration*/) const
307 {
308 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
309 "has not been implemented yet!");
310 }
311
315 template <class Evaluation>
316 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
317 const Evaluation& temperature,
318 const Evaluation& pressure,
319 const Evaluation& /*oilSaturation*/,
320 const Evaluation& /*maxOilSaturation*/) const
321 {
322 //TODO support VAPPARS
323 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
324 }
325
329 template <class Evaluation>
330 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
331 const Evaluation& temperature,
332 const Evaluation& pressure,
333 const Evaluation& saltConcentration) const
334 {
335 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
336 pressure, saltConcentration);
337 return rsSat_(regionIdx, temperature, pressure, salinity);
338 }
339
343 template <class Evaluation>
344 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
345 const Evaluation& temperature,
346 const Evaluation& pressure) const
347 {
348 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
349 }
350
351 Scalar oilReferenceDensity(unsigned regionIdx) const
352 { return brineReferenceDensity_[regionIdx]; }
353
354 Scalar waterReferenceDensity(unsigned regionIdx) const
355 { return brineReferenceDensity_[regionIdx]; }
356
357 Scalar gasReferenceDensity(unsigned regionIdx) const
358 { return h2ReferenceDensity_[regionIdx]; }
359
360 Scalar salinity(unsigned regionIdx) const
361 { return salinity_[regionIdx]; }
362
366 template <class Evaluation>
367 Evaluation diffusionCoefficient(const Evaluation& temperature,
368 const Evaluation& pressure,
369 unsigned /*compIdx*/) const
370 {
371 // Diffusion coefficient of H2 in pure water according to
372 // Ferrell & Himmelbau, AIChE Journal, 13(4), 1967 (Eq. 23)
373 // Some intermediate calculations and definitions
374
375 // molar volume at normal boiling point (20.271 K and 1 atm) in cm2/mol
376 const Scalar vm = 28.45;
377 const Scalar sigma = 2.96 * 1e-8; // Lennard-Jones 6-12 potential in cm (1 Å = 1e-8 cm)
378 const Scalar avogadro = 6.022e23; // Avogrado's number in mol^-1
379 const Scalar alpha = sigma / pow((vm / avogadro), 1 / 3); // Eq. (19)
380 const Scalar lambda = 1.729; // quantum parameter [-]
381 const Evaluation mu_pure = H2O::liquidViscosity(temperature, pressure, extrapolate) * 1e3; // [cP]
382 const Evaluation mu_brine = Brine::liquidViscosity(temperature, pressure,
383 Evaluation(salinity_[0])) * 1e3;
384
385 // Diffusion coeff in pure water in cm2/s
386 const Evaluation D_pure = ((4.8e-7 * temperature) / pow(mu_pure, alpha)) *
387 pow((1 + pow(lambda, 2)) / vm, 0.6);
388
389 // Diffusion coefficient in brine using Ratcliff and Holdcroft,
390 // J. G. Trans. Inst. Chem. Eng, 1963. OBS: Value
391 // of n is noted as the recommended single value according to
392 // Akita, Ind. Eng. Chem. Fundam., 1981.
393 const Evaluation log_D_brine = log10(D_pure) - 0.637 * log10(mu_brine / mu_pure);
394
395 return pow(Evaluation(10), log_D_brine) * 1e-4; // convert from cm2/s to m2/s
396 }
397
398private:
406 template <class LhsEval>
407 LhsEval density_(unsigned regionIdx,
408 const LhsEval& temperature,
409 const LhsEval& pressure,
410 const LhsEval& Rs,
411 const LhsEval& salinity) const
412 {
413 // convert Rs to mole fraction (via mass fraction)
414 LhsEval xlH2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx), salinity);
415
416 // calculate the density of solution
417 LhsEval result = liquidDensity_(temperature,
418 pressure,
419 xlH2,
420 salinity);
421
422 Valgrind::CheckDefined(result);
423 return result;
424 }
425
434 template <class LhsEval>
435 LhsEval liquidDensity_(const LhsEval& T,
436 const LhsEval& pl,
437 const LhsEval& xlH2,
438 const LhsEval& salinity) const
439 {
440 // check input variables
441 Valgrind::CheckDefined(T);
442 Valgrind::CheckDefined(pl);
443 Valgrind::CheckDefined(xlH2);
444
445 // check if pressure and temperature is valid
446 if (!extrapolate && T < 273.15) {
447 const std::string msg =
448 "Liquid density for Brine and H2 is only "
449 "defined above 273.15K (is " +
450 std::to_string(getValue(T)) + "K)";
451 throw NumericalProblem(msg);
452 }
453 if (!extrapolate && pl >= 2.5e8) {
454 const std::string msg =
455 "Liquid density for Brine and H2 is only "
456 "defined below 250MPa (is " +
457 std::to_string(getValue(pl)) + "Pa)";
458 throw NumericalProblem(msg);
459 }
460
461 // calculate individual contribution to density
462 const LhsEval& rho_brine = Brine::liquidDensity(T, pl, salinity, extrapolate);
463 const LhsEval& rho_pure = H2O::liquidDensity(T, pl, extrapolate);
464 const LhsEval& rho_lH2 = liquidDensityWaterH2_(T, pl, xlH2);
465 const LhsEval& contribH2 = rho_lH2 - rho_pure;
466
467 return rho_brine + contribH2;
468 }
469
479 template <class LhsEval>
480 LhsEval liquidDensityWaterH2_(const LhsEval& temperature,
481 const LhsEval& pl,
482 const LhsEval& xlH2) const
483 {
484 // molar masses
485 Scalar M_H2 = H2::molarMass();
486 Scalar M_H2O = H2O::molarMass();
487
488 // density of pure water
489 const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl, extrapolate);
490
491 // (apparent) molar volume of H2, Eq. (14) in Li et al. (2018)
492 const LhsEval& A1 = 51.1904 - 0.208062 * temperature +
493 3.4427e-4 * temperature * temperature;
494 const LhsEval& A2 = -0.022;
495 // pressure in [MPa] and Vphi in [m3/mol] (from [cm3/mol])
496 const LhsEval& V_phi = (A1 + A2 * (pl / 1e6)) / 1e6;
497
498 // density of solution, Eq. (19) in Garcia (2001)
499 const LhsEval xlH2O = 1.0 - xlH2;
500 const LhsEval& M_T = M_H2O * xlH2O + M_H2 * xlH2;
501 const LhsEval& rho_aq = 1 / (xlH2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
502
503 return rho_aq;
504 }
505
513 template <class LhsEval>
514 LhsEval convertRsToXoG_(const LhsEval& Rs, unsigned regionIdx) const
515 {
516 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
517 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
518
519 const LhsEval& rho_oG = Rs*rho_gRef;
520 return rho_oG/(rho_oRef + rho_oG);
521 }
522
529 template <class LhsEval>
530 LhsEval convertXoGToxoG_(const LhsEval& XoG, const LhsEval& salinity) const
531 {
532 Scalar M_H2 = H2::molarMass();
533 LhsEval M_Brine = Brine::molarMass(salinity);
534 return XoG*M_Brine / (M_H2*(1 - XoG) + XoG*M_Brine);
535 }
536
543 template <class LhsEval>
544 LhsEval convertxoGToXoG(const LhsEval& xoG, const LhsEval& salinity) const
545 {
546 Scalar M_H2 = H2::molarMass();
547 LhsEval M_Brine = Brine::molarMass(salinity);
548
549 return xoG*M_H2 / (xoG*(M_H2 - M_Brine) + M_Brine);
550 }
551
559 template <class LhsEval>
560 LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) const
561 {
562 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
563 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
564
565 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
566 }
567
575 template <class LhsEval>
576 LhsEval rsSat_(unsigned regionIdx,
577 const LhsEval& temperature,
578 const LhsEval& pressure,
579 const LhsEval& salinity) const
580 {
581 // Return Rs=0.0 if dissolution is disabled
582 if (!enableDissolution_)
583 return 0.0;
584
585 // calulate the equilibrium composition for the given temperature and pressure
586 LhsEval xlH2 = BinaryCoeffBrineH2::calculateMoleFractions(temperature, pressure,
587 salinity, extrapolate);
588
589 // normalize the phase compositions
590 xlH2 = max(0.0, min(1.0, xlH2));
591
592 return convertXoGToRs(convertxoGToXoG(xlH2, salinity), regionIdx);
593 }
594
595 template <class LhsEval>
596 static LhsEval liquidEnthalpyBrineH2_(const LhsEval& T,
597 const LhsEval& p,
598 const LhsEval& salinity,
599 const LhsEval& X_H2_w)
600 {
601 /* X_H2_w : mass fraction of H2 in brine */
602 /*NOTE: The heat of dissolution of H2 in brine is not included at the moment!*/
603
604 /*Numerical coefficents from PALLISER*/
605 static constexpr Scalar f[] = {
606 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
607 };
608
609 /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
610 static constexpr Scalar a[4][3] = {
611 { 9633.6, -4080.0, +286.49 },
612 { +166.58, +68.577, -4.6856 },
613 { -0.90963, -0.36524, +0.249667E-1 },
614 { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
615 };
616
617 LhsEval theta, h_NaCl;
618 LhsEval h_ls1, d_h;
619 LhsEval delta_h;
620 LhsEval hg, hw;
621
622 // Temperature in Celsius
623 theta = T - 273.15;
624
625 // Regularization
626 Scalar scalarTheta = scalarValue(theta);
627 Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
628
629 LhsEval S = salinity;
630 if (S > S_lSAT) {
631 S = S_lSAT;
632 }
633
634 hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
635
636 /*DAUBERT and DANNER*/
637 /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
638 +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
639
640 LhsEval m = 1E3/58.44 * S/(1-S);
641 d_h = 0;
642
643 for (int i = 0; i <=3 ; ++i) {
644 for (int j = 0; j <= 2; ++j) {
645 d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * pow(m, j);
646 }
647 }
648 /* heat of dissolution for halite according to Michaelides 1971 */
649 delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
650
651 /* Enthalpy of brine without H2 */
652 h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
653
654 /* enthalpy contribution of H2 gas (kJ/kg) */
655 hg = H2::gasEnthalpy(T, p, extrapolate) / 1E3;
656
657 /* Enthalpy of brine with dissolved H2 */
658 return (h_ls1 - X_H2_w*hw + hg*X_H2_w)*1E3; /*J/kg*/
659 }
660
661 template <class LhsEval>
662 const LhsEval salinityFromConcentration(unsigned regionIdx,
663 const LhsEval&T,
664 const LhsEval& P,
665 const LhsEval& saltConcentration) const
666 {
667 if (enableSaltConcentration_) {
668 return saltConcentration/H2O::liquidDensity(T, P, true);
669 }
670
671 return salinity(regionIdx);
672 }
673
674 std::vector<Scalar> brineReferenceDensity_{};
675 std::vector<Scalar> h2ReferenceDensity_{};
676 std::vector<Scalar> salinity_{};
677 bool enableDissolution_ = true;
678 bool enableSaltConcentration_ = false;
679}; // end class BrineH2Pvt
680
681} // end namespace Opm
682
683#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.
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 viscosity of pure water.
Definition BrineDynamic.hpp:341
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:264
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a H2-Brine sy...
Definition BrineH2Pvt.hpp:54
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:209
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:254
void setEnableSaltConcentration(bool yesno)
Specify whether the PVT model should consider salt concentration from the fluidstate or a fixed salin...
Definition BrineH2Pvt.hpp:114
void initEnd()
Finish initializing the oil phase PVT properties.
Definition BrineH2Pvt.hpp:95
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:303
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:237
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:150
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:181
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:316
void setEnableDissolvedGas(bool yesno)
Specify whether the PVT model should consider that the H2 component can dissolve in the brine phase.
Definition BrineH2Pvt.hpp:105
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:168
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition BrineH2Pvt.hpp:120
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:269
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned) const
Diffusion coefficient of H2 in water.
Definition BrineH2Pvt.hpp:367
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:344
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:220
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:288
Scalar hVap(unsigned) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition BrineH2Pvt.hpp:126
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:195
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:330
static Scalar molarMass()
The molar mass in of the component.
Definition Component.hpp:93
Definition EclipseState.hpp:63
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:264
Definition Schedule.hpp:101
A simple version of pure water with density from Hu et al.
Definition SimpleHuDuanH2O.hpp:65
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:313
static OPM_HOST_DEVICE Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate)
The dynamic viscosity of pure water.
Definition SimpleHuDuanH2O.hpp:358
static OPM_HOST_DEVICE Scalar molarMass()
The molar mass in of water.
Definition SimpleHuDuanH2O.hpp:102
static OPM_HOST_DEVICE Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of liquid water .
Definition SimpleHuDuanH2O.hpp:201
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30