3.2-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
28
33
35
36#include <dumux/io/name.hh>
37
38namespace Dumux {
39namespace FluidSystems {
40
46template< class Scalar, class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> >
47class Brine : public Base< Scalar, Brine<Scalar, H2OType>>
48{
51
52public:
54 using H2O = H2OType;
56
57 static const int numPhases = 1;
58 static const int numComponents = 2;
59
60 static constexpr int phase0Idx = 0;
61 static constexpr int liquidPhaseIdx = phase0Idx;
62
63 static constexpr int H2OIdx = 0;
64 static constexpr int NaClIdx = 1;
65 static constexpr int comp0Idx = H2OIdx;
66 static constexpr int comp1Idx = NaClIdx;
67
72 static const std::string phaseName(int phaseIdx = liquidPhaseIdx)
73 {
74 assert(phaseIdx == liquidPhaseIdx);
75 return IOName::liquidPhase();
76 }
77
82 static constexpr bool isMiscible()
83 {
84 return false;
85 }
86
91 static constexpr bool isGas(int phaseIdx = liquidPhaseIdx)
92 {
93 assert(phaseIdx == liquidPhaseIdx);
94 return false;
95 }
96
111 static bool isIdealMixture(int phaseIdx = liquidPhaseIdx)
112 {
113 assert(phaseIdx == liquidPhaseIdx);
114 return true;
115 }
116
126 static bool isCompressible(int phaseIdx = liquidPhaseIdx)
127 {
128 assert(phaseIdx == liquidPhaseIdx);
129 return H2O::liquidIsCompressible();
130 }
131
139 static bool isIdealGas(int phaseIdx = liquidPhaseIdx)
140 {
141 assert(phaseIdx == liquidPhaseIdx);
142 return false; /*we're a liquid!*/
143 }
144
145 /****************************************
146 * Component related static parameters
147 ****************************************/
148
153 static std::string componentName(int compIdx)
154 {
155 switch (compIdx)
156 {
157 case H2OIdx: return H2O::name();
158 case NaClIdx: return NaCl::name();
159 };
160 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
161 }
162
167 static Scalar molarMass(int compIdx)
168 {
169 switch (compIdx)
170 {
171 case H2OIdx: return H2O::molarMass();
172 case NaClIdx: return NaCl::molarMass();
173 };
174 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
175 }
176
177 /****************************************
178 * thermodynamic relations
179 ****************************************/
185 static void init()
186 {
187 init(/*tempMin=*/273.15,
188 /*tempMax=*/623.15,
189 /*numTemp=*/100,
190 /*pMin=*/-10.,
191 /*pMax=*/20e6,
192 /*numP=*/200);
193 }
194
206 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
207 Scalar pressMin, Scalar pressMax, unsigned nPress)
208 {
209
210 if (H2O::isTabulated)
211 {
212 std::cout << "Initializing tables for the H2O fluid properties ("
213 << nTemp*nPress
214 << " entries).\n";
215
216 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
217 }
218 }
219
220 using Base::density;
225 template <class FluidState>
226 static Scalar density(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
227 {
228 assert(phaseIdx == liquidPhaseIdx);
229 const Scalar temperature = fluidState.temperature(phaseIdx);
230 const Scalar pressure = fluidState.pressure(phaseIdx);
231 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
232
233 using std::max;
234 const Scalar TempC = temperature - 273.15;
235 const Scalar pMPa = pressure/1.0E6;
236 const Scalar salinity = max(0.0, xNaCl);
237
238 const Scalar rhow = H2O::liquidDensity(temperature, pressure);
239 const Scalar density = rhow + 1000*salinity*(0.668 +
240 0.44*salinity +
241 1.0E-6*(300*pMPa -
242 2400*pMPa*salinity +
243 TempC*(80.0 +
244 3*TempC -
245 3300*salinity -
246 13*pMPa +
247 47*pMPa*salinity)));
248 assert(density > 0.0);
249 return density;
250 }
251
260 template <class FluidState>
261 static Scalar fugacityCoefficient(const FluidState &fluidState,
262 int phaseIdx,
263 int compIdx)
264 {
265 assert(0 <= phaseIdx && phaseIdx < numPhases);
266 assert(0 <= compIdx && compIdx < numComponents);
267
268 if (phaseIdx == compIdx)
269 // We could calculate the real fugacity coefficient of
270 // the component in the fluid. Probably that's not worth
271 // the effort, since the fugacity coefficient of the other
272 // component is infinite anyway...
273 return 1.0;
274 return std::numeric_limits<Scalar>::infinity();
275 }
276
277 using Base::viscosity;
281 template <class FluidState>
282 static Scalar viscosity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
283 {
284 assert(phaseIdx == liquidPhaseIdx);
285 const Scalar temperature = fluidState.temperature(phaseIdx);
286 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
287
288 using std::pow;
289 using std::exp;
290 using std::max;
291 const Scalar T = max(temperature, 275.0);
292 const Scalar salinity = max(0.0, xNaCl);
293
294 const Scalar T_C = T - 273.15;
295 const Scalar A = ((0.42*pow((pow(salinity, 0.8)-0.17), 2)) + 0.045)*pow(T_C, 0.8);
296 const Scalar mu_brine = 0.1 + (0.333*salinity) + (1.65+(91.9*salinity*salinity*salinity))*exp(-A); // [cP]
297 assert(mu_brine > 0.0);
298 return mu_brine/1000.0; // [Pa·s]
299 }
300
307 template <class FluidState>
308 static Scalar vaporPressure(const FluidState& fluidState, int compIdx)
309 {
310 if (compIdx == H2OIdx)
311 {
312 // simplified version of Eq 2.29 in Vishal Jambhekar's Promo
313 const Scalar temperature = fluidState.temperature(H2OIdx);
314 return H2O::vaporPressure(temperature)/fluidState.massFraction(phase0Idx, H2OIdx);
315 }
316 else if (compIdx == NaClIdx)
317 DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)");
318 else
319 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
320 }
321
322 using Base::enthalpy;
336 template <class FluidState>
337 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
338 {
339 /*Numerical coefficients from PALLISER*/
340 static const Scalar f[] = { 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10 };
341
342 /*Numerical coefficients from MICHAELIDES for the enthalpy of brine*/
343 static const Scalar a[4][3] = { { +9633.6, -4080.0, +286.49 },
344 { +166.58, +68.577, -4.6856 },
345 { -0.90963, -0.36524, +0.249667E-1 },
346 { +0.17965E-2, +0.71924E-3, -0.4900E-4 } };
347
348 const Scalar p = fluidState.pressure(phaseIdx);
349 const Scalar T = fluidState.temperature(phaseIdx);
350 const Scalar theta = T - 273.15;
351 const Scalar salSat = f[0] + f[1]*theta + f[2]*theta*theta + f[3]*theta*theta*theta;
352
353 /*Regularization*/
354 using std::min;
355 using std::max;
356 const Scalar xNaCl = fluidState.massFraction(phase0Idx, NaClIdx);
357 const Scalar salinity = min(max(xNaCl,0.0), salSat);
358
359 const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3; /* kJ/kg */
360
361 /*DAUBERT and DANNER*/
362 const Scalar h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
363 + ((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* U [kJ/kg] */
364
365 const Scalar m = (1E3/58.44)*(salinity/(1-salinity));
366
367 using std::pow;
368 Scalar d_h = 0;
369 for (int i = 0; i<=3; i++)
370 for (int j=0; j<=2; j++)
371 d_h = d_h + a[i][j] * pow(theta, i) * pow(m, j);
372
373 /* heat of dissolution for halite according to Michaelides 1971 */
374 const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
375
376 /* Enthalpy of brine without any dissolved gas */
377 const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */
378 return h_ls1*1E3; /*J/kg*/
379 }
380
388 template <class FluidState>
389 static Scalar componentEnthalpy(const FluidState &fluidState,
390 int phaseIdx,
391 int componentIdx)
392 {
393 const Scalar T = fluidState.temperature(liquidPhaseIdx);
394 const Scalar p = fluidState.pressure(liquidPhaseIdx);
395
396 if (phaseIdx == liquidPhaseIdx)
397 {
398 if (componentIdx == H2OIdx)
399 return H2O::liquidEnthalpy(T, p);
400 else if (componentIdx == NaClIdx)
401 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NaCl is not implemented.");
402 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
403 }
404 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
405 }
406
407 using Base::molarDensity;
417 template <class FluidState>
418 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
419 {
420 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
421 }
422
428 template <class FluidState>
429 static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
430 {
431 DUNE_THROW(Dune::NotImplemented, "FluidSystems::Brine::diffusionCoefficient()");
432 }
433
444 template <class FluidState>
445 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
446 int phaseIdx,
447 int compIIdx,
448 int compJIdx)
449 {
450 if (phaseIdx == liquidPhaseIdx)
451 {
452 if (compIIdx > compJIdx)
453 {
454 using std::swap;
455 swap(compIIdx, compJIdx);
456 }
458 // http://webserver.dmt.upm.es/~isidoro/dat1/Mass%20diffusivity%20data.htm
459 // The link above was given as a reference in brine_air fluid system.
460 // Doesn't work anymore though...
461 if (compJIdx == NaClIdx)
462 return 0.12e-9;
463 else
464 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
465 << compIIdx << " and " << compJIdx
466 << " in phase " << phaseIdx);
467 }
468
469 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
470 }
471
482 template <class FluidState>
483 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
484 {
485 if (phaseIdx == liquidPhaseIdx)
486 return H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
487 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
488 }
489
490 using Base::heatCapacity;
500 template <class FluidState>
501 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
502 {
503 if (phaseIdx == liquidPhaseIdx)
504 return H2O::liquidHeatCapacity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
505 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
506 }
507};
508
509} // end namespace FluidSystems
510} // end namespace Dumux
511
512#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
Material properties of pure water .
Material properties of pure salt .
Tabulates all thermodynamic properties of a given untabulated chemical species.
A central place for various physical constants occuring in some equations.
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:48
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/brine.hh:483
static Scalar viscosity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the viscosity of the phase.
Definition: fluidsystems/brine.hh:282
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:389
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:139
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the diffusion coefficient of a component in a phase.
Definition: fluidsystems/brine.hh:429
static constexpr int phase0Idx
Index of the first (and only) phase.
Definition: fluidsystems/brine.hh:60
static constexpr bool isGas(int phaseIdx=liquidPhaseIdx)
Return whether a phase is gaseous.
Definition: fluidsystems/brine.hh:91
static constexpr int comp0Idx
index of the first component
Definition: fluidsystems/brine.hh:65
static const int numComponents
Number of components in the fluid system (H2O, NaCl)
Definition: fluidsystems/brine.hh:58
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:206
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:261
static constexpr int comp1Idx
index of the second component
Definition: fluidsystems/brine.hh:66
static constexpr int H2OIdx
index of the water component
Definition: fluidsystems/brine.hh:63
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: fluidsystems/brine.hh:501
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: fluidsystems/brine.hh:167
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: fluidsystems/brine.hh:337
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: fluidsystems/brine.hh:308
static const std::string phaseName(int phaseIdx=liquidPhaseIdx)
Return the human readable name of the phase.
Definition: fluidsystems/brine.hh:72
static constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:64
H2OType H2O
export the involved components
Definition: fluidsystems/brine.hh:54
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: fluidsystems/brine.hh:82
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
The molar density of the fluid phase in .
Definition: fluidsystems/brine.hh:418
static constexpr int liquidPhaseIdx
The one considered phase is liquid.
Definition: fluidsystems/brine.hh:61
static bool isCompressible(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: fluidsystems/brine.hh:126
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: fluidsystems/brine.hh:153
static Scalar density(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the phase density [kg/m^3].
Definition: fluidsystems/brine.hh:226
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:111
static const int numPhases
Number of phases in the fluid system.
Definition: fluidsystems/brine.hh:57
static void init()
Initialize the fluid system's static parameters generically.
Definition: fluidsystems/brine.hh:185
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:445
Fluid system base class.