3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fluidsystems/brine.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_BRINE_FLUID_SYSTEM_HH
25#define DUMUX_BRINE_FLUID_SYSTEM_HH
26
27#include <dune/common/math.hh>
28
30#include <dumux/io/name.hh>
36
37namespace Dumux {
38namespace FluidSystems {
39
45template< class Scalar, class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> >
46class Brine : public Base< Scalar, Brine<Scalar, H2OType>>
47{
50
51public:
53 using H2O = H2OType;
55
56 static const int numPhases = 1;
57 static const int numComponents = 2;
58
59 static constexpr int phase0Idx = 0;
60 static constexpr int liquidPhaseIdx = phase0Idx;
61
62 static constexpr int H2OIdx = 0;
63 static constexpr int NaClIdx = 1;
64 static constexpr int comp0Idx = H2OIdx;
65 static constexpr int comp1Idx = NaClIdx;
66
71 static const std::string phaseName(int phaseIdx = liquidPhaseIdx)
72 {
73 assert(phaseIdx == liquidPhaseIdx);
74 return IOName::liquidPhase();
75 }
76
81 static constexpr bool isMiscible()
82 {
83 return false;
84 }
85
90 static constexpr bool isGas(int phaseIdx = liquidPhaseIdx)
91 {
92 assert(phaseIdx == liquidPhaseIdx);
93 return false;
94 }
95
110 static bool isIdealMixture(int phaseIdx = liquidPhaseIdx)
111 {
112 assert(phaseIdx == liquidPhaseIdx);
113 return true;
114 }
115
125 static bool isCompressible(int phaseIdx = liquidPhaseIdx)
126 {
127 assert(phaseIdx == liquidPhaseIdx);
128 return H2O::liquidIsCompressible();
129 }
130
138 static bool isIdealGas(int phaseIdx = liquidPhaseIdx)
139 {
140 assert(phaseIdx == liquidPhaseIdx);
141 return false; /*we're a liquid!*/
142 }
143
144 /****************************************
145 * Component related static parameters
146 ****************************************/
147
152 static std::string componentName(int compIdx)
153 {
154 switch (compIdx)
155 {
156 case H2OIdx: return H2O::name();
157 case NaClIdx: return NaCl::name();
158 };
159 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
160 }
161
166 static Scalar molarMass(int compIdx)
167 {
168 switch (compIdx)
169 {
170 case H2OIdx: return H2O::molarMass();
171 case NaClIdx: return NaCl::molarMass();
172 };
173 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
174 }
175
176 /****************************************
177 * thermodynamic relations
178 ****************************************/
184 static void init()
185 {
186 init(/*tempMin=*/273.15,
187 /*tempMax=*/623.15,
188 /*numTemp=*/100,
189 /*pMin=*/-10.,
190 /*pMax=*/20e6,
191 /*numP=*/200);
192 }
193
205 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
206 Scalar pressMin, Scalar pressMax, unsigned nPress)
207 {
208
209 if (H2O::isTabulated)
210 {
211 std::cout << "Initializing tables for the H2O fluid properties ("
212 << nTemp*nPress
213 << " entries).\n";
214
215 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
216 }
217 }
218
219 using Base::density;
229 template <class FluidState>
230 static Scalar density(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
231 {
232 assert(phaseIdx == liquidPhaseIdx);
233 const Scalar temperature = fluidState.temperature(phaseIdx);
234 const Scalar pressure = fluidState.pressure(phaseIdx);
235 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
236
237 using std::max;
238 const Scalar TempC = temperature - 273.15;
239 const Scalar pMPa = pressure/1.0E6;
240 const Scalar salinity = max(0.0, xNaCl);
241
242 const Scalar rhow = H2O::liquidDensity(temperature, pressure);
243 const Scalar density = rhow + 1000*salinity*(0.668 +
244 0.44*salinity +
245 1.0E-6*(300*pMPa -
246 2400*pMPa*salinity +
247 TempC*(80.0 +
248 3*TempC -
249 3300*salinity -
250 13*pMPa +
251 47*pMPa*salinity)));
252 assert(density > 0.0);
253 return density;
254 }
255
264 template <class FluidState>
265 static Scalar fugacityCoefficient(const FluidState &fluidState,
266 int phaseIdx,
267 int compIdx)
268 {
269 assert(0 <= phaseIdx && phaseIdx < numPhases);
270 assert(0 <= compIdx && compIdx < numComponents);
271
272 if (phaseIdx == compIdx)
273 // We could calculate the real fugacity coefficient of
274 // the component in the fluid. Probably that's not worth
275 // the effort, since the fugacity coefficient of the other
276 // component is infinite anyway...
277 return 1.0;
278 return std::numeric_limits<Scalar>::infinity();
279 }
280
281 using Base::viscosity;
291 template <class FluidState>
292 static Scalar viscosity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
293 {
294 assert(phaseIdx == liquidPhaseIdx);
295 const Scalar temperature = fluidState.temperature(phaseIdx);
296 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
297
298 using std::pow;
299 using Dune::power;
300 using std::exp;
301 using std::max;
302 const Scalar T = max(temperature, 275.0);
303 const Scalar salinity = max(0.0, xNaCl);
304
305 const Scalar T_C = T - 273.15;
306 const Scalar A = ((0.42*power((pow(salinity, 0.8)-0.17), 2)) + 0.045)*pow(T_C, 0.8);
307 const Scalar mu_brine = 0.1 + (0.333*salinity) + (1.65+(91.9*salinity*salinity*salinity))*exp(-A); // [cP]
308 assert(mu_brine > 0.0);
309 return mu_brine/1000.0; // [Pa·s]
310 }
311
322 template <class FluidState>
323 static Scalar vaporPressure(const FluidState& fluidState, int compIdx)
324 {
325 if (compIdx == H2OIdx)
326 {
327 const Scalar temperature = fluidState.temperature(H2OIdx);
328 // Raoult's law, see Thomas Fetzer's Dissertation Eq. 2.11.
329 return H2O::vaporPressure(temperature)*fluidState.moleFraction(phase0Idx, H2OIdx);
330 }
331 else if (compIdx == NaClIdx)
332 DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)");
333 else
334 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
335 }
336
337 using Base::enthalpy;
351 template <class FluidState>
352 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
353 {
354 //use private enthalpy function to recycle it for the heat capacity calculation
355 return enthalpy_(fluidState.pressure(phaseIdx),
356 fluidState.temperature(phaseIdx),
357 fluidState.massFraction(phase0Idx, NaClIdx)); /*J/kg*/
358 }
359
367 template <class FluidState>
368 static Scalar componentEnthalpy(const FluidState &fluidState,
369 int phaseIdx,
370 int componentIdx)
371 {
372 const Scalar T = fluidState.temperature(liquidPhaseIdx);
373 const Scalar p = fluidState.pressure(liquidPhaseIdx);
374
375 if (phaseIdx == liquidPhaseIdx)
376 {
377 if (componentIdx == H2OIdx)
378 return H2O::liquidEnthalpy(T, p);
379 else if (componentIdx == NaClIdx)
380 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NaCl is not implemented.");
381 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
382 }
383 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
384 }
385
386 using Base::molarDensity;
396 template <class FluidState>
397 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
398 {
399 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
400 }
401
407 template <class FluidState>
408 static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
409 {
410 DUNE_THROW(Dune::NotImplemented, "FluidSystems::Brine::diffusionCoefficient()");
411 }
412
428 template <class FluidState>
429 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
430 int phaseIdx,
431 int compIIdx,
432 int compJIdx)
433 {
434 if (phaseIdx == liquidPhaseIdx)
435 {
436 if (compIIdx > compJIdx)
437 {
438 using std::swap;
439 swap(compIIdx, compJIdx);
440 }
441
442 if (compJIdx == NaClIdx)
443 return 1.54e-9;
444 else
445 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
446 << compIIdx << " and " << compJIdx
447 << " in phase " << phaseIdx);
448 }
449
450 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
451 }
452
461 template <class FluidState>
462 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
463 {
464 if (phaseIdx == liquidPhaseIdx)
465 {
466 Scalar tempC = fluidState.temperature(phaseIdx)-273.15;
467 Scalar m = fluidState.moleFraction(phaseIdx, NaClIdx)/(molarMass(H2OIdx)*(1- fluidState.moleFraction(phaseIdx, NaClIdx))); // molality of NaCl
468 Scalar S = 5844.3 * m / (1000 + 58.443 *m);
469 Scalar contribNaClFactor = 1.0 - (2.3434e-3 - 7.924e-6*tempC + 3.924e-8*tempC*tempC)*S + (1.06e-5 - 2.0e-8*tempC + 1.2e-10*tempC*tempC)*S*S;
470 return contribNaClFactor * H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
471 }
472 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
473 }
474
475 using Base::heatCapacity;
483 template <class FluidState>
484 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
485 {
486 if (phaseIdx == liquidPhaseIdx){
487 const Scalar eps = fluidState.temperature(phaseIdx)*1e-8;
488 //calculate heat capacity from the difference in enthalpy with temperature at constant pressure.
489 return (enthalpy_(fluidState.pressure(phaseIdx),
490 fluidState.temperature(phaseIdx) +eps ,
491 fluidState.massFraction(phase0Idx, NaClIdx))
492 - enthalpy_(fluidState.pressure(phaseIdx),
493 fluidState.temperature(phaseIdx),
494 fluidState.massFraction(phase0Idx, NaClIdx)))/eps; /*J/kg*/
495 }
496 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
497 }
498
499private:
510 static Scalar enthalpy_(const Scalar p, const Scalar T, const Scalar xNaCl)
511 {
512 /*Numerical coefficients from PALLISER*/
513 static const Scalar f[] = { 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10 };
514
515 /*Numerical coefficients from MICHAELIDES for the enthalpy of brine*/
516 static const Scalar a[4][3] = { { +9633.6, -4080.0, +286.49 },
517 { +166.58, +68.577, -4.6856 },
518 { -0.90963, -0.36524, +0.249667E-1 },
519 { +0.17965E-2, +0.71924E-3, -0.4900E-4 } };
520
521 const Scalar theta = T - 273.15;
522 const Scalar salSat = f[0] + f[1]*theta + f[2]*theta*theta + f[3]*theta*theta*theta;
523
524 /*Regularization*/
525 using std::min;
526 using std::max;
527 const Scalar salinity = min(max(xNaCl,0.0), salSat);
528
529 const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3; /* kJ/kg */
530
531 /*component enthalpy of soluted NaCl after DAUBERT and DANNER*/
532 const Scalar h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
533 + ((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* U [kJ/kg] */
534
535 const Scalar m = (1E3/58.44)*(salinity/(1-salinity));
536
537 using Dune::power;
538 Scalar d_h = 0;
539 for (int i = 0; i<=3; i++)
540 for (int j=0; j<=2; j++)
541 d_h = d_h + a[i][j] * power(theta, i) * power(m, j);
542
543 /* heat of dissolution for halite according to Michaelides 1981 */
544 const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
545
546 /* Enthalpy of brine without any dissolved gas */
547 const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */
548 return h_ls1*1E3; /*J/kg*/
549 }
550};
551
552} // end namespace FluidSystems
553} // end namespace Dumux
554
555#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
A central place for various physical constants occurring in some equations.
Material properties of pure water .
Tabulates all thermodynamic properties of a given untabulated chemical species.
Material properties of pure salt .
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:119
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
A class for the NaCl properties.
Definition: nacl.hh:47
static std::string name()
A human readable name for the NaCl.
Definition: nacl.hh:52
static constexpr Scalar molarMass()
The molar mass of NaCl in .
Definition: nacl.hh:60
Fluid system base class.
Definition: fluidsystems/base.hh:45
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculate the density of a fluid phase.
Definition: fluidsystems/base.hh:134
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/base.hh:390
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/base.hh:197
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: fluidsystems/base.hh:278
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/base.hh:326
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: fluidsystems/base.hh:363
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/base.hh:160
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: fluidsystems/base.hh:236
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: fluidsystems/base.hh:424
A compositional single phase fluid system consisting of two components, which are H2O and NaCl.
Definition: fluidsystems/brine.hh:47
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/brine.hh:462
static Scalar viscosity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the viscosity of the phase.
Definition: fluidsystems/brine.hh:292
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: fluidsystems/brine.hh:368
static bool isIdealGas(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: fluidsystems/brine.hh:138
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the diffusion coefficient of a component in a phase.
Definition: fluidsystems/brine.hh:408
static constexpr int phase0Idx
Index of the first (and only) phase.
Definition: fluidsystems/brine.hh:59
static constexpr bool isGas(int phaseIdx=liquidPhaseIdx)
Return whether a phase is gaseous.
Definition: fluidsystems/brine.hh:90
static constexpr int comp0Idx
index of the first component
Definition: fluidsystems/brine.hh:64
static const int numComponents
Number of components in the fluid system (H2O, NaCl)
Definition: fluidsystems/brine.hh:57
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: fluidsystems/brine.hh:205
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/brine.hh:265
static constexpr int comp1Idx
index of the second component
Definition: fluidsystems/brine.hh:65
static constexpr int H2OIdx
index of the water component
Definition: fluidsystems/brine.hh:62
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: fluidsystems/brine.hh:484
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: fluidsystems/brine.hh:166
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: fluidsystems/brine.hh:352
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: fluidsystems/brine.hh:323
static const std::string phaseName(int phaseIdx=liquidPhaseIdx)
Return the human readable name of the phase.
Definition: fluidsystems/brine.hh:71
static constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:63
H2OType H2O
export the involved components
Definition: fluidsystems/brine.hh:53
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: fluidsystems/brine.hh:81
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
The molar density of the fluid phase in .
Definition: fluidsystems/brine.hh:397
static constexpr int liquidPhaseIdx
The one considered phase is liquid.
Definition: fluidsystems/brine.hh:60
static bool isCompressible(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: fluidsystems/brine.hh:125
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: fluidsystems/brine.hh:152
static Scalar density(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the phase density [kg/m^3].
Definition: fluidsystems/brine.hh:230
static bool isIdealMixture(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: fluidsystems/brine.hh:110
static const int numPhases
Number of phases in the fluid system.
Definition: fluidsystems/brine.hh:56
static void init()
Initialize the fluid system's static parameters generically.
Definition: fluidsystems/brine.hh:184
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/brine.hh:429
Fluid system base class.