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