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
45class EclipseState;
46class Schedule;
47
51template <class Scalar>
52class BrineH2Pvt
53{
54 static const bool extrapolate = true;
55public:
56 using H2O = SimpleHuDuanH2O<Scalar>;
58 using H2 = ::Opm::H2<Scalar>;
59
60 // The binary coefficients for brine and H2 used by this fluid system
61 using BinaryCoeffBrineH2 = BinaryCoeff::Brine_H2<Scalar, H2O, H2>;
62
63 explicit BrineH2Pvt() = default;
64
65 explicit BrineH2Pvt(const std::vector<Scalar>& salinity,
66 Scalar T_ref = 288.71, //(273.15 + 15.56)
67 Scalar P_ref = 101325);
68
73 void initFromState(const EclipseState& eclState, const Schedule&);
74
75 void setNumRegions(std::size_t numRegions);
76
77 void setVapPars(const Scalar, const Scalar)
78 {
79 }
80
84 void setReferenceDensities(unsigned regionIdx,
85 Scalar rhoRefBrine,
86 Scalar rhoRefH2,
87 Scalar /*rhoRefWater*/);
88
92 void initEnd()
93 {
94 }
95
102 void setEnableDissolvedGas(bool yesno)
103 { enableDissolution_ = yesno; }
104
112 { enableSaltConcentration_ = yesno; }
113
117 unsigned numRegions() const
118 { return brineReferenceDensity_.size(); }
119
123 Scalar hVap(unsigned /*regionIdx*/) const
124 { return 0.0; }
125
126 template <class Evaluation>
127 Evaluation internalEnergy(unsigned regionIdx,
128 const Evaluation& temperature,
129 const Evaluation& pressure,
130 const Evaluation& Rs,
131 const Evaluation& saltConcentration) const
132 {
133 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
134 pressure, saltConcentration);
135 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
136 return liquidEnthalpyBrineH2_(temperature,
137 pressure,
138 salinity,
139 xlH2)
140 - pressure / density_(regionIdx, temperature, pressure, Rs, salinity);
141 }
142
146 template <class Evaluation>
147 Evaluation internalEnergy(unsigned regionIdx,
148 const Evaluation& temperature,
149 const Evaluation& pressure,
150 const Evaluation& Rs) const
151 {
152 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
153 return liquidEnthalpyBrineH2_(temperature,
154 pressure,
155 Evaluation(salinity_[regionIdx]),
156 xlH2)
157 - pressure / density_(regionIdx, temperature, pressure,
158 Rs, Evaluation(salinity_[regionIdx]));
159 }
160
164 template <class Evaluation>
165 Evaluation viscosity(unsigned regionIdx,
166 const Evaluation& temperature,
167 const Evaluation& pressure,
168 const Evaluation& /*Rs*/) const
169 {
170 //TODO: The viscosity does not yet depend on the composition
171 return saturatedViscosity(regionIdx, temperature, pressure);
172 }
173
177 template <class Evaluation>
178 Evaluation saturatedViscosity(unsigned regionIdx,
179 const Evaluation& temperature,
180 const Evaluation& pressure,
181 const Evaluation& saltConcentration) const
182 {
183 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
184 pressure, saltConcentration);
185 return Brine::liquidViscosity(temperature, pressure, salinity);
186 }
187
191 template <class Evaluation>
192 Evaluation viscosity(unsigned regionIdx,
193 const Evaluation& temperature,
194 const Evaluation& pressure,
195 const Evaluation& /*Rsw*/,
196 const Evaluation& saltConcentration) const
197 {
198 //TODO: The viscosity does not yet depend on the composition
199 return saturatedViscosity(regionIdx, temperature, pressure, saltConcentration);
200 }
201
205 template <class Evaluation>
206 Evaluation saturatedViscosity(unsigned regionIdx,
207 const Evaluation& temperature,
208 const Evaluation& pressure) const
209 {
210 return Brine::liquidViscosity(temperature, pressure, Evaluation(salinity_[regionIdx]));
211 }
212
216 template <class Evaluation>
217 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
218 const Evaluation& temperature,
219 const Evaluation& pressure,
220 const Evaluation& saltconcentration) const
221 {
222 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
223 pressure, saltconcentration);
224 Evaluation rsSat = rsSat_(regionIdx, temperature, pressure, salinity);
225 return (1.0 - convertRsToXoG_(rsSat,regionIdx))
226 * density_(regionIdx, temperature, pressure, rsSat, salinity)
227 / brineReferenceDensity_[regionIdx];
228 }
229
233 template <class Evaluation>
234 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
235 const Evaluation& temperature,
236 const Evaluation& pressure,
237 const Evaluation& Rs,
238 const Evaluation& saltConcentration) const
239 {
240 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
241 pressure, saltConcentration);
242 return (1.0 - convertRsToXoG_(Rs,regionIdx))
243 * density_(regionIdx, temperature, pressure, Rs, salinity)
244 / brineReferenceDensity_[regionIdx];
245 }
246
250 template <class Evaluation>
251 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
252 const Evaluation& temperature,
253 const Evaluation& pressure,
254 const Evaluation& Rs) const
255 {
256 return (1.0 - convertRsToXoG_(Rs, regionIdx))
257 * density_(regionIdx, temperature, pressure, Rs, Evaluation(salinity_[regionIdx]))
258 / brineReferenceDensity_[regionIdx];
259 }
260
264 template <class FluidState, class LhsEval = typename FluidState::ValueType>
265 std::pair<LhsEval, LhsEval>
266 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
267 {
268 // Deal with the possibility that we are in a two-phase H2STORE with OIL and GAS as phases.
269 const bool waterIsActive = fluidState.phaseIsActive(FluidState::waterPhaseIdx);
270 const int myPhaseIdx = waterIsActive ? FluidState::waterPhaseIdx : FluidState::oilPhaseIdx;
271 const LhsEval& Rsw = waterIsActive ? decay<LhsEval>(fluidState.Rsw()) : decay<LhsEval>(fluidState.Rs());
272
273 const LhsEval& T = decay<LhsEval>(fluidState.temperature(myPhaseIdx));
274 const LhsEval& p = decay<LhsEval>(fluidState.pressure(myPhaseIdx));
275 const LhsEval& saltConcentration
276 = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
277 // TODO: The viscosity does not yet depend on the composition
278 return { this->inverseFormationVolumeFactor(regionIdx, T, p, Rsw, saltConcentration) ,
279 this->saturatedViscosity(regionIdx, T, p, saltConcentration) };
280 }
281
286 template <class Evaluation>
287 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
288 const Evaluation& temperature,
289 const Evaluation& pressure) const
290 {
291 Evaluation rsSat = rsSat_(regionIdx, temperature,
292 pressure, Evaluation(salinity_[regionIdx]));
293 return (1.0 - convertRsToXoG_(rsSat, regionIdx))
294 * density_(regionIdx, temperature, pressure, rsSat,
295 Evaluation(salinity_[regionIdx]))
296 / brineReferenceDensity_[regionIdx];
297 }
298
305 template <class Evaluation>
306 Evaluation saturationPressure(unsigned /*regionIdx*/,
307 const Evaluation& /*temperature*/,
308 const Evaluation& /*Rs*/) const
309 {
310 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
311 "has not been implemented yet!");
312 }
313
320 template <class Evaluation>
321 Evaluation saturationPressure(unsigned /*regionIdx*/,
322 const Evaluation& /*temperature*/,
323 const Evaluation& /*Rs*/,
324 const Evaluation& /*saltConcentration*/) const
325 {
326 throw std::runtime_error("Saturation pressure for the Brine-H2 PVT module "
327 "has not been implemented yet!");
328 }
329
333 template <class Evaluation>
334 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
335 const Evaluation& temperature,
336 const Evaluation& pressure,
337 const Evaluation& /*oilSaturation*/,
338 const Evaluation& /*maxOilSaturation*/) const
339 {
340 //TODO support VAPPARS
341 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
342 }
343
347 template <class Evaluation>
348 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
349 const Evaluation& temperature,
350 const Evaluation& pressure,
351 const Evaluation& saltConcentration) const
352 {
353 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
354 pressure, saltConcentration);
355 return rsSat_(regionIdx, temperature, pressure, salinity);
356 }
357
361 template <class Evaluation>
362 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
363 const Evaluation& temperature,
364 const Evaluation& pressure) const
365 {
366 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
367 }
368
369 Scalar oilReferenceDensity(unsigned regionIdx) const
370 { return brineReferenceDensity_[regionIdx]; }
371
372 Scalar waterReferenceDensity(unsigned regionIdx) const
373 { return brineReferenceDensity_[regionIdx]; }
374
375 Scalar gasReferenceDensity(unsigned regionIdx) const
376 { return h2ReferenceDensity_[regionIdx]; }
377
378 Scalar salinity(unsigned regionIdx) const
379 { return salinity_[regionIdx]; }
380
384 template <class Evaluation>
385 Evaluation diffusionCoefficient(const Evaluation& temperature,
386 const Evaluation& pressure,
387 unsigned /*compIdx*/,
388 unsigned regionIdx = 0) const
389 {
390 // Diffusion coefficient of H2 in pure water according to
391 // Ferrell & Himmelbau, AIChE Journal, 13(4), 1967 (Eq. 23)
392 // Some intermediate calculations and definitions
393
394 // molar volume at normal boiling point (20.271 K and 1 atm) in cm2/mol
395 const Scalar vm = 28.45;
396 const Scalar sigma = 2.96 * 1e-8; // Lennard-Jones 6-12 potential in cm (1 Å = 1e-8 cm)
397 const Scalar avogadro = 6.022e23; // Avogrado's number in mol^-1
398 const Scalar alpha = sigma / pow((vm / avogadro), 1 / 3); // Eq. (19)
399 const Scalar lambda = 1.729; // quantum parameter [-]
400 const Evaluation mu_pure = H2O::liquidViscosity(temperature, pressure, extrapolate) * 1e3; // [cP]
401 const Evaluation mu_brine = Brine::liquidViscosity(temperature, pressure,
402 Evaluation(salinity_[regionIdx])) * 1e3;
403
404 // Diffusion coeff in pure water in cm2/s
405 const Evaluation D_pure = ((4.8e-7 * temperature) / pow(mu_pure, alpha)) *
406 pow((1 + pow(lambda, 2)) / vm, 0.6);
407
408 // Diffusion coefficient in brine using Ratcliff and Holdcroft,
409 // J. G. Trans. Inst. Chem. Eng, 1963. OBS: Value
410 // of n is noted as the recommended single value according to
411 // Akita, Ind. Eng. Chem. Fundam., 1981.
412 const Evaluation log_D_brine = log10(D_pure) - 0.637 * log10(mu_brine / mu_pure);
413
414 return pow(Evaluation(10), log_D_brine) * 1e-4; // convert from cm2/s to m2/s
415 }
416
417private:
425 template <class LhsEval>
426 LhsEval density_(unsigned regionIdx,
427 const LhsEval& temperature,
428 const LhsEval& pressure,
429 const LhsEval& Rs,
430 const LhsEval& salinity) const
431 {
432 // convert Rs to mole fraction (via mass fraction)
433 LhsEval xlH2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx), salinity);
434
435 // calculate the density of solution
436 LhsEval result = liquidDensity_(temperature,
437 pressure,
438 xlH2,
439 salinity);
440
442 return result;
443 }
444
453 template <class LhsEval>
454 LhsEval liquidDensity_(const LhsEval& T,
455 const LhsEval& pl,
456 const LhsEval& xlH2,
457 const LhsEval& salinity) const
458 {
459 // check input variables
460 Valgrind::CheckDefined(T);
461 Valgrind::CheckDefined(pl);
462 Valgrind::CheckDefined(xlH2);
463
464 // check if pressure and temperature is valid
465 if (!extrapolate && T < 273.15) {
466 const std::string msg =
467 "Liquid density for Brine and H2 is only "
468 "defined above 273.15K (is " +
469 std::to_string(getValue(T)) + "K)";
470 throw NumericalProblem(msg);
471 }
472 if (!extrapolate && pl >= 2.5e8) {
473 const std::string msg =
474 "Liquid density for Brine and H2 is only "
475 "defined below 250MPa (is " +
476 std::to_string(getValue(pl)) + "Pa)";
477 throw NumericalProblem(msg);
478 }
479
480 // calculate individual contribution to density
481 const LhsEval& rho_pure = H2O::liquidDensity(T, pl, extrapolate);
482 const LhsEval& rho_brine = Brine::liquidDensity(T, pl, salinity, rho_pure);
483 const LhsEval& rho_lH2 = liquidDensityWaterH2_(T, pl, xlH2);
484 const LhsEval& contribH2 = rho_lH2 - rho_pure;
485
486 return rho_brine + contribH2;
487 }
488
498 template <class LhsEval>
499 LhsEval liquidDensityWaterH2_(const LhsEval& temperature,
500 const LhsEval& pl,
501 const LhsEval& xlH2) const
502 {
503 // molar masses
504 Scalar M_H2 = H2::molarMass();
505 Scalar M_H2O = H2O::molarMass();
506
507 // density of pure water
508 const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl, extrapolate);
509
510 // (apparent) molar volume of H2, Eq. (14) in Li et al. (2018)
511 const LhsEval& A1 = 51.1904 - 0.208062 * temperature +
512 3.4427e-4 * temperature * temperature;
513 const LhsEval& A2 = -0.022;
514 // pressure in [MPa] and Vphi in [m3/mol] (from [cm3/mol])
515 const LhsEval& V_phi = (A1 + A2 * (pl / 1e6)) / 1e6;
516
517 // density of solution, Eq. (19) in Garcia (2001)
518 const LhsEval xlH2O = 1.0 - xlH2;
519 const LhsEval& M_T = M_H2O * xlH2O + M_H2 * xlH2;
520 const LhsEval& rho_aq = 1 / (xlH2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
521
522 return rho_aq;
523 }
524
532 template <class LhsEval>
533 LhsEval convertRsToXoG_(const LhsEval& Rs, unsigned regionIdx) const
534 {
535 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
536 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
537
538 const LhsEval& rho_oG = Rs*rho_gRef;
539 return rho_oG/(rho_oRef + rho_oG);
540 }
541
548 template <class LhsEval>
549 LhsEval convertXoGToxoG_(const LhsEval& XoG, const LhsEval& salinity) const
550 {
551 Scalar M_H2 = H2::molarMass();
552 LhsEval M_Brine = Brine::molarMass(salinity);
553 return XoG*M_Brine / (M_H2*(1 - XoG) + XoG*M_Brine);
554 }
555
562 template <class LhsEval>
563 LhsEval convertxoGToXoG(const LhsEval& xoG, const LhsEval& salinity) const
564 {
565 Scalar M_H2 = H2::molarMass();
566 LhsEval M_Brine = Brine::molarMass(salinity);
567
568 return xoG*M_H2 / (xoG*(M_H2 - M_Brine) + M_Brine);
569 }
570
578 template <class LhsEval>
579 LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) const
580 {
581 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
582 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
583
584 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
585 }
586
594 template <class LhsEval>
595 LhsEval rsSat_(unsigned regionIdx,
596 const LhsEval& temperature,
597 const LhsEval& pressure,
598 const LhsEval& salinity) const
599 {
600 // Return Rs=0.0 if dissolution is disabled
601 if (!enableDissolution_)
602 return 0.0;
603
604 // calulate the equilibrium composition for the given temperature and pressure
605 LhsEval xlH2 = BinaryCoeffBrineH2::calculateMoleFractions(temperature, pressure,
606 salinity, extrapolate);
607
608 // normalize the phase compositions
609 xlH2 = max(0.0, min(1.0, xlH2));
610
611 return convertXoGToRs(convertxoGToXoG(xlH2, salinity), regionIdx);
612 }
613
614 template <class LhsEval>
615 static LhsEval liquidEnthalpyBrineH2_(const LhsEval& T,
616 const LhsEval& p,
617 const LhsEval& salinity,
618 const LhsEval& X_H2_w)
619 {
620 /* X_H2_w : mass fraction of H2 in brine */
621 /*NOTE: The heat of dissolution of H2 in brine is not included at the moment!*/
622
623 /*Numerical coefficents from PALLISER*/
624 static constexpr Scalar f[] = {
625 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
626 };
627
628 /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
629 static constexpr Scalar a[4][3] = {
630 { 9633.6, -4080.0, +286.49 },
631 { +166.58, +68.577, -4.6856 },
632 { -0.90963, -0.36524, +0.249667E-1 },
633 { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
634 };
635
636 LhsEval theta, h_NaCl;
637 LhsEval h_ls1, d_h;
638 LhsEval delta_h;
639 LhsEval hg, hw;
640
641 // Temperature in Celsius
642 theta = T - 273.15;
643
644 // Regularization
645 Scalar scalarTheta = scalarValue(theta);
646 Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
647
648 LhsEval S = salinity;
649 if (S > S_lSAT) {
650 S = S_lSAT;
651 }
652
653 hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
654
655 /*DAUBERT and DANNER*/
656 /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
657 +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
658
659 LhsEval m = 1E3/58.44 * S/(1-S);
660 d_h = 0;
661
662 for (int i = 0; i <=3 ; ++i) {
663 for (int j = 0; j <= 2; ++j) {
664 d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * pow(m, j);
665 }
666 }
667 /* heat of dissolution for halite according to Michaelides 1971 */
668 delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
669
670 /* Enthalpy of brine without H2 */
671 h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
672
673 /* enthalpy contribution of H2 gas (kJ/kg) */
674 hg = H2::gasEnthalpy(T, p, extrapolate) / 1E3;
675
676 /* Enthalpy of brine with dissolved H2 */
677 return (h_ls1 - X_H2_w*hw + hg*X_H2_w)*1E3; /*J/kg*/
678 }
679
680 template <class LhsEval>
681 const LhsEval salinityFromConcentration(unsigned regionIdx,
682 const LhsEval&T,
683 const LhsEval& P,
684 const LhsEval& saltConcentration) const
685 {
686 if (enableSaltConcentration_) {
687 // Convert concentration [kg/m³] to mass fraction [kg_salt/kg_solution].
688 // First approximation using pure water density
689 const LhsEval rho_w = H2O::liquidDensity(T, P, true);
690 const LhsEval S_approx = saltConcentration / rho_w;
691 // Improved estimate using Batzle-Wang brine density
692 const LhsEval rho_brine = Brine::liquidDensity(T, P, S_approx, rho_w);
693 return saltConcentration / rho_brine;
694 }
695
696 return salinity(regionIdx);
697 }
698
699 std::vector<Scalar> brineReferenceDensity_{};
700 std::vector<Scalar> h2ReferenceDensity_{};
701 std::vector<Scalar> salinity_{};
702 bool enableDissolution_ = true;
703 bool enableSaltConcentration_ = false;
704}; // end class BrineH2Pvt
705
706} // end namespace Opm
707
708#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:206
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:251
void setEnableSaltConcentration(bool yesno)
Specify whether the PVT model should consider salt concentration from the fluidstate or a fixed salin...
Definition BrineH2Pvt.hpp:111
void initEnd()
Finish initializing the oil phase PVT properties.
Definition BrineH2Pvt.hpp:92
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:321
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned, unsigned regionIdx=0) const
Diffusion coefficient of H2 in water.
Definition BrineH2Pvt.hpp:385
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:234
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:147
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:178
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:334
void setEnableDissolvedGas(bool yesno)
Specify whether the PVT model should consider that the H2 component can dissolve in the brine phase.
Definition BrineH2Pvt.hpp:102
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:165
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:266
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition BrineH2Pvt.hpp:117
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:287
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:362
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefBrine, Scalar rhoRefH2, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition BrineH2Pvt.cpp:126
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:217
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:306
Scalar hVap(unsigned) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition BrineH2Pvt.hpp:123
void initFromState(const EclipseState &eclState, const Schedule &)
Initialize the parameters for Brine-H2 system using an ECL deck.
Definition BrineH2Pvt.cpp:55
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:192
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:348
static Scalar molarMass()
Definition Component.hpp:93
Definition EclipseState.hpp:66
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