https://mooseframework.inl.gov
PorousFlowBrineCO2.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "PorousFlowBrineCO2.h"
11 #include "BrineFluidProperties.h"
13 #include "MathUtils.h"
14 #include "Conversion.h"
15 
16 registerMooseObject("PorousFlowApp", PorousFlowBrineCO2);
17 
20 {
22  params.addRequiredParam<UserObjectName>("brine_fp", "The name of the user object for brine");
23  params.addRequiredParam<UserObjectName>("co2_fp", "The name of the user object for CO2");
24  params.addParam<unsigned int>("salt_component", 2, "The component number of salt");
25  params.addClassDescription("Fluid state class for brine and CO2");
26  return params;
27 }
28 
31  _salt_component(getParam<unsigned int>("salt_component")),
32  _brine_fp(getUserObject<BrineFluidProperties>("brine_fp")),
33  _co2_fp(getUserObject<SinglePhaseFluidProperties>("co2_fp")),
34  _Mh2o(_brine_fp.molarMassH2O()),
35  _invMh2o(1.0 / _Mh2o),
36  _Mco2(_co2_fp.molarMass()),
37  _Mnacl(_brine_fp.molarMassNaCl()),
38  _Rbar(_R * 10.0),
39  _Tlower(372.15),
40  _Tupper(382.15),
41  _Zmin(1.0e-4),
42  _co2_henry(_co2_fp.henryCoefficients())
43 {
44  // Check that the correct FluidProperties UserObjects have been provided
45  if (_co2_fp.fluidName() != "co2")
46  paramError("co2_fp", "A valid CO2 FluidProperties UserObject must be provided");
47 
48  if (_brine_fp.fluidName() != "brine")
49  paramError("brine_fp", "A valid Brine FluidProperties UserObject must be provided");
50 
51  // Set the number of phases and components, and their indexes
52  _num_phases = 2;
53  _num_components = 3;
56 
57  // Check that _aqueous_phase_number is <= total number of phases
59  paramError("liquid_phase_number",
60  "This value is larger than the possible number of phases ",
61  _num_phases);
62 
63  // Check that _aqueous_fluid_component is <= total number of fluid components
65  paramError("liquid_fluid_component",
66  "This value is larger than the possible number of fluid components",
68 
69  // Check that the salt component index is not identical to the liquid fluid component
71  paramError(
72  "salt_component",
73  "The value provided must be different from the value entered in liquid_fluid_component");
74 
75  // Check that _salt_component is <= total number of fluid components
77  paramError("salt_component",
78  "The value provided is larger than the possible number of fluid components",
80 
82 }
83 
84 std::string
86 {
87  return "brine-co2";
88 }
89 
90 void
92  Real temperature,
93  Real Xnacl,
94  Real Z,
95  unsigned int qp,
96  std::vector<FluidStateProperties> & fsp) const
97 {
98  // Make AD versions of primary variables then call AD thermophysicalProperties()
99  ADReal p = pressure;
100  Moose::derivInsert(p.derivatives(), _pidx, 1.0);
101  ADReal T = temperature;
102  Moose::derivInsert(T.derivatives(), _Tidx, 1.0);
103  ADReal Zco2 = Z;
104  Moose::derivInsert(Zco2.derivatives(), _Zidx, 1.0);
105  ADReal X = Xnacl;
106  Moose::derivInsert(X.derivatives(), _Xidx, 1.0);
107 
108  thermophysicalProperties(p, T, X, Zco2, qp, fsp);
109 }
110 
111 void
113  const ADReal & temperature,
114  const ADReal & Xnacl,
115  const ADReal & Z,
116  unsigned int qp,
117  std::vector<FluidStateProperties> & fsp) const
118 {
121 
122  // Check whether the input temperature is within the region of validity
123  checkVariables(pressure.value(), temperature.value());
124 
125  // Clear all of the FluidStateProperties data
127 
128  FluidStatePhaseEnum phase_state;
129  massFractions(pressure, temperature, Xnacl, Z, phase_state, fsp);
130 
131  switch (phase_state)
132  {
134  {
135  // Set the gas saturations
136  gas.saturation = 1.0;
137 
138  // Calculate gas properties
140 
141  break;
142  }
143 
145  {
146  // Calculate the liquid properties
147  const ADReal liquid_pressure = pressure - _pc.capillaryPressure(1.0, qp);
148  liquidProperties(liquid_pressure, temperature, Xnacl, fsp);
149 
150  break;
151  }
152 
154  {
155  // Calculate the gas and liquid properties in the two phase region
156  twoPhaseProperties(pressure, temperature, Xnacl, Z, qp, fsp);
157 
158  break;
159  }
160  }
161 
162  // Liquid saturations can now be set
163  liquid.saturation = 1.0 - gas.saturation;
164 
165  // Save pressures to FluidStateProperties object
166  gas.pressure = pressure;
167  liquid.pressure = pressure - _pc.capillaryPressure(liquid.saturation, qp);
168 }
169 
170 void
172  const ADReal & temperature,
173  const ADReal & Xnacl,
174  const ADReal & Z,
175  FluidStatePhaseEnum & phase_state,
176  std::vector<FluidStateProperties> & fsp) const
177 {
180 
181  ADReal Xco2 = 0.0;
182  ADReal Yh2o = 0.0;
183  ADReal Yco2 = 0.0;
184 
185  // If the amount of CO2 is less than the smallest solubility, then all CO2 will
186  // be dissolved, and the equilibrium mass fractions do not need to be computed
187  if (Z.value() < _Zmin)
188  phase_state = FluidStatePhaseEnum::LIQUID;
189 
190  else
191  {
192  // Equilibrium mass fraction of CO2 in liquid and H2O in gas phases
193  equilibriumMassFractions(pressure, temperature, Xnacl, Xco2, Yh2o);
194 
195  Yco2 = 1.0 - Yh2o;
196 
197  // Determine which phases are present based on the value of z
198  phaseState(Z.value(), Xco2.value(), Yco2.value(), phase_state);
199  }
200 
201  // The equilibrium mass fractions calculated above are only correct in the two phase
202  // state. If only liquid or gas phases are present, the mass fractions are given by
203  // the total mass fraction z
204  ADReal Xh2o = 0.0;
205 
206  switch (phase_state)
207  {
209  {
210  Xco2 = Z;
211  Yco2 = 0.0;
212  Xh2o = 1.0 - Z;
213  Yh2o = 0.0;
214  break;
215  }
216 
218  {
219  Xco2 = 0.0;
220  Yco2 = Z;
221  Yh2o = 1.0 - Z;
222  break;
223  }
224 
226  {
227  // Keep equilibrium mass fractions
228  Xh2o = 1.0 - Xco2;
229  break;
230  }
231  }
232 
233  // Save the mass fractions in the FluidStateProperties object
235  liquid.mass_fraction[_gas_fluid_component] = Xco2;
236  liquid.mass_fraction[_salt_component] = Xnacl;
239 }
240 
241 void
243  const ADReal & temperature,
244  std::vector<FluidStateProperties> & fsp) const
245 {
247 
248  // Gas density, viscosity and enthalpy are approximated with pure CO2 - no correction due
249  // to the small amount of water vapor is made
250  ADReal co2_density, co2_viscosity;
251  _co2_fp.rho_mu_from_p_T(pressure, temperature, co2_density, co2_viscosity);
252 
253  ADReal co2_enthalpy = _co2_fp.h_from_p_T(pressure, temperature);
254 
255  // Save the values to the FluidStateProperties object. Note that derivatives wrt z are 0
256  gas.density = co2_density;
257  gas.viscosity = co2_viscosity;
258  gas.enthalpy = co2_enthalpy;
259 
260  mooseAssert(gas.density.value() > 0.0, "Gas density must be greater than zero");
261  gas.internal_energy = gas.enthalpy - pressure / gas.density;
262 }
263 
264 void
266  const ADReal & temperature,
267  const ADReal & Xnacl,
268  std::vector<FluidStateProperties> & fsp) const
269 {
271 
272  // The liquid density includes the density increase due to dissolved CO2
273  const ADReal brine_density = _brine_fp.rho_from_p_T_X(pressure, temperature, Xnacl);
274 
275  // Mass fraction of CO2 in liquid phase
276  const ADReal Xco2 = liquid.mass_fraction[_gas_fluid_component];
277 
278  // The liquid density
279  const ADReal co2_partial_density = partialDensityCO2(temperature);
280 
281  const ADReal liquid_density = 1.0 / (Xco2 / co2_partial_density + (1.0 - Xco2) / brine_density);
282 
283  // Assume that liquid viscosity is just the brine viscosity
284  const ADReal liquid_viscosity = _brine_fp.mu_from_p_T_X(pressure, temperature, Xnacl);
285 
286  // Liquid enthalpy (including contribution due to the enthalpy of dissolution)
287  const ADReal brine_enthalpy = _brine_fp.h_from_p_T_X(pressure, temperature, Xnacl);
288 
289  // Enthalpy of CO2
290  const ADReal co2_enthalpy = _co2_fp.h_from_p_T(pressure, temperature);
291 
292  // Enthalpy of dissolution
294 
295  const ADReal liquid_enthalpy = (1.0 - Xco2) * brine_enthalpy + Xco2 * (co2_enthalpy + hdis);
296 
297  // Save the values to the FluidStateProperties object
298  liquid.density = liquid_density;
299  liquid.viscosity = liquid_viscosity;
300  liquid.enthalpy = liquid_enthalpy;
301 
302  mooseAssert(liquid.density.value() > 0.0, "Liquid density must be greater than zero");
303  liquid.internal_energy = liquid.enthalpy - pressure / liquid.density;
304 }
305 
306 ADReal
308  const ADReal & temperature,
309  const ADReal & Xnacl,
310  const ADReal & Z,
311  std::vector<FluidStateProperties> & fsp) const
312 {
313  auto & gas = fsp[_gas_phase_number];
314  auto & liquid = fsp[_aqueous_fluid_component];
315 
316  // Approximate liquid density as saturation isn't known yet, by using the gas
317  // pressure rather than the liquid pressure. This does result in a small error
318  // in the calculated saturation, but this is below the error associated with
319  // the correlations. A more accurate saturation could be found iteraviely,
320  // at the cost of increased computational expense
321 
322  // Gas density
323  const ADReal gas_density = _co2_fp.rho_from_p_T(pressure, temperature);
324 
325  // Approximate liquid density as saturation isn't known yet
326  const ADReal brine_density = _brine_fp.rho_from_p_T_X(pressure, temperature, Xnacl);
327 
328  // Mass fraction of CO2 in liquid phase
329  const ADReal Xco2 = liquid.mass_fraction[_gas_fluid_component];
330 
331  // The liquid density
332  const ADReal co2_partial_density = partialDensityCO2(temperature);
333 
334  const ADReal liquid_density = 1.0 / (Xco2 / co2_partial_density + (1.0 - Xco2) / brine_density);
335 
336  const ADReal Yco2 = gas.mass_fraction[_gas_fluid_component];
337 
338  // Set mass equilibrium constants used in the calculation of vapor mass fraction
339  const ADReal K0 = Yco2 / Xco2;
340  const ADReal K1 = (1.0 - Yco2) / (1.0 - Xco2);
341  const ADReal vapor_mass_fraction = vaporMassFraction(Z, K0, K1);
342 
343  // The gas saturation in the two phase case
344  const ADReal saturation = vapor_mass_fraction * liquid_density /
345  (gas_density + vapor_mass_fraction * (liquid_density - gas_density));
346 
347  return saturation;
348 }
349 
350 void
352  const ADReal & temperature,
353  const ADReal & Xnacl,
354  const ADReal & Z,
355  unsigned int qp,
356  std::vector<FluidStateProperties> & fsp) const
357 {
358  auto & gas = fsp[_gas_phase_number];
359 
360  // Calculate all of the gas phase properties, as these don't depend on saturation
362 
363  // The gas saturation in the two phase case
364  gas.saturation = saturation(pressure, temperature, Xnacl, Z, fsp);
365 
366  // The liquid pressure and properties can now be calculated
367  const ADReal liquid_pressure = pressure - _pc.capillaryPressure(1.0 - gas.saturation, qp);
368  liquidProperties(liquid_pressure, temperature, Xnacl, fsp);
369 }
370 
371 void
373  const ADReal & temperature,
374  const ADReal & Xnacl,
375  ADReal & Xco2,
376  ADReal & Yh2o) const
377 {
378  // Mole fractions at equilibrium
379  ADReal xco2, yh2o;
380  equilibriumMoleFractions(pressure, temperature, Xnacl, xco2, yh2o);
381 
382  // The mass fraction of H2O in gas (assume no salt in gas phase) and derivatives
383  // wrt p, T, and X
384  Yh2o = yh2o * _Mh2o / (yh2o * _Mh2o + (1.0 - yh2o) * _Mco2);
385 
386  // NaCl molality (mol/kg)
387  const ADReal mnacl = Xnacl / (1.0 - Xnacl) / _Mnacl;
388 
389  // The molality of CO2 in 1kg of H2O
390  const ADReal mco2 = xco2 * (2.0 * mnacl + _invMh2o) / (1.0 - xco2);
391  // The mass fraction of CO2 in brine is then
392  const ADReal denominator = (1.0 + mnacl * _Mnacl + mco2 * _Mco2);
393  Xco2 = mco2 * _Mco2 / denominator;
394 }
395 
396 void
398  const ADReal & temperature,
399  const ADReal & co2_density,
400  ADReal & fco2,
401  ADReal & fh2o) const
402 {
403  if (temperature.value() > 373.15)
404  mooseError(name(),
405  ": fugacityCoefficientsLowTemp() is not valid for T > 373.15K. Use "
406  "fugacityCoefficientsHighTemp() instead");
407 
408  // Need pressure in bar
409  const ADReal pbar = pressure * 1.0e-5;
410 
411  // Molar volume in cm^3/mol
412  const ADReal V = _Mco2 / co2_density * 1.0e6;
413 
414  // Redlich-Kwong parameters
415  const ADReal aCO2 = 7.54e7 - 4.13e4 * temperature;
416  const Real bCO2 = 27.8;
417  const Real aCO2H2O = 7.89e7;
418  const Real bH2O = 18.18;
419 
420  const ADReal t15 = std::pow(temperature, 1.5);
421 
422  // The fugacity coefficients for H2O and CO2
423  auto lnPhi = [V, aCO2, bCO2, t15, this](ADReal a, ADReal b)
424  {
425  return std::log(V / (V - bCO2)) + b / (V - bCO2) -
426  2.0 * a / (_Rbar * t15 * bCO2) * std::log((V + bCO2) / V) +
427  aCO2 * b / (_Rbar * t15 * bCO2 * bCO2) * (std::log((V + bCO2) / V) - bCO2 / (V + bCO2));
428  };
429 
430  const ADReal lnPhiH2O = lnPhi(aCO2H2O, bH2O) - std::log(pbar * V / (_Rbar * temperature));
431  const ADReal lnPhiCO2 = lnPhi(aCO2, bCO2) - std::log(pbar * V / (_Rbar * temperature));
432 
433  fh2o = std::exp(lnPhiH2O);
434  fco2 = std::exp(lnPhiCO2);
435 }
436 
437 void
439  const ADReal & temperature,
440  const ADReal & co2_density,
441  const ADReal & xco2,
442  const ADReal & yh2o,
443  ADReal & fco2,
444  ADReal & fh2o) const
445 {
446  if (temperature.value() <= 373.15)
447  mooseError(name(),
448  ": fugacityCoefficientsHighTemp() is not valid for T <= 373.15K. Use "
449  "fugacityCoefficientsLowTemp() instead");
450 
451  fh2o = fugacityCoefficientH2OHighTemp(pressure, temperature, co2_density, xco2, yh2o);
452  fco2 = fugacityCoefficientCO2HighTemp(pressure, temperature, co2_density, xco2, yh2o);
453 }
454 
455 ADReal
457  const ADReal & temperature,
458  const ADReal & co2_density,
459  const ADReal & xco2,
460  const ADReal & yh2o) const
461 {
462  // Need pressure in bar
463  const ADReal pbar = pressure * 1.0e-5;
464  // Molar volume in cm^3/mol
465  const ADReal V = _Mco2 / co2_density * 1.0e6;
466 
467  // Redlich-Kwong parameters
468  const ADReal yco2 = 1.0 - yh2o;
469  const ADReal xh2o = 1.0 - xco2;
470 
471  const ADReal aCO2 = 8.008e7 - 4.984e4 * temperature;
472  const ADReal aH2O = 1.337e8 - 1.4e4 * temperature;
473  const Real bCO2 = 28.25;
474  const Real bH2O = 15.7;
475  const ADReal KH2OCO2 = 1.427e-2 - 4.037e-4 * temperature;
476  const ADReal KCO2H2O = 0.4228 - 7.422e-4 * temperature;
477  const ADReal kH2OCO2 = KH2OCO2 * yh2o + KCO2H2O * yco2;
478  const ADReal kCO2H2O = KCO2H2O * yh2o + KH2OCO2 * yco2;
479 
480  const ADReal aH2OCO2 = std::sqrt(aCO2 * aH2O) * (1.0 - kH2OCO2);
481  const ADReal aCO2H2O = std::sqrt(aCO2 * aH2O) * (1.0 - kCO2H2O);
482 
483  const ADReal amix = yh2o * yh2o * aH2O + yh2o * yco2 * (aH2OCO2 + aCO2H2O) + yco2 * yco2 * aCO2;
484  const ADReal bmix = yh2o * bH2O + yco2 * bCO2;
485 
486  const ADReal t15 = std::pow(temperature, 1.5);
487 
488  ADReal lnPhiH2O = bH2O / bmix * (pbar * V / (_Rbar * temperature) - 1.0) -
489  std::log(pbar * (V - bmix) / (_Rbar * temperature));
490  ADReal term3 = (2.0 * yh2o * aH2O + yco2 * (aH2OCO2 + aCO2H2O) -
491  yh2o * yco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O) * (yh2o - yco2) +
492  xh2o * xco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O)) /
493  amix;
494  term3 -= bH2O / bmix;
495  term3 *= amix / (bmix * _Rbar * t15) * std::log(V / (V + bmix));
496  lnPhiH2O += term3;
497 
498  return std::exp(lnPhiH2O);
499 }
500 
501 ADReal
503  const ADReal & temperature,
504  const ADReal & co2_density,
505  const ADReal & xco2,
506  const ADReal & yh2o) const
507 {
508  // Need pressure in bar
509  const ADReal pbar = pressure * 1.0e-5;
510  // Molar volume in cm^3/mol
511  const ADReal V = _Mco2 / co2_density * 1.0e6;
512 
513  // Redlich-Kwong parameters
514  const ADReal yco2 = 1.0 - yh2o;
515  const ADReal xh2o = 1.0 - xco2;
516 
517  const ADReal aCO2 = 8.008e7 - 4.984e4 * temperature;
518  const ADReal aH2O = 1.337e8 - 1.4e4 * temperature;
519  const Real bCO2 = 28.25;
520  const Real bH2O = 15.7;
521  const ADReal KH2OCO2 = 1.427e-2 - 4.037e-4 * temperature;
522  const ADReal KCO2H2O = 0.4228 - 7.422e-4 * temperature;
523  const ADReal kH2OCO2 = KH2OCO2 * yh2o + KCO2H2O * yco2;
524  const ADReal kCO2H2O = KCO2H2O * yh2o + KH2OCO2 * yco2;
525 
526  const ADReal aH2OCO2 = std::sqrt(aCO2 * aH2O) * (1.0 - kH2OCO2);
527  const ADReal aCO2H2O = std::sqrt(aCO2 * aH2O) * (1.0 - kCO2H2O);
528 
529  const ADReal amix = yh2o * yh2o * aH2O + yh2o * yco2 * (aH2OCO2 + aCO2H2O) + yco2 * yco2 * aCO2;
530  const ADReal bmix = yh2o * bH2O + yco2 * bCO2;
531 
532  const ADReal t15 = std::pow(temperature, 1.5);
533 
534  ADReal lnPhiCO2 = bCO2 / bmix * (pbar * V / (_Rbar * temperature) - 1.0) -
535  std::log(pbar * (V - bmix) / (_Rbar * temperature));
536 
537  ADReal term3 = (2.0 * yco2 * aCO2 + yh2o * (aH2OCO2 + aCO2H2O) -
538  yh2o * yco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O) * (yh2o - yco2) +
539  xh2o * xco2 * std::sqrt(aH2O * aCO2) * (kCO2H2O - kH2OCO2)) /
540  amix;
541 
542  lnPhiCO2 += (term3 - bCO2 / bmix) * amix / (bmix * _Rbar * t15) * std::log(V / (V + bmix));
543 
544  return std::exp(lnPhiCO2);
545 }
546 
547 ADReal
549 {
550  if (temperature.value() <= 373.15)
551  return 1.0;
552  else
553  {
554  const ADReal Tref = temperature - 373.15;
555  const ADReal xh2o = 1.0 - xco2;
556  const ADReal Am = -3.084e-2 * Tref + 1.927e-5 * Tref * Tref;
557 
558  return std::exp((Am - 2.0 * Am * xh2o) * xco2 * xco2);
559  }
560 }
561 
562 ADReal
564 {
565  if (temperature.value() <= 373.15)
566  return 1.0;
567  else
568  {
569  const ADReal Tref = temperature - 373.15;
570  const ADReal xh2o = 1.0 - xco2;
571  const ADReal Am = -3.084e-2 * Tref + 1.927e-5 * Tref * Tref;
572 
573  return std::exp(2.0 * Am * xco2 * xh2o * xh2o);
574  }
575 }
576 
577 ADReal
579  const ADReal & temperature,
580  const ADReal & Xnacl) const
581 {
582  // Need pressure in bar
583  const ADReal pbar = pressure * 1.0e-5;
584  // Need NaCl molality (mol/kg)
585  const ADReal mnacl = Xnacl / (1.0 - Xnacl) / _Mnacl;
586 
587  const ADReal lambda = -0.411370585 + 6.07632013e-4 * temperature + 97.5347708 / temperature -
588  0.0237622469 * pbar / temperature +
589  0.0170656236 * pbar / (630.0 - temperature) +
590  1.41335834e-5 * temperature * std::log(pbar);
591 
592  const ADReal xi = 3.36389723e-4 - 1.9829898e-5 * temperature +
593  2.12220830e-3 * pbar / temperature -
594  5.24873303e-3 * pbar / (630.0 - temperature);
595 
596  return std::exp(2.0 * lambda * mnacl + xi * mnacl * mnacl);
597 }
598 
599 ADReal
601  const ADReal & Xnacl) const
602 {
603  // Need NaCl molality (mol/kg)
604  const ADReal mnacl = Xnacl / (1.0 - Xnacl) / _Mnacl;
605 
606  const ADReal T2 = temperature * temperature;
607  const ADReal T3 = temperature * T2;
608 
609  const ADReal lambda = 2.217e-4 * temperature + 1.074 / temperature + 2648.0 / T2;
610  const ADReal xi = 1.3e-5 * temperature - 20.12 / temperature + 5259.0 / T2;
611 
612  return (1.0 - mnacl / _invMh2o) * std::exp(2.0 * lambda * mnacl + xi * mnacl * mnacl);
613 }
614 
615 ADReal
617 {
618  // Uses temperature in Celsius
619  const ADReal Tc = temperature - _T_c2k;
620  const ADReal Tc2 = Tc * Tc;
621  const ADReal Tc3 = Tc2 * Tc;
622  const ADReal Tc4 = Tc3 * Tc;
623 
624  ADReal logK0H2O;
625 
626  if (Tc <= 99.0)
627  logK0H2O = -2.209 + 3.097e-2 * Tc - 1.098e-4 * Tc2 + 2.048e-7 * Tc3;
628 
629  else if (Tc > 99.0 && Tc < 109.0)
630  {
631  const ADReal Tint = (Tc - 99.0) / 10.0;
632  const ADReal Tint2 = Tint * Tint;
633  logK0H2O = -0.0204026 + 0.0152513 * Tint + 0.417565 * Tint2 - 0.278636 * Tint * Tint2;
634  }
635 
636  else // 109 <= Tc <= 300
637  logK0H2O = -2.1077 + 2.8127e-2 * Tc - 8.4298e-5 * Tc2 + 1.4969e-7 * Tc3 - 1.1812e-10 * Tc4;
638 
639  return std::pow(10.0, logK0H2O);
640 }
641 
642 ADReal
644 {
645  // Uses temperature in Celsius
646  const ADReal Tc = temperature - _T_c2k;
647  const ADReal Tc2 = Tc * Tc;
648  const ADReal Tc3 = Tc2 * Tc;
649 
650  ADReal logK0CO2;
651 
652  if (Tc <= 99.0)
653  logK0CO2 = 1.189 + 1.304e-2 * Tc - 5.446e-5 * Tc2;
654 
655  else if (Tc > 99.0 && Tc < 109.0)
656  {
657  const ADReal Tint = (Tc - 99.0) / 10.0;
658  const ADReal Tint2 = Tint * Tint;
659  logK0CO2 = 1.9462 + 2.25692e-2 * Tint - 9.49577e-3 * Tint2 - 6.77721e-3 * Tint * Tint2;
660  }
661 
662  else // 109 <= Tc <= 300
663  logK0CO2 = 1.668 + 3.992e-3 * Tc - 1.156e-5 * Tc2 + 1.593e-9 * Tc3;
664 
665  return std::pow(10.0, logK0CO2);
666 }
667 
668 void
670  const ADReal & temperature,
671  const ADReal & Xnacl,
672  ADReal & xco2,
673  ADReal & yh2o) const
674 {
675  if (temperature.value() <= _Tlower)
676  {
678  }
679  else if (temperature.value() > _Tlower && temperature.value() < _Tupper)
680  {
681  // Cubic polynomial in this regime
682  const Real Tint = (temperature.value() - _Tlower) / 10.0;
683 
684  // Equilibrium mole fractions and derivatives at the lower temperature
685  ADReal Tlower = _Tlower;
686  Moose::derivInsert(Tlower.derivatives(), _Tidx, 1.0);
687 
688  ADReal xco2_lower, yh2o_lower;
689  equilibriumMoleFractionsLowTemp(pressure, Tlower, Xnacl, xco2_lower, yh2o_lower);
690 
691  const Real dxco2_dT_lower = xco2_lower.derivatives()[_Tidx];
692  const Real dyh2o_dT_lower = yh2o_lower.derivatives()[_Tidx];
693 
694  // Equilibrium mole fractions and derivatives at the upper temperature
695  Real xco2_upper, yh2o_upper;
696  Real co2_density_upper = _co2_fp.rho_from_p_T(pressure.value(), _Tupper);
697 
699  pressure.value(), _Tupper, Xnacl.value(), co2_density_upper, xco2_upper, yh2o_upper);
700 
701  Real A, dA_dp, dA_dT, B, dB_dp, dB_dT, dB_dX;
702  funcABHighTemp(pressure.value(),
703  _Tupper,
704  Xnacl.value(),
705  co2_density_upper,
706  xco2_upper,
707  yh2o_upper,
708  A,
709  dA_dp,
710  dA_dT,
711  B,
712  dB_dp,
713  dB_dT,
714  dB_dX);
715 
716  const Real dyh2o_dT_upper =
717  ((1.0 - B) * dA_dT + (A - 1.0) * A * dB_dT) / (1.0 - A * B) / (1.0 - A * B);
718  const Real dxco2_dT_upper = dB_dT * (1.0 - yh2o_upper) - B * dyh2o_dT_upper;
719 
720  // The mole fractions in this regime are then found by interpolation
721  Real xco2r, yh2or, dxco2_dT, dyh2o_dT;
723  Tint, xco2_lower.value(), dxco2_dT_lower, xco2_upper, dxco2_dT_upper, xco2r, dxco2_dT);
725  Tint, yh2o_lower.value(), dyh2o_dT_lower, yh2o_upper, dyh2o_dT_upper, yh2or, dyh2o_dT);
726 
727  xco2 = xco2r;
728  Moose::derivInsert(xco2.derivatives(), _pidx, xco2_lower.derivatives()[_pidx]);
729  Moose::derivInsert(xco2.derivatives(), _Tidx, dxco2_dT);
730  Moose::derivInsert(xco2.derivatives(), _Xidx, xco2_lower.derivatives()[_Xidx]);
731 
732  yh2o = yh2or;
733  Moose::derivInsert(yh2o.derivatives(), _pidx, yh2o_lower.derivatives()[_pidx]);
734  Moose::derivInsert(yh2o.derivatives(), _Tidx, dyh2o_dT);
735  Moose::derivInsert(yh2o.derivatives(), _Xidx, yh2o_lower.derivatives()[_Xidx]);
736  }
737  else
738  {
739  // CO2 density and derivatives wrt pressure and temperature
740  const Real co2_density = _co2_fp.rho_from_p_T(pressure.value(), temperature.value());
741 
742  // Equilibrium mole fractions solved using iteration in this regime
743  Real xco2r, yh2or;
745  pressure.value(), temperature.value(), Xnacl.value(), co2_density, xco2r, yh2or);
746 
747  // Can use these in funcABHighTemp() to get derivatives analytically rather than by iteration
748  Real A, dA_dp, dA_dT, B, dB_dp, dB_dT, dB_dX;
749  funcABHighTemp(pressure.value(),
750  temperature.value(),
751  Xnacl.value(),
752  co2_density,
753  xco2r,
754  yh2or,
755  A,
756  dA_dp,
757  dA_dT,
758  B,
759  dB_dp,
760  dB_dT,
761  dB_dX);
762 
763  const Real dyh2o_dp =
764  ((1.0 - B) * dA_dp + (A - 1.0) * A * dB_dp) / (1.0 - A * B) / (1.0 - A * B);
765  const Real dxco2_dp = dB_dp * (1.0 - yh2or) - B * dyh2o_dp;
766 
767  const Real dyh2o_dT =
768  ((1.0 - B) * dA_dT + (A - 1.0) * A * dB_dT) / (1.0 - A * B) / (1.0 - A * B);
769  const Real dxco2_dT = dB_dT * (1.0 - yh2or) - B * dyh2o_dT;
770 
771  const Real dyh2o_dX = ((A - 1.0) * A * dB_dX) / (1.0 - A * B) / (1.0 - A * B);
772  const Real dxco2_dX = dB_dX * (1.0 - yh2or) - B * dyh2o_dX;
773 
774  xco2 = xco2r;
775  Moose::derivInsert(xco2.derivatives(), _pidx, dxco2_dp);
776  Moose::derivInsert(xco2.derivatives(), _Tidx, dxco2_dT);
777  Moose::derivInsert(xco2.derivatives(), _Xidx, dxco2_dX);
778 
779  yh2o = yh2or;
780  Moose::derivInsert(yh2o.derivatives(), _pidx, dyh2o_dp);
781  Moose::derivInsert(yh2o.derivatives(), _Tidx, dyh2o_dT);
782  Moose::derivInsert(yh2o.derivatives(), _Xidx, dyh2o_dX);
783  }
784 }
785 
786 void
788  const ADReal & temperature,
789  const ADReal & Xnacl,
790  ADReal & xco2,
791  ADReal & yh2o) const
792 {
793  if (temperature.value() > 373.15)
794  mooseError(name(),
795  ": equilibriumMoleFractionsLowTemp() is not valid for T > 373.15K. Use "
796  "equilibriumMoleFractions() instead");
797 
798  // CO2 density and derivatives wrt pressure and temperature
799  const ADReal co2_density = _co2_fp.rho_from_p_T(pressure, temperature);
800 
801  // Assume infinite dilution (yh20 = 0 and xco2 = 0) in low temperature regime
802  ADReal A, B;
803  funcABLowTemp(pressure, temperature, co2_density, A, B);
804 
805  // As the activity coefficient for CO2 in brine used in this regime isn't a 'true'
806  // activity coefficient, we instead calculate the molality of CO2 in water, then
807  // correct it for brine, and then calculate the mole fractions.
808  // The mole fraction in pure water is
809  const ADReal yh2ow = (1.0 - B) / (1.0 / A - B);
810  const ADReal xco2w = B * (1.0 - yh2ow);
811 
812  // Molality of CO2 in pure water
813  const ADReal mco2w = xco2w * _invMh2o / (1.0 - xco2w);
814  // Molality of CO2 in brine is then calculated using gamma
815  const ADReal gamma = activityCoefficient(pressure, temperature, Xnacl);
816  const ADReal mco2 = mco2w / gamma;
817 
818  // Need NaCl molality (mol/kg)
819  const ADReal mnacl = Xnacl / (1.0 - Xnacl) / _Mnacl;
820 
821  // Mole fractions of CO2 and H2O in liquid and gas phases
822  const ADReal total_moles = 2.0 * mnacl + _invMh2o + mco2;
823  xco2 = mco2 / total_moles;
824  yh2o = A * (1.0 - xco2 - 2.0 * mnacl / total_moles);
825 }
826 
827 void
829  const ADReal & temperature,
830  const ADReal & co2_density,
831  ADReal & A,
832  ADReal & B) const
833 {
834  if (temperature.value() > 373.15)
835  mooseError(name(),
836  ": funcABLowTemp() is not valid for T > 373.15K. Use funcABHighTemp() instead");
837 
838  // Pressure in bar
839  const ADReal pbar = pressure * 1.0e-5;
840 
841  // Reference pressure and partial molar volumes
842  const Real pref = 1.0;
843  const Real vCO2 = 32.6;
844  const Real vH2O = 18.1;
845 
846  const ADReal delta_pbar = pbar - pref;
847  const ADReal Rt = _Rbar * temperature;
848 
849  // Equilibrium constants
852 
853  // Fugacity coefficients
854  ADReal phiH2O, phiCO2;
855  fugacityCoefficientsLowTemp(pressure, temperature, co2_density, phiCO2, phiH2O);
856 
857  A = K0H2O / (phiH2O * pbar) * std::exp(delta_pbar * vH2O / Rt);
858  B = phiCO2 * pbar / (_invMh2o * K0CO2) * std::exp(-delta_pbar * vCO2 / Rt);
859 }
860 
861 void
863  Real temperature,
864  Real Xnacl,
865  Real co2_density,
866  Real xco2,
867  Real yh2o,
868  Real & A,
869  Real & B) const
870 {
871  if (temperature <= 373.15)
872  mooseError(name(),
873  ": funcABHighTemp() is not valid for T <= 373.15K. Use funcABLowTemp() instead");
874 
875  // Pressure in bar
876  const Real pbar = pressure * 1.0e-5;
877  // Temperature in C
878  const Real Tc = temperature - _T_c2k;
879 
880  // Reference pressure and partial molar volumes
881  const Real pref = -1.9906e-1 + 2.0471e-3 * Tc + 1.0152e-4 * Tc * Tc - 1.4234e-6 * Tc * Tc * Tc +
882  1.4168e-8 * Tc * Tc * Tc * Tc;
883  const Real vCO2 = 32.6 + 3.413e-2 * (Tc - 100.0);
884  const Real vH2O = 18.1 + 3.137e-2 * (Tc - 100.0);
885 
886  const Real delta_pbar = pbar - pref;
887  const Real Rt = _Rbar * temperature;
888 
889  // Equilibrium constants
890  // Use dummy ADReal temperature as derivatives aren't required
891  const ADReal T = temperature;
892  Real K0H2O = equilibriumConstantH2O(T).value();
893  Real K0CO2 = equilibriumConstantCO2(T).value();
894 
895  // Fugacity coefficients
896  // Use dummy ADReal variables as derivatives aren't required
897  const ADReal p = pressure;
898  const ADReal rhoco2 = co2_density;
899  const ADReal x = xco2;
900  const ADReal y = yh2o;
901 
902  ADReal phiH2O, phiCO2;
903  fugacityCoefficientsHighTemp(p, T, rhoco2, x, y, phiCO2, phiH2O);
904 
905  // Activity coefficients
906  const Real gammaH2O = activityCoefficientH2O(T, x).value();
907  const Real gammaCO2 = activityCoefficientCO2(T, x).value();
908 
909  // Activity coefficient for CO2 in brine
910  // Use dummy ADReal Xnacl as derivatives aren't required
911  const ADReal X = Xnacl;
912  const Real gamma = activityCoefficientHighTemp(T, X).value();
913 
914  A = K0H2O * gammaH2O / (phiH2O.value() * pbar) * std::exp(delta_pbar * vH2O / Rt);
915  B = phiCO2.value() * pbar / (_invMh2o * K0CO2 * gamma * gammaCO2) *
916  std::exp(-delta_pbar * vCO2 / Rt);
917 }
918 
919 void
921  Real temperature,
922  Real Xnacl,
923  Real co2_density,
924  Real xco2,
925  Real yh2o,
926  Real & A,
927  Real & dA_dp,
928  Real & dA_dT,
929  Real & B,
930  Real & dB_dp,
931  Real & dB_dT,
932  Real & dB_dX) const
933 {
934  funcABHighTemp(pressure, temperature, Xnacl, co2_density, xco2, yh2o, A, B);
935 
936  // Use finite differences for derivatives in the high temperature regime
937  const Real dp = 1.0e-2;
938  const Real dT = 1.0e-6;
939  const Real dX = 1.0e-8;
940 
941  Real A2, B2;
942  funcABHighTemp(pressure + dp, temperature, Xnacl, co2_density, xco2, yh2o, A2, B2);
943  dA_dp = (A2 - A) / dp;
944  dB_dp = (B2 - B) / dp;
945 
946  funcABHighTemp(pressure, temperature + dT, Xnacl, co2_density, xco2, yh2o, A2, B2);
947  dA_dT = (A2 - A) / dT;
948  dB_dT = (B2 - B) / dT;
949 
950  funcABHighTemp(pressure, temperature, Xnacl + dX, co2_density, xco2, yh2o, A2, B2);
951  dB_dX = (B2 - B) / dX;
952 }
953 
954 void
956  Real pressure, Real temperature, Real Xnacl, Real co2_density, Real & xco2, Real & yh2o) const
957 {
958  // Initial guess for yh2o and xco2 (from Spycher and Pruess (2010))
960  Real x = 0.009;
961 
962  // Need salt mass fraction in molality
963  const Real mnacl = Xnacl / (1.0 - Xnacl) / _Mnacl;
964 
965  // If y > 1, then just use y = 1, x = 0 (only a gas phase)
966  if (y >= 1.0)
967  {
968  y = 1.0;
969  x = 0.0;
970  }
971  else
972  {
973  // Residual function for Newton-Raphson
974  auto fy = [mnacl, this](Real y, Real A, Real B) {
975  return y -
976  (1.0 - B) * _invMh2o / ((1.0 / A - B) * (2.0 * mnacl + _invMh2o) + 2.0 * mnacl * B);
977  };
978 
979  // Derivative of fy wrt y
980  auto dfy = [mnacl, this](Real A, Real B, Real dA, Real dB)
981  {
982  const Real denominator = (1.0 / A - B) * (2.0 * mnacl + _invMh2o) + 2.0 * mnacl * B;
983  return 1.0 + _invMh2o * dB / denominator +
984  (1.0 - B) * _invMh2o *
985  (2.0 * mnacl * dB - (2.0 * mnacl + _invMh2o) * (dB + dA / A / A)) / denominator /
986  denominator;
987  };
988 
989  Real A, B;
990  Real dA, dB;
991  const Real dy = 1.0e-8;
992 
993  // Solve for yh2o using Newton-Raphson method
994  unsigned int iter = 0;
995  const Real tol = 1.0e-12;
996  const unsigned int max_its = 10;
997  funcABHighTemp(pressure, temperature, Xnacl, co2_density, x, y, A, B);
998 
999  while (std::abs(fy(y, A, B)) > tol)
1000  {
1001  funcABHighTemp(pressure, temperature, Xnacl, co2_density, x, y, A, B);
1002  // Finite difference derivatives of A and B wrt y
1003  funcABHighTemp(pressure, temperature, Xnacl, co2_density, x, y + dy, dA, dB);
1004  dA = (dA - A) / dy;
1005  dB = (dB - B) / dy;
1006 
1007  y = y - fy(y, A, B) / dfy(A, B, dA, dB);
1008 
1009  x = B * (1.0 - y);
1010 
1011  // Break if not converged and just use the value
1012  if (iter > max_its)
1013  break;
1014  }
1015  }
1016 
1017  yh2o = y;
1018  xco2 = x;
1019 }
1020 
1021 ADReal
1023 {
1024  // This correlation uses temperature in C
1025  const ADReal Tc = temperature - _T_c2k;
1026  // The parial molar volume
1027  const ADReal V = 37.51 - 9.585e-2 * Tc + 8.74e-4 * Tc * Tc - 5.044e-7 * Tc * Tc * Tc;
1028 
1029  return 1.0e6 * _Mco2 / V;
1030 }
1031 
1032 Real
1034  Real pressure, Real temperature, Real Xnacl, Real saturation, unsigned int qp) const
1035 {
1036  // Check whether the input pressure and temperature are within the region of validity
1038 
1039  // As we do not require derivatives, we can simply ignore their initialisation
1040  const ADReal p = pressure;
1041  const ADReal T = temperature;
1042  const ADReal X = Xnacl;
1043 
1044  // FluidStateProperties data structure
1045  std::vector<FluidStateProperties> fsp(_num_phases, FluidStateProperties(_num_components));
1048 
1049  // Calculate equilibrium mass fractions in the two-phase state
1050  ADReal Xco2, Yh2o;
1051  equilibriumMassFractions(p, T, X, Xco2, Yh2o);
1052 
1053  // Save the mass fractions in the FluidStateMassFractions object
1054  const ADReal Yco2 = 1.0 - Yh2o;
1055  liquid.mass_fraction[_aqueous_fluid_component] = 1.0 - Xco2;
1056  liquid.mass_fraction[_gas_fluid_component] = Xco2;
1058  gas.mass_fraction[_gas_fluid_component] = Yco2;
1059 
1060  // Gas properties
1062 
1063  // Liquid properties
1064  const ADReal liquid_saturation = 1.0 - saturation;
1065  const ADReal liquid_pressure = p - _pc.capillaryPressure(liquid_saturation, qp);
1066  liquidProperties(liquid_pressure, T, X, fsp);
1067 
1068  // The total mass fraction of ncg (z) can now be calculated
1069  const ADReal Z = (saturation * gas.density * Yco2 + liquid_saturation * liquid.density * Xco2) /
1070  (saturation * gas.density + liquid_saturation * liquid.density);
1071 
1072  return Z.value();
1073 }
1074 
1075 ADReal
1077 {
1078  // Henry's constant for dissolution in water
1080 
1081  // The correction to salt is obtained through the salting out coefficient
1082  const std::vector<Real> b{1.19784e-1, -7.17823e-4, 4.93854e-6, -1.03826e-8, 1.08233e-11};
1083 
1084  // Need temperature in Celsius
1085  const ADReal Tc = temperature - _T_c2k;
1086 
1087  ADReal kb = 0.0;
1088  for (unsigned int i = 0; i < b.size(); ++i)
1089  kb += b[i] * std::pow(Tc, i);
1090 
1091  // Need salt mass fraction in molality
1092  const ADReal xmol = Xnacl / (1.0 - Xnacl) / _Mnacl;
1093 
1094  // Henry's constant and its derivative wrt temperature and salt mass fraction
1095  return Kh_h2o * std::pow(10.0, xmol * kb);
1096 }
1097 
1098 ADReal
1100 {
1101  // Henry's constant
1102  const ADReal Kh = henryConstant(temperature, Xnacl);
1103 
1104  ADReal hdis = -_R * temperature * temperature * Kh.derivatives()[_Tidx] / Kh / _Mco2;
1105 
1106  // Derivative of enthalpy of dissolution wrt temperature and xnacl requires the second
1107  // derivatives of Henry's constant. For simplicity, approximate these numerically
1108  const Real dT = temperature.value() * 1.0e-8;
1109  const ADReal T2 = temperature + dT;
1110  const ADReal Kh2 = henryConstant(T2, Xnacl);
1111 
1112  const Real dhdis_dT =
1113  (-_R * T2 * T2 * Kh2.derivatives()[_Tidx] / Kh2 / _Mco2 - hdis).value() / dT;
1114 
1115  const Real dX = Xnacl.value() * 1.0e-8;
1116  const ADReal X2 = Xnacl + dX;
1117  const ADReal Kh3 = henryConstant(temperature, X2);
1118 
1119  const Real dhdis_dX =
1120  (-_R * temperature * temperature * Kh3.derivatives()[_Tidx] / Kh3 / _Mco2 - hdis).value() /
1121  dX;
1122 
1123  hdis.derivatives() = temperature.derivatives() * dhdis_dT + Xnacl.derivatives() * dhdis_dX;
1124 
1125  return hdis;
1126 }
1127 
1128 ADReal
1130 {
1131  // Linear fit to model of Duan and Sun (2003) (in kJ/mol)
1132  const ADReal delta_h = -58.3533 + 0.134519 * temperature;
1133 
1134  // Convert to J/kg
1135  return delta_h * 1000.0 / _Mco2;
1136 }
1137 
1138 void
1140  Real temperature, Real f0, Real df0, Real f1, Real df1, Real & value, Real & deriv) const
1141 {
1142  // Coefficients of cubic polynomial
1143  const Real dT = _Tupper - _Tlower;
1144 
1145  const Real a = f0;
1146  const Real b = df0 * dT;
1147  const Real c = 3.0 * (f1 - f0) - (2.0 * df0 + df1) * dT;
1148  const Real d = 2.0 * (f0 - f1) + (df0 + df1) * dT;
1149 
1150  const Real t2 = temperature * temperature;
1151  const Real t3 = temperature * t2;
1152 
1153  value = a + b * temperature + c * t2 + d * t3;
1154  deriv = b + 2.0 * c * temperature + 3.0 * d * t2;
1155 }
1156 
1157 void
1159 {
1160  // The calculation of mass fractions is valid from 12C <= T <= 300C, and
1161  // pressure less than 60 MPa
1162  if (temperature < 285.15 || temperature > 573.15)
1163  mooseException(name() + ": temperature " + Moose::stringify(temperature) +
1164  " is outside range 285.15 K <= T <= 573.15 K");
1165 
1166  if (pressure > 6.0e7)
1167  mooseException(name() + ": pressure " + Moose::stringify(pressure) +
1168  " must be less than 60 MPa");
1169 }
void liquidProperties(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, std::vector< FluidStateProperties > &fsp) const
Thermophysical properties of the liquid state.
const unsigned int _aqueous_fluid_component
Fluid component number of the aqueous component.
const Real _Tlower
Temperature below which the Spycher, Pruess & Ennis-King (2003) model is used (K) ...
const unsigned int _aqueous_phase_number
Phase number of the aqueous phase.
unsigned int _num_phases
Number of phases.
const unsigned int _salt_component
Salt component index.
const Real _Mh2o
Molar mass of water (kg/mol)
ADReal activityCoefficient(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl) const
Activity coefficient for CO2 in brine.
const Real _T_c2k
Conversion from C to K.
FPADReal h_from_p_T_X(const FPADReal &pressure, const FPADReal &temperature, const FPADReal &xnacl) const
static InputParameters validParams()
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
virtual Real rho_from_p_T_X(Real pressure, Real temperature, Real xnacl) const override
ADReal activityCoefficientHighTemp(const ADReal &temperature, const ADReal &Xnacl) const
Activity coefficient for CO2 in brine used in the elevated temperature formulation.
const Real _Mnacl
Molar mass of NaCL.
ADReal equilibriumConstantH2O(const ADReal &temperature) const
Equilibrium constant for H2O For temperatures 12C <= T <= 99C, uses Spycher, Pruess and Ennis-King (2...
const Real _invMh2o
Inverse of molar mass of H2O (mol/kg)
const double tol
void fugacityCoefficientsHighTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &co2_density, const ADReal &xco2, const ADReal &yh2o, ADReal &fco2, ADReal &fh2o) const
Fugacity coefficients for H2O and CO2 at elevated temperatures (100C < T <= 300C).
unsigned int _gas_fluid_component
Fluid component number of the gas phase.
virtual std::string fluidName() const override
Fluid name.
const std::vector< Real > _co2_henry
Henry&#39;s coefficeients for CO2.
const unsigned int _Tidx
Index of derivative wrt temperature.
virtual Real totalMassFraction(Real pressure, Real temperature, Real Xnacl, Real saturation, unsigned int qp) const override
Total mass fraction of fluid component summed over all phases in the two-phase state for a specified ...
void funcABLowTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &co2_density, ADReal &A, ADReal &B) const
unsigned int _gas_phase_number
Phase number of the gas phase.
void clearFluidStateProperties(std::vector< FluidStateProperties > &fsp) const
Clears the contents of the FluidStateProperties data structure.
unsigned int _num_components
Number of components.
ADReal fugacityCoefficientH2OHighTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &co2_density, const ADReal &xco2, const ADReal &yh2o) const
const Real _Mco2
Molar mass of CO2 (kg/mol)
Compositional flash routines for miscible multiphase flow classes with multiple fluid components...
ADReal enthalpyOfDissolutionGas(const ADReal &temperature, const ADReal &Xnacl) const
Enthalpy of dissolution of gas phase CO2 in brine calculated using Henry&#39;s constant From Himmelblau...
const std::vector< double > y
static const std::string temperature
Definition: NS.h:59
virtual Real mu_from_p_T_X(Real pressure, Real temperature, Real xnacl) const override
Brine (NaCl in H2O) fluid properties as a function of pressure (Pa), temperature (K) and NaCl mass fr...
DualNumber< Real, DNDerivativeType, true > ADReal
void solveEquilibriumMoleFractionHighTemp(Real pressure, Real temperature, Real Xnacl, Real co2_density, Real &xco2, Real &yh2o) const
Function to solve for yh2o and xco2 iteratively in the elevated temperature regime (T > 100C) ...
virtual const std::string & name() const
void addRequiredParam(const std::string &name, const std::string &doc_string)
const BrineFluidProperties & _brine_fp
Fluid properties UserObject for water.
Real vaporPressure(Real temperature, Real xnacl) const
Brine vapour pressure From Haas, Physical properties of the coexisting phases and thermochemical prop...
virtual std::string fluidStateName() const override
Name of FluidState.
const SinglePhaseFluidProperties & _co2_fp
Fluid properties UserObject for the CO2.
const unsigned int _pidx
Index of derivative wrt pressure.
ADReal partialDensityCO2(const ADReal &temperature) const
Partial density of dissolved CO2 From Garcia, Density of aqueous solutions of CO2, LBNL-49023 (2001)
void funcABHighTemp(Real pressure, Real temperature, Real Xnacl, Real co2_density, Real xco2, Real yh2o, Real &A, Real &B) const
The function A (Eq.
virtual Real capillaryPressure(Real saturation, unsigned qp=0) const
Capillary pressure is calculated as a function of true saturation.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
Real deriv(unsigned n, unsigned alpha, unsigned beta, Real x)
const std::vector< double > x
void equilibriumMassFractions(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, ADReal &Xco2, ADReal &Yh2o) const
Mass fractions of CO2 in brine and water vapor in CO2 at equilibrium.
void thermophysicalProperties(Real pressure, Real temperature, Real Xnacl, Real Z, unsigned int qp, std::vector< FluidStateProperties > &fsp) const override
Determines the complete thermophysical state of the system for a given set of primary variables...
ADReal enthalpyOfDissolution(const ADReal &temperature) const
Enthalpy of dissolution of CO2 in brine calculated using linear fit to model of Duan and Sun...
void equilibriumMoleFractionsLowTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, ADReal &xco2, ADReal &yh2o) const
Mole fractions of CO2 in brine and water vapor in CO2 at equilibrium in the low temperature regime (T...
const Real _R
Universal gas constant (J/mol/K)
ADReal activityCoefficientH2O(const ADReal &temperature, const ADReal &xco2) const
Activity coefficient of H2O Eq.
AD data structure to pass calculated thermophysical properties.
Common class for single phase fluid properties.
FluidStatePhaseEnum
Phase state enum.
registerMooseObject("PorousFlowApp", PorousFlowBrineCO2)
void fugacityCoefficientsLowTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &co2_density, ADReal &fco2, ADReal &fh2o) const
Fugacity coefficients for H2O and CO2 for T <= 100C Eq.
void paramError(const std::string &param, Args... args) const
ADReal activityCoefficientCO2(const ADReal &temperature, const ADReal &xco2) const
Activity coefficient of CO2 Eq.
ADReal fugacityCoefficientCO2HighTemp(const ADReal &pressure, const ADReal &temperature, const ADReal &co2_density, const ADReal &xco2, const ADReal &yh2o) const
std::string stringify(const T &t)
virtual void rho_mu_from_p_T(Real p, Real T, Real &rho, Real &mu) const
Combined methods.
void equilibriumMoleFractions(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, ADReal &xco2, ADReal &yh2o) const
Mole fractions of CO2 in brine and water vapor in CO2 at equilibrium.
void twoPhaseProperties(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, const ADReal &Z, unsigned int qp, std::vector< FluidStateProperties > &fsp) const
Gas and liquid properties in the two-phase region.
static const std::string Z
Definition: NS.h:169
const Real _Rbar
Molar gas constant in bar cm^3 /(K mol)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const PorousFlowCapillaryPressure & _pc
Capillary pressure UserObject.
void phaseState(Real Zi, Real Xi, Real Yi, FluidStatePhaseEnum &phase_state) const
Determines the phase state gven the total mass fraction and equilibrium mass fractions.
Real vaporMassFraction(Real Z0, Real K0, Real K1) const
Solves Rachford-Rice equation to provide vapor mass fraction.
void massFractions(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, const ADReal &Z, FluidStatePhaseEnum &phase_state, std::vector< FluidStateProperties > &fsp) const
Mass fractions of CO2 and H2O in both phases, as well as derivatives wrt PorousFlow variables...
void gasProperties(const ADReal &pressure, const ADReal &temperature, std::vector< FluidStateProperties > &fsp) const
Thermophysical properties of the gaseous state.
Specialized class for brine and CO2 including calculation of mutual solubility of the two fluids usin...
static const std::string pressure
Definition: NS.h:56
const unsigned int _Zidx
Index of derivative wrt total mass fraction Z.
void mooseError(Args &&... args) const
e e e e s T T T T T rho v v T e p T T virtual T std::string fluidName() const
Fluid name.
void addClassDescription(const std::string &doc_string)
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
const unsigned int _Xidx
Index of derivative wrt salt mass fraction X.
void smoothCubicInterpolation(Real temperature, Real f0, Real df0, Real f1, Real df1, Real &value, Real &deriv) const
Cubic function to smoothly interpolate between the low temperature and elevated temperature models fo...
FluidStateProperties _empty_fsp
Empty FluidStateProperties object.
Real henryConstant(Real temperature, const std::vector< Real > &coeffs) const
IAPWS formulation of Henry&#39;s law constant for dissolution in water (implemented in water FluidPropert...
virtual void checkVariables(Real pressure, Real temperature) const
Check the input variables.
const Real _Zmin
Minimum Z - below this value all CO2 will be dissolved.
MooseUnits pow(const MooseUnits &, int)
PorousFlowBrineCO2(const InputParameters &parameters)
ADReal henryConstant(const ADReal &temperature, const ADReal &Xnacl) const
Henry&#39;s constant of dissolution of gas phase CO2 in brine.
ADReal equilibriumConstantCO2(const ADReal &temperature) const
Equilibrium constant for CO2 For temperatures 12C <= T <= 99C, uses Spycher, Pruess and Ennis-King (2...
void ErrorVector unsigned int
std::vector< ADReal > mass_fraction
ADReal saturation(const ADReal &pressure, const ADReal &temperature, const ADReal &Xnacl, const ADReal &Z, std::vector< FluidStateProperties > &fsp) const
Gas saturation in the two-phase region.
const Real _Tupper
Temperature above which the Spycher & Pruess (2010) model is used (K)