27#ifndef OPM_BRINE_H2_PVT_HPP
28#define OPM_BRINE_H2_PVT_HPP
53template <
class Scalar>
56 static const bool extrapolate =
true;
65 explicit BrineH2Pvt() =
default;
67 explicit BrineH2Pvt(
const std::vector<Scalar>& salinity,
68 Scalar T_ref = 288.71,
69 Scalar P_ref = 101325);
81 void setVapPars(
const Scalar,
const Scalar)
107 { enableDissolution_ = yesno; }
116 { enableSaltConcentration_ = yesno; }
122 {
return brineReferenceDensity_.size(); }
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
137 const Evaluation
salinity = salinityFromConcentration(regionIdx, temperature,
138 pressure, saltConcentration);
139 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
140 return liquidEnthalpyBrineH2_(temperature,
144 - pressure / density_(regionIdx, temperature, pressure, Rs,
salinity);
150 template <
class Evaluation>
152 const Evaluation& temperature,
153 const Evaluation& pressure,
154 const Evaluation& Rs)
const
156 const Evaluation xlH2 = convertRsToXoG_(Rs,regionIdx);
157 return liquidEnthalpyBrineH2_(temperature,
159 Evaluation(salinity_[regionIdx]),
161 - pressure / density_(regionIdx, temperature, pressure,
162 Rs, Evaluation(salinity_[regionIdx]));
168 template <
class Evaluation>
170 const Evaluation& temperature,
171 const Evaluation& pressure,
172 const Evaluation& )
const
181 template <
class Evaluation>
183 const Evaluation& temperature,
184 const Evaluation& pressure,
185 const Evaluation& saltConcentration)
const
187 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
188 pressure, saltConcentration);
195 template <
class Evaluation>
197 const Evaluation& temperature,
198 const Evaluation& pressure,
200 const Evaluation& saltConcentration)
const
209 template <
class Evaluation>
211 const Evaluation& temperature,
212 const Evaluation& pressure)
const
220 template <
class Evaluation>
222 const Evaluation& temperature,
223 const Evaluation& pressure,
224 const Evaluation& saltconcentration)
const
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];
237 template <
class Evaluation>
239 const Evaluation& temperature,
240 const Evaluation& pressure,
241 const Evaluation& Rs,
242 const Evaluation& saltConcentration)
const
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];
254 template <
class Evaluation>
256 const Evaluation& temperature,
257 const Evaluation& pressure,
258 const Evaluation& Rs)
const
260 return (1.0 - convertRsToXoG_(Rs, regionIdx))
261 * density_(regionIdx, temperature, pressure, Rs, Evaluation(salinity_[regionIdx]))
262 / brineReferenceDensity_[regionIdx];
268 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
269 std::pair<LhsEval, LhsEval>
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());
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);
290 template <
class Evaluation>
292 const Evaluation& temperature,
293 const Evaluation& pressure)
const
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];
309 template <
class Evaluation>
312 const Evaluation& )
const
314 throw std::runtime_error(
"Saturation pressure for the Brine-H2 PVT module "
315 "has not been implemented yet!");
324 template <
class Evaluation>
328 const Evaluation& )
const
330 throw std::runtime_error(
"Saturation pressure for the Brine-H2 PVT module "
331 "has not been implemented yet!");
337 template <
class Evaluation>
339 const Evaluation& temperature,
340 const Evaluation& pressure,
342 const Evaluation& )
const
345 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
351 template <
class Evaluation>
353 const Evaluation& temperature,
354 const Evaluation& pressure,
355 const Evaluation& saltConcentration)
const
357 const Evaluation salinity = salinityFromConcentration(regionIdx, temperature,
358 pressure, saltConcentration);
359 return rsSat_(regionIdx, temperature, pressure, salinity);
365 template <
class Evaluation>
367 const Evaluation& temperature,
368 const Evaluation& pressure)
const
370 return rsSat_(regionIdx, temperature, pressure, Evaluation(salinity_[regionIdx]));
373 Scalar oilReferenceDensity(
unsigned regionIdx)
const
374 {
return brineReferenceDensity_[regionIdx]; }
376 Scalar waterReferenceDensity(
unsigned regionIdx)
const
377 {
return brineReferenceDensity_[regionIdx]; }
379 Scalar gasReferenceDensity(
unsigned regionIdx)
const
380 {
return h2ReferenceDensity_[regionIdx]; }
382 Scalar
salinity(
unsigned regionIdx)
const
383 {
return salinity_[regionIdx]; }
388 template <
class Evaluation>
390 const Evaluation& pressure,
398 const Scalar vm = 28.45;
399 const Scalar sigma = 2.96 * 1e-8;
400 const Scalar avogadro = 6.022e23;
401 const Scalar alpha = sigma / pow((vm / avogadro), 1 / 3);
402 const Scalar lambda = 1.729;
405 Evaluation(salinity_[0])) * 1e3;
408 const Evaluation D_pure = ((4.8e-7 * temperature) / pow(mu_pure, alpha)) *
409 pow((1 + pow(lambda, 2)) / vm, 0.6);
415 const Evaluation log_D_brine = log10(D_pure) - 0.637 * log10(mu_brine / mu_pure);
417 return pow(Evaluation(10), log_D_brine) * 1e-4;
428 template <
class LhsEval>
429 LhsEval density_(
unsigned regionIdx,
430 const LhsEval& temperature,
431 const LhsEval& pressure,
436 LhsEval xlH2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx),
salinity);
439 LhsEval result = liquidDensity_(temperature,
456 template <
class LhsEval>
457 LhsEval liquidDensity_(
const LhsEval& T,
463 Valgrind::CheckDefined(T);
464 Valgrind::CheckDefined(pl);
465 Valgrind::CheckDefined(xlH2);
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);
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);
486 const LhsEval& rho_lH2 = liquidDensityWaterH2_(T, pl, xlH2);
487 const LhsEval& contribH2 = rho_lH2 - rho_pure;
489 return rho_brine + contribH2;
501 template <
class LhsEval>
502 LhsEval liquidDensityWaterH2_(
const LhsEval& temperature,
504 const LhsEval& xlH2)
const
514 const LhsEval& A1 = 51.1904 - 0.208062 * temperature +
515 3.4427e-4 * temperature * temperature;
516 const LhsEval& A2 = -0.022;
518 const LhsEval& V_phi = (A1 + A2 * (pl / 1e6)) / 1e6;
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));
535 template <
class LhsEval>
536 LhsEval convertRsToXoG_(
const LhsEval& Rs,
unsigned regionIdx)
const
538 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
539 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
541 const LhsEval& rho_oG = Rs*rho_gRef;
542 return rho_oG/(rho_oRef + rho_oG);
551 template <
class LhsEval>
552 LhsEval convertXoGToxoG_(
const LhsEval& XoG,
const LhsEval& salinity)
const
556 return XoG*M_Brine / (M_H2*(1 - XoG) + XoG*M_Brine);
565 template <
class LhsEval>
566 LhsEval convertxoGToXoG(
const LhsEval& xoG,
const LhsEval& salinity)
const
571 return xoG*M_H2 / (xoG*(M_H2 - M_Brine) + M_Brine);
581 template <
class LhsEval>
582 LhsEval convertXoGToRs(
const LhsEval& XoG,
unsigned regionIdx)
const
584 Scalar rho_oRef = brineReferenceDensity_[regionIdx];
585 Scalar rho_gRef = h2ReferenceDensity_[regionIdx];
587 return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
597 template <
class LhsEval>
598 LhsEval rsSat_(
unsigned regionIdx,
599 const LhsEval& temperature,
600 const LhsEval& pressure,
601 const LhsEval& salinity)
const
604 if (!enableDissolution_)
609 salinity, extrapolate);
612 xlH2 = max(0.0, min(1.0, xlH2));
614 return convertXoGToRs(convertxoGToXoG(xlH2, salinity), regionIdx);
617 template <
class LhsEval>
618 static LhsEval liquidEnthalpyBrineH2_(
const LhsEval& T,
620 const LhsEval& salinity,
621 const LhsEval& X_H2_w)
627 static constexpr Scalar f[] = {
628 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
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 }
639 LhsEval theta, h_NaCl;
648 Scalar scalarTheta = scalarValue(theta);
649 Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
651 LhsEval S = salinity;
659 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;
662 LhsEval m = 1E3/58.44 * S/(1-S);
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);
671 delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
674 h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h;
680 return (h_ls1 - X_H2_w*hw + hg*X_H2_w)*1E3;
683 template <
class LhsEval>
684 const LhsEval salinityFromConcentration(
unsigned regionIdx,
687 const LhsEval& saltConcentration)
const
689 if (enableSaltConcentration_) {
696 std::vector<Scalar> brineReferenceDensity_{};
697 std::vector<Scalar> h2ReferenceDensity_{};
698 std::vector<Scalar> salinity_{};
699 bool enableDissolution_ =
true;
700 bool enableSaltConcentration_ =
false;
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.
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