My Project
BrineCo2Pvt.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 */
27 #ifndef OPM_BRINE_CO2_PVT_HPP
28 #define OPM_BRINE_CO2_PVT_HPP
29 
31 
40 #include <opm/material/components/co2tables.inc>
41 
42 
43 #if HAVE_ECL_INPUT
44 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
45 #include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
46 #include <opm/parser/eclipse/EclipseState/Tables/TableManager.hpp>
47 #endif
48 
49 #include <vector>
50 
51 namespace Opm {
56 template <class Scalar>
58 {
59  typedef std::vector<std::pair<Scalar, Scalar> > SamplingPoints;
60  static const bool extrapolate = true;
61  //typedef H2O<Scalar> H2O_IAPWS;
62  //typedef Brine<Scalar, H2O_IAPWS> Brine_IAPWS;
63  //typedef TabulatedComponent<Scalar, H2O_IAPWS> H2O_Tabulated;
64  //typedef TabulatedComponent<Scalar, Brine_IAPWS> Brine_Tabulated;
65 
66  //typedef H2O_Tabulated H2O;
67  //typedef Brine_Tabulated Brine;
68 
69 
70 public:
72  typedef ::Opm::Brine<Scalar, H2O> Brine;
73  typedef ::Opm::CO2<Scalar, CO2Tables> CO2;
74 
76 
79 
80  explicit BrineCo2Pvt() = default;
81  BrineCo2Pvt(const std::vector<Scalar>& brineReferenceDensity,
82  const std::vector<Scalar>& co2ReferenceDensity,
83  const std::vector<Scalar>& salinity)
84  : brineReferenceDensity_(brineReferenceDensity),
85  co2ReferenceDensity_(co2ReferenceDensity),
86  salinity_(salinity)
87  {
88  }
89 #if HAVE_ECL_INPUT
94  void initFromState(const EclipseState& eclState, const Schedule&)
95  {
96  if( !eclState.getTableManager().getDensityTable().empty()) {
97  std::cerr << "WARNING: CO2STOR is enabled but DENSITY is in the deck. \n" <<
98  "The surface density is computed based on CO2-BRINE PVT at standard conditions (STCOND) and DENSITY is ignored " << std::endl;
99  }
100 
101  if( eclState.getTableManager().hasTables("PVDO") || !eclState.getTableManager().getPvtoTables().empty()) {
102  std::cerr << "WARNING: CO2STOR is enabled but PVDO or PVTO is in the deck. \n" <<
103  "BRINE PVT properties are computed based on the Hu et al. pvt model and PVDO/PVTO input is ignored. " << std::endl;
104  }
105 
106  // We only supported single pvt region for the co2-brine module
107  size_t numRegions = 1;
108  setNumRegions(numRegions);
109  size_t regionIdx = 0;
110  // Currently we only support constant salinity
111  const Scalar molality = eclState.getTableManager().salinity(); // mol/kg
112  const Scalar MmNaCl = 58e-3; // molar mass of NaCl [kg/mol]
113  // convert to mass fraction
114  Brine::salinity = 1 / ( 1 + 1 / (molality*MmNaCl)); //
115  salinity_[regionIdx] = Brine::salinity;
116  // set the surface conditions using the STCOND keyword
117  Scalar T_ref = eclState.getTableManager().stCond().temperature;
118  Scalar P_ref = eclState.getTableManager().stCond().pressure;
119 
120  brineReferenceDensity_[regionIdx] = Brine::liquidDensity(T_ref, P_ref, extrapolate);
121  co2ReferenceDensity_[regionIdx] = CO2::gasDensity(T_ref, P_ref, extrapolate);
122  }
123 #endif
124 
125  void setNumRegions(size_t numRegions)
126  {
127  brineReferenceDensity_.resize(numRegions);
128  co2ReferenceDensity_.resize(numRegions);
129  salinity_.resize(numRegions);
130  }
131 
132 
136  void setReferenceDensities(unsigned regionIdx,
137  Scalar rhoRefBrine,
138  Scalar rhoRefCO2,
139  Scalar /*rhoRefWater*/)
140  {
141  brineReferenceDensity_[regionIdx] = rhoRefBrine;
142  co2ReferenceDensity_[regionIdx] = rhoRefCO2;
143  }
144 
145 
149  void initEnd()
150  {
151 
152  }
153 
157  unsigned numRegions() const
158  { return brineReferenceDensity_.size(); }
159 
163  template <class Evaluation>
164  Evaluation internalEnergy(unsigned regionIdx,
165  const Evaluation& temperature,
166  const Evaluation& pressure,
167  const Evaluation& Rs) const
168  {
169 
170  const Evaluation xlCO2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx));
171  return (liquidEnthalpyBrineCO2_(temperature,
172  pressure,
173  salinity_[regionIdx],
174  xlCO2)
175  - pressure / density_(regionIdx, temperature, pressure, Rs));
176  }
177 
181  template <class Evaluation>
182  Evaluation viscosity(unsigned regionIdx,
183  const Evaluation& temperature,
184  const Evaluation& pressure,
185  const Evaluation& /*Rs*/) const
186  {
187  //TODO: The viscosity does not yet depend on the composition
188  return saturatedViscosity(regionIdx, temperature, pressure);
189  }
190 
194  template <class Evaluation>
195  Evaluation saturatedViscosity(unsigned /*regionIdx*/,
196  const Evaluation& temperature,
197  const Evaluation& pressure) const
198  {
199  return Brine::liquidViscosity(temperature, pressure);
200  }
201 
205  template <class Evaluation>
206  Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
207  const Evaluation& temperature,
208  const Evaluation& pressure,
209  const Evaluation& Rs) const
210  {
211  return density_(regionIdx, temperature, pressure, Rs)/brineReferenceDensity_[regionIdx];
212  }
213 
217  template <class Evaluation>
218  Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
219  const Evaluation& temperature,
220  const Evaluation& pressure) const
221  {
222  Evaluation rsSat = rsSat_(regionIdx, temperature, pressure);
223  return density_(regionIdx, temperature, pressure, rsSat)/brineReferenceDensity_[regionIdx];
224  }
225 
232  template <class Evaluation>
233  Evaluation saturationPressure(unsigned /*regionIdx*/,
234  const Evaluation& /*temperature*/,
235  const Evaluation& /*Rs*/) const
236  {
237  throw std::runtime_error("Requested the saturation pressure for the brine-co2 pvt module. Not yet implemented.");
238  }
239 
243  template <class Evaluation>
244  Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
245  const Evaluation& temperature,
246  const Evaluation& pressure,
247  const Evaluation& /*oilSaturation*/,
248  const Evaluation& /*maxOilSaturation*/) const
249  {
250  //TODO support VAPPARS
251  return rsSat_(regionIdx, temperature, pressure);
252  }
253 
257  template <class Evaluation>
258  Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
259  const Evaluation& temperature,
260  const Evaluation& pressure) const
261  {
262  return rsSat_(regionIdx, temperature, pressure);
263  }
264 
265  const Scalar oilReferenceDensity(unsigned regionIdx) const
266  { return brineReferenceDensity_[regionIdx]; }
267 
268  const Scalar gasReferenceDensity(unsigned regionIdx) const
269  { return co2ReferenceDensity_[regionIdx]; }
270 
271  const Scalar salinity(unsigned regionIdx) const
272  { return salinity_[regionIdx]; }
273 
274  bool operator==(const BrineCo2Pvt<Scalar>& data) const
275  {
276  return co2ReferenceDensity_ == data.co2ReferenceDensity_ &&
277  brineReferenceDensity_ == data.brineReferenceDensity_;
278  }
279 
280  template <class Evaluation>
281  Evaluation diffusionCoefficient(const Evaluation& temperature,
282  const Evaluation& pressure,
283  unsigned /*compIdx*/) const
284  {
285  //Diffusion coefficient of CO2 in pure water according to (McLachlan and Danckwerts, 1972)
286  const Evaluation log_D_H20 = -4.1764 + 712.52 / temperature - 2.5907e5 / (temperature*temperature);
287 
288  //Diffusion coefficient of CO2 in the brine phase modified following (Ratcliff and Holdcroft,1963 and Al-Rawajfeh, 2004)
289  const Evaluation& mu_H20 = H2O::liquidViscosity(temperature, pressure, extrapolate); // Water viscosity
290  const Evaluation& mu_Brine = Brine::liquidViscosity(temperature, pressure); // Brine viscosity
291  const Evaluation log_D_Brine = log_D_H20 - 0.87*log10(mu_Brine / mu_H20);
292 
293  return pow(Evaluation(10), log_D_Brine) * 1e-4; // convert from cm2/s to m2/s
294  }
295 
296 private:
297  std::vector<Scalar> brineReferenceDensity_;
298  std::vector<Scalar> co2ReferenceDensity_;
299  std::vector<Scalar> salinity_;
300 
301  template <class LhsEval>
302  LhsEval density_(unsigned regionIdx,
303  const LhsEval& temperature,
304  const LhsEval& pressure,
305  const LhsEval& Rs) const
306  {
307  LhsEval xlCO2 = convertXoGToxoG_(convertRsToXoG_(Rs,regionIdx));
308  LhsEval result = liquidDensity_(temperature,
309  pressure,
310  xlCO2);
311 
312  Valgrind::CheckDefined(result);
313  return result;
314  }
315 
316 
317  template <class LhsEval>
318  LhsEval liquidDensity_(const LhsEval& T,
319  const LhsEval& pl,
320  const LhsEval& xlCO2) const
321  {
322  Valgrind::CheckDefined(T);
323  Valgrind::CheckDefined(pl);
324  Valgrind::CheckDefined(xlCO2);
325 
326  if(!extrapolate && T < 273.15) {
327  std::ostringstream oss;
328  oss << "Liquid density for Brine and CO2 is only "
329  "defined above 273.15K (is "<<T<<"K)";
330  throw NumericalIssue(oss.str());
331  }
332  if(!extrapolate && pl >= 2.5e8) {
333  std::ostringstream oss;
334  oss << "Liquid density for Brine and CO2 is only "
335  "defined below 250MPa (is "<<pl<<"Pa)";
336  throw NumericalIssue(oss.str());
337  }
338 
339  const LhsEval& rho_brine = Brine::liquidDensity(T, pl, extrapolate);
340  const LhsEval& rho_pure = H2O::liquidDensity(T, pl, extrapolate);
341  const LhsEval& rho_lCO2 = liquidDensityWaterCO2_(T, pl, xlCO2);
342  const LhsEval& contribCO2 = rho_lCO2 - rho_pure;
343 
344  return rho_brine + contribCO2;
345  }
346 
347  template <class LhsEval>
348  LhsEval liquidDensityWaterCO2_(const LhsEval& temperature,
349  const LhsEval& pl,
350  const LhsEval& xlCO2) const
351  {
352  Scalar M_CO2 = CO2::molarMass();
353  Scalar M_H2O = H2O::molarMass();
354 
355  const LhsEval& tempC = temperature - 273.15; /* tempC : temperature in °C */
356  const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl, extrapolate);
357  // calculate the mole fraction of CO2 in the liquid. note that xlH2O is available
358  // as a function parameter, but in the case of a pure gas phase the value of M_T
359  // for the virtual liquid phase can become very large
360  const LhsEval xlH2O = 1.0 - xlCO2;
361  const LhsEval& M_T = M_H2O * xlH2O + M_CO2 * xlCO2;
362  const LhsEval& V_phi =
363  (37.51 +
364  tempC*(-9.585e-2 +
365  tempC*(8.74e-4 -
366  tempC*5.044e-7))) / 1.0e6;
367  return 1/ (xlCO2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
368  }
369 
374  template <class LhsEval>
375  LhsEval convertRsToXoG_(const LhsEval& Rs, unsigned regionIdx) const
376  {
377  Scalar rho_oRef = brineReferenceDensity_[regionIdx];
378  Scalar rho_gRef = co2ReferenceDensity_[regionIdx];
379 
380  const LhsEval& rho_oG = Rs*rho_gRef;
381  return rho_oG/(rho_oRef + rho_oG);
382  }
383 
384 
388  template <class LhsEval>
389  LhsEval convertXoGToxoG_(const LhsEval& XoG) const
390  {
391  Scalar M_CO2 = CO2::molarMass();
392  Scalar M_Brine = Brine::molarMass();
393  return XoG*M_Brine / (M_CO2*(1 - XoG) + XoG*M_Brine);
394  }
395 
396 
400  template <class LhsEval>
401  LhsEval convertxoGToXoG(const LhsEval& xoG) const
402  {
403  Scalar M_CO2 = CO2::molarMass();
404  Scalar M_Brine = Brine::molarMass();
405 
406  return xoG*M_CO2 / (xoG*(M_CO2 - M_Brine) + M_Brine);
407  }
408 
409 
414  template <class LhsEval>
415  LhsEval convertXoGToRs(const LhsEval& XoG, unsigned regionIdx) const
416  {
417  Scalar rho_oRef = brineReferenceDensity_[regionIdx];
418  Scalar rho_gRef = co2ReferenceDensity_[regionIdx];
419 
420  return XoG/(1.0 - XoG)*(rho_oRef/rho_gRef);
421  }
422 
423 
424  template <class LhsEval>
425  LhsEval rsSat_(unsigned regionIdx,
426  const LhsEval& temperature,
427  const LhsEval& pressure) const
428  {
429  // calulate the equilibrium composition for the given
430  // temperature and pressure.
431  LhsEval xgH2O;
432  LhsEval xlCO2;
434  pressure,
435  salinity_[regionIdx],
436  /*knownPhaseIdx=*/-1,
437  xlCO2,
438  xgH2O);
439 
440  // normalize the phase compositions
441  xlCO2 = max(0.0, min(1.0, xlCO2));
442 
443  return convertXoGToRs(convertxoGToXoG(xlCO2), regionIdx);
444  }
445 
446  template <class LhsEval>
447  static LhsEval liquidEnthalpyBrineCO2_(const LhsEval& T,
448  const LhsEval& p,
449  Scalar S, // salinity
450  const LhsEval& X_CO2_w)
451  {
452  /* X_CO2_w : mass fraction of CO2 in brine */
453 
454  /* same function as enthalpy_brine, only extended by CO2 content */
455 
456  /*Numerical coefficents from PALLISER*/
457  static Scalar f[] = {
458  2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
459  };
460 
461  /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
462  static Scalar a[4][3] = {
463  { 9633.6, -4080.0, +286.49 },
464  { +166.58, +68.577, -4.6856 },
465  { -0.90963, -0.36524, +0.249667E-1 },
466  { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
467  };
468 
469  LhsEval theta, h_NaCl;
470  LhsEval h_ls1, d_h;
471  LhsEval delta_h;
472  LhsEval delta_hCO2, hg, hw;
473 
474  theta = T - 273.15;
475 
476  // Regularization
477  Scalar scalarTheta = scalarValue(theta);
478  Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
479  if (S > S_lSAT)
480  S = S_lSAT;
481 
482  hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
483 
484  /*DAUBERT and DANNER*/
485  /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
486  +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
487 
488  Scalar m = 1E3/58.44 * S/(1-S);
489  int i = 0;
490  int j = 0;
491  d_h = 0;
492 
493  for (i = 0; i<=3; i++) {
494  for (j=0; j<=2; j++) {
495  d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * std::pow(m, j);
496  }
497  }
498  /* heat of dissolution for halite according to Michaelides 1971 */
499  delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
500 
501  /* Enthalpy of brine without CO2 */
502  h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
503 
504  /* heat of dissolution for CO2 according to Fig. 6 in Duan and Sun 2003. (kJ/kg)
505  In the relevant temperature ranges CO2 dissolution is
506  exothermal */
507  delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
508 
509  /* enthalpy contribution of CO2 (kJ/kg) */
510  hg = CO2::gasEnthalpy(T, p, extrapolate)/1E3 + delta_hCO2;
511 
512  /* Enthalpy of brine with dissolved CO2 */
513  return (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1E3; /*J/kg*/
514  }
515 
516 };
517 
518 } // namespace Opm
519 
520 #endif
A class for the brine fluid properties.
Binary coefficients for brine and CO2.
A class for the CO2 fluid properties.
A central place for various physical constants occuring in some equations.
Binary coefficients for water and CO2.
A simple version of pure water with density from Hu et al.
Implements a linearly interpolated scalar function that depends on one variable.
A generic class which tabulates all thermodynamic properties of a given component.
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Binary coefficients for brine and CO2.
Definition: Brine_CO2.hpp:41
static void calculateMoleFractions(const Evaluation &temperature, const Evaluation &pg, Scalar salinity, const int knownPhaseIdx, Evaluation &xlCO2, Evaluation &ygH2O)
Returns the mol (!) fraction of CO2 in the liquid phase and the mol_ (!) fraction of H2O in the gas p...
Definition: Brine_CO2.hpp:96
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
Definition: BrineCo2Pvt.hpp:58
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition: BrineCo2Pvt.hpp:164
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &, const Evaluation &) const
Returns the gas dissoluiton factor [m^3/m^3] of the liquid phase.
Definition: BrineCo2Pvt.hpp:244
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the formation volume factor [-] of the fluid phase.
Definition: BrineCo2Pvt.hpp:206
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: BrineCo2Pvt.hpp:157
Evaluation saturationPressure(unsigned, const Evaluation &, const Evaluation &) const
Returns the saturation pressure of the brine phase [Pa] depending on its mass fraction of the gas com...
Definition: BrineCo2Pvt.hpp:233
void initEnd()
Finish initializing the oil phase PVT properties.
Definition: BrineCo2Pvt.hpp:149
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of brine saturated with CO2 at a given pressure.
Definition: BrineCo2Pvt.hpp:218
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefBrine, Scalar rhoRefCO2, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition: BrineCo2Pvt.hpp:136
BinaryCoeff::Brine_CO2< Scalar, H2O, CO2 > BinaryCoeffBrineCO2
The binary coefficients for brine and CO2 used by this fluid system.
Definition: BrineCo2Pvt.hpp:78
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition: BrineCo2Pvt.hpp:182
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns thegas dissoluiton factor [m^3/m^3] of the liquid phase.
Definition: BrineCo2Pvt.hpp:258
Evaluation saturatedViscosity(unsigned, const Evaluation &temperature, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of oil saturated gas at given pressure.
Definition: BrineCo2Pvt.hpp:195
A class for the brine fluid properties.
Definition: Brine.hpp:46
static Scalar molarMass()
The molar mass in of the component.
Definition: Brine.hpp:80
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
The density of the liquid component at a given pressure in and temperature in .
Definition: Brine.hpp:250
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &)
The dynamic viscosity of pure water.
Definition: Brine.hpp:327
static Scalar salinity
The mass fraction of salt assumed to be in the brine.
Definition: Brine.hpp:49
A class for the CO2 fluid properties.
Definition: CO2.hpp:53
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 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
A simple version of pure water with density from Hu et al.
Definition: SimpleHuDuanH2O.hpp:70
static Evaluation liquidEnthalpy(const Evaluation &temperature, const Evaluation &)
Specific enthalpy of liquid water .
Definition: SimpleHuDuanH2O.hpp:198
static Evaluation liquidViscosity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate)
The dynamic viscosity of pure water.
Definition: SimpleHuDuanH2O.hpp:350
static Evaluation liquidDensity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate)
The density of pure water at a given pressure and temperature .
Definition: SimpleHuDuanH2O.hpp:309
static Scalar molarMass()
The molar mass in of water.
Definition: SimpleHuDuanH2O.hpp:104
Implements a linearly interpolated scalar function that depends on one variable.
Definition: Tabulated1DFunction.hpp:47