My Project
BrineCO2FluidSystem.hpp
Go to the documentation of this file.
1 // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=4 sw=4 sts=4:
3 /*
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 2 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 
19  Consult the COPYING file in the top-level source directory of this
20  module for the precise wording of the license and the list of
21  copyright holders.
22 */
28 #ifndef OPM_BRINE_CO2_SYSTEM_HPP
29 #define OPM_BRINE_CO2_SYSTEM_HPP
30 
31 #include "BaseFluidSystem.hpp"
32 #include "NullParameterCache.hpp"
33 
35 
41 
45 
46 #include <iostream>
47 
48 namespace Opm {
49 
57 template <class Scalar, class CO2Tables>
59  : public BaseFluidSystem<Scalar, BrineCO2FluidSystem<Scalar, CO2Tables> >
60 {
61  typedef ::Opm::H2O<Scalar> H2O_IAPWS;
62  typedef ::Opm::Brine<Scalar, H2O_IAPWS> Brine_IAPWS;
65 
66  typedef H2O_Tabulated H2O;
67 
68 public:
69  template <class Evaluation>
70  struct ParameterCache : public NullParameterCache<Evaluation>
71  {};
72 
76  typedef ::Opm::CO2<Scalar, CO2Tables> CO2;
77 
80 
81  /****************************************
82  * Fluid phase related static parameters
83  ****************************************/
84 
86  static const int numPhases = 2;
87 
89  static const int liquidPhaseIdx = 0;
91  static const int gasPhaseIdx = 1;
92 
96  static const char* phaseName(unsigned phaseIdx)
97  {
98  static const char* name[] = {
99  "liquid",
100  "gas"
101  };
102 
103  assert(phaseIdx < numPhases);
104  return name[phaseIdx];
105  }
106 
110  static bool isLiquid(unsigned phaseIdx)
111  {
112  assert(phaseIdx < numPhases);
113 
114  return phaseIdx != gasPhaseIdx;
115  }
116 
120  static bool isIdealGas(unsigned phaseIdx)
121  {
122  assert(phaseIdx < numPhases);
123 
124  if (phaseIdx == gasPhaseIdx)
125  return CO2::gasIsIdeal();
126  return false;
127  }
128 
132  static bool isIdealMixture([[maybe_unused]] unsigned phaseIdx)
133  {
134  assert(phaseIdx < numPhases);
135 
136  return true;
137  }
138 
142  static bool isCompressible([[maybe_unused]] unsigned phaseIdx)
143  {
144  assert(phaseIdx < numPhases);
145 
146  return true;
147  }
148 
149  /****************************************
150  * Component related static parameters
151  ****************************************/
153  static const int numComponents = 2;
154 
156  static const int BrineIdx = 0;
158  static const int CO2Idx = 1;
159 
163  static const char* componentName(unsigned compIdx)
164  {
165  static const char* name[] = {
166  Brine::name(),
167  CO2::name(),
168  };
169 
170  assert(compIdx < numComponents);
171  return name[compIdx];
172  }
173 
177  static Scalar molarMass(unsigned compIdx)
178  {
179  assert(compIdx < numComponents);
180  return (compIdx==BrineIdx)
181  ? Brine::molarMass()
182  : CO2::molarMass();
183  }
184 
185  /****************************************
186  * thermodynamic relations
187  ****************************************/
188 
192  static void init()
193  {
194  init(/*startTemp=*/273.15, /*endTemp=*/623.15, /*tempSteps=*/50,
195  /*startPressure=*/1e4, /*endPressure=*/40e6, /*pressureSteps=*/50);
196  }
197 
209  static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
210  Scalar pressMin, Scalar pressMax, unsigned nPress)
211  {
212  if (H2O::isTabulated) {
213  H2O_Tabulated::init(tempMin, tempMax, nTemp,
214  pressMin, pressMax, nPress);
215  }
216 
217  // set the salinity of brine to the one used by the CO2 tables
218  Brine_IAPWS::salinity = CO2Tables::brineSalinity;
219 
220  if (Brine::isTabulated) {
221  Brine_Tabulated::init(tempMin, tempMax, nTemp,
222  pressMin, pressMax, nPress);
223  }
224  }
225 
229  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
230  static LhsEval density(const FluidState& fluidState,
231  const ParameterCache<ParamCacheEval>& /*paramCache*/,
232  unsigned phaseIdx)
233  {
234  assert(phaseIdx < numPhases);
235 
236  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
237  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
238 
239  if (phaseIdx == liquidPhaseIdx) {
240  // use normalized composition for to calculate the density
241  // (the relations don't seem to take non-normalized
242  // compositions too well...)
243  LhsEval xlBrine = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, BrineIdx))));
244  LhsEval xlCO2 = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, CO2Idx))));
245  LhsEval sumx = xlBrine + xlCO2;
246  xlBrine /= sumx;
247  xlCO2 /= sumx;
248 
249  LhsEval result = liquidDensity_(temperature,
250  pressure,
251  xlBrine,
252  xlCO2);
253 
254  Valgrind::CheckDefined(result);
255  return result;
256  }
257 
258  assert(phaseIdx == gasPhaseIdx);
259 
260  // use normalized composition for to calculate the density
261  // (the relations don't seem to take non-normalized
262  // compositions too well...)
263  LhsEval xgBrine = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, BrineIdx))));
264  LhsEval xgCO2 = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, CO2Idx))));
265  LhsEval sumx = xgBrine + xgCO2;
266  xgBrine /= sumx;
267  xgCO2 /= sumx;
268 
269  LhsEval result = gasDensity_(temperature,
270  pressure,
271  xgBrine,
272  xgCO2);
273  Valgrind::CheckDefined(result);
274  return result;
275  }
276 
280  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
281  static LhsEval viscosity(const FluidState& fluidState,
282  const ParameterCache<ParamCacheEval>& /*paramCache*/,
283  unsigned phaseIdx)
284  {
285  assert(phaseIdx < numPhases);
286 
287  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
288  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
289 
290  if (phaseIdx == liquidPhaseIdx) {
291  // assume pure brine for the liquid phase. TODO: viscosity
292  // of mixture
293  LhsEval result = Brine::liquidViscosity(temperature, pressure);
294  Valgrind::CheckDefined(result);
295  return result;
296  }
297 
298  assert(phaseIdx == gasPhaseIdx);
299  LhsEval result = CO2::gasViscosity(temperature, pressure);
300  Valgrind::CheckDefined(result);
301  return result;
302  }
303 
307  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
308  static LhsEval fugacityCoefficient(const FluidState& fluidState,
309  const ParameterCache<ParamCacheEval>& /*paramCache*/,
310  unsigned phaseIdx,
311  unsigned compIdx)
312  {
313  assert(phaseIdx < numPhases);
314  assert(compIdx < numComponents);
315 
316  if (phaseIdx == gasPhaseIdx)
317  // use the fugacity coefficients of an ideal gas. the
318  // actual value of the fugacity is not relevant, as long
319  // as the relative fluid compositions are observed,
320  return 1.0;
321 
322  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
323  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
324  assert(temperature > 0);
325  assert(pressure > 0);
326 
327  // calulate the equilibrium composition for the given
328  // temperature and pressure. TODO: calculateMoleFractions()
329  // could use some cleanup.
330  LhsEval xlH2O, xgH2O;
331  LhsEval xlCO2, xgCO2;
333  pressure,
335  /*knownPhaseIdx=*/-1,
336  xlCO2,
337  xgH2O);
338 
339  // normalize the phase compositions
340  xlCO2 = max(0.0, min(1.0, xlCO2));
341  xgH2O = max(0.0, min(1.0, xgH2O));
342 
343  xlH2O = 1.0 - xlCO2;
344  xgCO2 = 1.0 - xgH2O;
345 
346  if (compIdx == BrineIdx) {
347  Scalar phigH2O = 1.0;
348  return phigH2O * xgH2O / xlH2O;
349  }
350  else {
351  assert(compIdx == CO2Idx);
352 
353  Scalar phigCO2 = 1.0;
354  return phigCO2 * xgCO2 / xlCO2;
355  };
356  }
357 
361  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
362  static LhsEval diffusionCoefficient(const FluidState& fluidState,
363  const ParameterCache<ParamCacheEval>& /*paramCache*/,
364  unsigned phaseIdx,
365  unsigned /*compIdx*/)
366  {
367  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
368  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
369  if (phaseIdx == liquidPhaseIdx)
370  return BinaryCoeffBrineCO2::liquidDiffCoeff(temperature, pressure);
371 
372  assert(phaseIdx == gasPhaseIdx);
373  return BinaryCoeffBrineCO2::gasDiffCoeff(temperature, pressure);
374  }
375 
379  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
380  static LhsEval enthalpy(const FluidState& fluidState,
381  const ParameterCache<ParamCacheEval>& /*paramCache*/,
382  unsigned phaseIdx)
383  {
384  assert(phaseIdx < numPhases);
385 
386  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
387  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
388 
389  if (phaseIdx == liquidPhaseIdx) {
390  const LhsEval& XlCO2 = decay<LhsEval>(fluidState.massFraction(phaseIdx, CO2Idx));
391  const LhsEval& result = liquidEnthalpyBrineCO2_(temperature,
392  pressure,
394  XlCO2);
395  Valgrind::CheckDefined(result);
396  return result;
397  }
398  else {
399  const LhsEval& XCO2 = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, CO2Idx));
400  const LhsEval& XBrine = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, BrineIdx));
401 
402  LhsEval result = 0;
403  result += XBrine * Brine::gasEnthalpy(temperature, pressure);
404  result += XCO2 * CO2::gasEnthalpy(temperature, pressure);
405  Valgrind::CheckDefined(result);
406  return result;
407  }
408  }
409 
413  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
414  static LhsEval thermalConductivity(const FluidState& /*fluidState*/,
415  const ParameterCache<ParamCacheEval>& /*paramCache*/,
416  unsigned phaseIdx)
417  {
418  // TODO way too simple!
419  if (phaseIdx == liquidPhaseIdx)
420  return 0.6; // conductivity of water[W / (m K ) ]
421 
422  // gas phase
423  return 0.025; // conductivity of air [W / (m K ) ]
424  }
425 
438  template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
439  static LhsEval heatCapacity(const FluidState& fluidState,
440  const ParameterCache<ParamCacheEval>& /*paramCache*/,
441  unsigned phaseIdx)
442  {
443  assert(phaseIdx < numPhases);
444 
445  const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
446  const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
447 
448  if(phaseIdx == liquidPhaseIdx)
449  return H2O::liquidHeatCapacity(temperature, pressure);
450  else
451  return CO2::gasHeatCapacity(temperature, pressure);
452  }
453 
454 private:
455  template <class LhsEval>
456  static LhsEval gasDensity_(const LhsEval& T,
457  const LhsEval& pg,
458  const LhsEval& xgH2O,
459  const LhsEval& xgCO2)
460  {
461  Valgrind::CheckDefined(T);
462  Valgrind::CheckDefined(pg);
463  Valgrind::CheckDefined(xgH2O);
464  Valgrind::CheckDefined(xgCO2);
465 
466  return CO2::gasDensity(T, pg);
467  }
468 
469  /***********************************************************************/
470  /* */
471  /* Total brine density with dissolved CO2 */
472  /* rho_{b,CO2} = rho_w + contribution(salt) + contribution(CO2) */
473  /* */
474  /***********************************************************************/
475  template <class LhsEval>
476  static LhsEval liquidDensity_(const LhsEval& T,
477  const LhsEval& pl,
478  const LhsEval& xlH2O,
479  const LhsEval& xlCO2)
480  {
481  Valgrind::CheckDefined(T);
482  Valgrind::CheckDefined(pl);
483  Valgrind::CheckDefined(xlH2O);
484  Valgrind::CheckDefined(xlCO2);
485 
486  if(T < 273.15) {
487  std::ostringstream oss;
488  oss << "Liquid density for Brine and CO2 is only "
489  "defined above 273.15K (is "<<T<<"K)";
490  throw NumericalIssue(oss.str());
491  }
492  if(pl >= 2.5e8) {
493  std::ostringstream oss;
494  oss << "Liquid density for Brine and CO2 is only "
495  "defined below 250MPa (is "<<pl<<"Pa)";
496  throw NumericalIssue(oss.str());
497  }
498 
499  const LhsEval& rho_brine = Brine::liquidDensity(T, pl);
500  const LhsEval& rho_pure = H2O::liquidDensity(T, pl);
501  const LhsEval& rho_lCO2 = liquidDensityWaterCO2_(T, pl, xlH2O, xlCO2);
502  const LhsEval& contribCO2 = rho_lCO2 - rho_pure;
503 
504  return rho_brine + contribCO2;
505  }
506 
507  template <class LhsEval>
508  static LhsEval liquidDensityWaterCO2_(const LhsEval& temperature,
509  const LhsEval& pl,
510  const LhsEval& /*xlH2O*/,
511  const LhsEval& xlCO2)
512  {
513  Scalar M_CO2 = CO2::molarMass();
514  Scalar M_H2O = H2O::molarMass();
515 
516  const LhsEval& tempC = temperature - 273.15; /* tempC : temperature in °C */
517  const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl);
518  // calculate the mole fraction of CO2 in the liquid. note that xlH2O is available
519  // as a function parameter, but in the case of a pure gas phase the value of M_T
520  // for the virtual liquid phase can become very large
521  const LhsEval xlH2O = 1.0 - xlCO2;
522  const LhsEval& M_T = M_H2O * xlH2O + M_CO2 * xlCO2;
523  const LhsEval& V_phi =
524  (37.51 +
525  tempC*(-9.585e-2 +
526  tempC*(8.74e-4 -
527  tempC*5.044e-7))) / 1.0e6;
528  return 1/ (xlCO2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
529  }
530 
531  template <class LhsEval>
532  static LhsEval liquidEnthalpyBrineCO2_(const LhsEval& T,
533  const LhsEval& p,
534  Scalar S, // salinity
535  const LhsEval& X_CO2_w)
536  {
537  /* X_CO2_w : mass fraction of CO2 in brine */
538 
539  /* same function as enthalpy_brine, only extended by CO2 content */
540 
541  /*Numerical coefficents from PALLISER*/
542  static Scalar f[] = {
543  2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
544  };
545 
546  /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
547  static Scalar a[4][3] = {
548  { 9633.6, -4080.0, +286.49 },
549  { +166.58, +68.577, -4.6856 },
550  { -0.90963, -0.36524, +0.249667E-1 },
551  { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
552  };
553 
554  LhsEval theta, h_NaCl;
555  LhsEval h_ls1, d_h;
556  LhsEval delta_h;
557  LhsEval delta_hCO2, hg, hw;
558 
559  theta = T - 273.15;
560 
561  // Regularization
562  Scalar scalarTheta = scalarValue(theta);
563  Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
564  if (S > S_lSAT)
565  S = S_lSAT;
566 
567  hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
568 
569  /*DAUBERT and DANNER*/
570  /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
571  +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
572 
573  Scalar m = 1E3/58.44 * S/(1-S);
574  int i = 0;
575  int j = 0;
576  d_h = 0;
577 
578  for (i = 0; i<=3; i++) {
579  for (j=0; j<=2; j++) {
580  d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * std::pow(m, j);
581  }
582  }
583  /* heat of dissolution for halite according to Michaelides 1971 */
584  delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
585 
586  /* Enthalpy of brine without CO2 */
587  h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
588 
589  /* heat of dissolution for CO2 according to Fig. 6 in Duan and Sun 2003. (kJ/kg)
590  In the relevant temperature ranges CO2 dissolution is
591  exothermal */
592  delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
593 
594  /* enthalpy contribution of CO2 (kJ/kg) */
595  hg = CO2::gasEnthalpy(T, p)/1E3 + delta_hCO2;
596 
597  /* Enthalpy of brine with dissolved CO2 */
598  return (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1E3; /*J/kg*/
599  }
600 };
601 
602 } // namespace Opm
603 
604 #endif
The base class for all fluid systems.
A class for the brine fluid properties.
Binary coefficients for brine and CO2.
A class for the CO2 fluid properties.
Binary coefficients for water and CO2.
Binary coefficients for water and nitrogen.
Relations valid for an ideal gas.
A parameter cache which does nothing.
A simplistic class representing the fluid properties.
A simple version of pure water.
A generic class which tabulates all thermodynamic properties of a given component.
The base class for all fluid systems.
Definition: BaseFluidSystem.hpp:44
Scalar Scalar
The type used for scalar quantities.
Definition: BaseFluidSystem.hpp:49
Binary coefficients for brine and CO2.
Definition: Brine_CO2.hpp:42
static void calculateMoleFractions(const Evaluation &temperature, const Evaluation &pg, Scalar salinity, const int knownPhaseIdx, Evaluation &xlCO2, Evaluation &ygH2O, bool extrapolate=false)
Returns the mol (!) fraction of CO2 in the liquid phase and the mol_ (!) fraction of H2O in the gas p...
Definition: Brine_CO2.hpp:97
static Evaluation gasDiffCoeff(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
Binary diffusion coefficent [m^2/s] of water in the CO2 phase.
Definition: Brine_CO2.hpp:56
static Evaluation liquidDiffCoeff(const Evaluation &, const Evaluation &)
Binary diffusion coefficent [m^2/s] of CO2 in the brine phase.
Definition: Brine_CO2.hpp:73
A two-phase fluid system with water and CO2.
Definition: BrineCO2FluidSystem.hpp:60
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: BrineCO2FluidSystem.hpp:120
static LhsEval enthalpy(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy [J/kg...
Definition: BrineCO2FluidSystem.hpp:380
static const char * phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition: BrineCO2FluidSystem.hpp:96
static LhsEval thermalConductivity(const FluidState &, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition: BrineCO2FluidSystem.hpp:414
static const int gasPhaseIdx
The index of the gas phase.
Definition: BrineCO2FluidSystem.hpp:91
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition: BrineCO2FluidSystem.hpp:177
static const char * componentName(unsigned compIdx)
Return the human readable name of a component.
Definition: BrineCO2FluidSystem.hpp:163
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition: BrineCO2FluidSystem.hpp:230
static const int numComponents
Number of chemical species in the fluid system.
Definition: BrineCO2FluidSystem.hpp:153
static LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition: BrineCO2FluidSystem.hpp:308
::Opm::CO2< Scalar, CO2Tables > CO2
The type of the component for pure CO2 used by the fluid system.
Definition: BrineCO2FluidSystem.hpp:76
static const int CO2Idx
The index of the CO2 component.
Definition: BrineCO2FluidSystem.hpp:158
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition: BrineCO2FluidSystem.hpp:281
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition: BrineCO2FluidSystem.hpp:110
BinaryCoeff::Brine_CO2< Scalar, H2O, CO2 > BinaryCoeffBrineCO2
The binary coefficients for brine and CO2 used by this fluid system.
Definition: BrineCO2FluidSystem.hpp:79
static bool isCompressible([[maybe_unused]] unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: BrineCO2FluidSystem.hpp:142
static const int numPhases
The number of phases considered by the fluid system.
Definition: BrineCO2FluidSystem.hpp:86
static const int BrineIdx
The index of the brine component.
Definition: BrineCO2FluidSystem.hpp:156
static bool isIdealMixture([[maybe_unused]] unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: BrineCO2FluidSystem.hpp:132
static void init()
Initialize the fluid system's static parameters.
Definition: BrineCO2FluidSystem.hpp:192
static const int liquidPhaseIdx
The index of the liquid phase.
Definition: BrineCO2FluidSystem.hpp:89
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the fluid system's static parameters using problem specific temperature and pressure range...
Definition: BrineCO2FluidSystem.hpp:209
Brine_Tabulated Brine
The type of the component for brine used by the fluid system.
Definition: BrineCO2FluidSystem.hpp:74
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition: BrineCO2FluidSystem.hpp:439
static LhsEval diffusionCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx, unsigned)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition: BrineCO2FluidSystem.hpp:362
A class for the brine fluid properties.
Definition: Brine.hpp:46
static Scalar salinity
The mass fraction of salt assumed to be in the brine.
Definition: Brine.hpp:49
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the component [J/kg] as a liquid.
Definition: CO2.hpp:255
static const char * name()
A human readable name for the CO2.
Definition: CO2.hpp:60
static Scalar molarMass()
The mass in [kg] of one mole of CO2.
Definition: CO2.hpp:66
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
Specific enthalpy of gaseous CO2 [J/kg].
Definition: CO2.hpp:164
static Evaluation gasViscosity(Evaluation temperature, const Evaluation &pressure, bool extrapolate=false)
The dynamic viscosity [Pa s] of CO2.
Definition: CO2.hpp:203
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition: CO2.hpp:157
static Evaluation gasDensity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
The density of CO2 at a given pressure and temperature [kg/m^3].
Definition: CO2.hpp:189
Material properties of pure water .
Definition: H2O.hpp:62
A parameter cache which does nothing.
Definition: NullParameterCache.hpp:40
A generic class which tabulates all thermodynamic properties of a given component.
Definition: TabulatedComponent.hpp:56
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the gas .
Definition: TabulatedComponent.hpp:282
static Evaluation liquidHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the liquid .
Definition: TabulatedComponent.hpp:333
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the tables.
Definition: TabulatedComponent.hpp:72
static Scalar molarMass()
The molar mass in of the component.
Definition: TabulatedComponent.hpp:221
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure)
The dynamic viscosity of liquid.
Definition: TabulatedComponent.hpp:478
static const char * name()
A human readable name for the component.
Definition: TabulatedComponent.hpp:215
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure)
The density of liquid at a given pressure and temperature .
Definition: TabulatedComponent.hpp:444
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &pressure)
Specific enthalpy of the liquid .
Definition: TabulatedComponent.hpp:299
Definition: BrineCO2FluidSystem.hpp:71