My Project
Loading...
Searching...
No Matches
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"
33
35
37
43
47
48#include <string_view>
49
50namespace Opm {
51
52// Silence compiler warnings about use of variables
53// that are instantiated in a different compilation unit.
54template<>
55const float CO2<float>::brineSalinity;
56template<>
57const double CO2<double>::brineSalinity;
58
66template <class Scalar>
68 : public BaseFluidSystem<Scalar, BrineCO2FluidSystem<Scalar> >
69{
70 typedef ::Opm::H2O<Scalar> H2O_IAPWS;
71 typedef ::Opm::Brine<Scalar, H2O_IAPWS> Brine_IAPWS;
74
75 typedef H2O_Tabulated H2O;
76
77public:
78 template <class Evaluation>
79 struct ParameterCache : public NullParameterCache<Evaluation>
80 {};
81
85 typedef ::Opm::CO2<Scalar> CO2;
86
89
90 /****************************************
91 * Fluid phase related static parameters
92 ****************************************/
93
95 static const int numPhases = 2;
96
98 static const int liquidPhaseIdx = 0;
100 static const int gasPhaseIdx = 1;
101
105 static std::string_view phaseName(unsigned phaseIdx)
106 {
107 static const std::string_view name[] = {
108 "liquid",
109 "gas"
110 };
111
112 assert(phaseIdx < numPhases);
113 return name[phaseIdx];
114 }
115
119 static bool isLiquid(unsigned phaseIdx)
120 {
121 assert(phaseIdx < numPhases);
122
123 return phaseIdx != gasPhaseIdx;
124 }
125
129 static bool isIdealGas(unsigned phaseIdx)
130 {
131 assert(phaseIdx < numPhases);
132
133 if (phaseIdx == gasPhaseIdx)
134 return CO2::gasIsIdeal();
135 return false;
136 }
137
141 static bool isIdealMixture([[maybe_unused]] unsigned phaseIdx)
142 {
143 assert(phaseIdx < numPhases);
144
145 return true;
146 }
147
151 static bool isCompressible([[maybe_unused]] unsigned phaseIdx)
152 {
153 assert(phaseIdx < numPhases);
154
155 return true;
156 }
157
158 /****************************************
159 * Component related static parameters
160 ****************************************/
162 static const int numComponents = 2;
163
165 static const int BrineIdx = 0;
167 static const int CO2Idx = 1;
168
172 static std::string_view componentName(unsigned compIdx)
173 {
174 static const std::string_view name[] = {
175 Brine::name(),
176 CO2::name(),
177 };
178
179 assert(compIdx < numComponents);
180 return name[compIdx];
181 }
182
186 static Scalar molarMass(unsigned compIdx)
187 {
188 assert(compIdx < numComponents);
189 return (compIdx==BrineIdx)
191 : CO2::molarMass();
192 }
193
194 /****************************************
195 * thermodynamic relations
196 ****************************************/
197
201 static void init()
202 {
203 init(/*startTemp=*/273.15, /*endTemp=*/623.15, /*tempSteps=*/50,
204 /*startPressure=*/1e4, /*endPressure=*/40e6, /*pressureSteps=*/50);
205 }
206
218 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
219 Scalar pressMin, Scalar pressMax, unsigned nPress)
220 {
221 if (H2O::isTabulated) {
222 H2O_Tabulated::init(tempMin, tempMax, nTemp,
223 pressMin, pressMax, nPress);
224 }
225
226 // set the salinity of brine to the one used by the CO2 tables
227 Brine_IAPWS::salinity = CO2::brineSalinity;
228
229 if (Brine::isTabulated) {
230 Brine_Tabulated::init(tempMin, tempMax, nTemp,
231 pressMin, pressMax, nPress);
232 }
233 }
234
238 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
239 static LhsEval density(const FluidState& fluidState,
240 const ParameterCache<ParamCacheEval>& /*paramCache*/,
241 unsigned phaseIdx)
242 {
243 assert(phaseIdx < numPhases);
244
245 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
246 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
247
248 if (phaseIdx == liquidPhaseIdx) {
249 // use normalized composition for to calculate the density
250 // (the relations don't seem to take non-normalized
251 // compositions too well...)
252 LhsEval xlBrine = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, BrineIdx))));
253 LhsEval xlCO2 = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(liquidPhaseIdx, CO2Idx))));
254 LhsEval sumx = xlBrine + xlCO2;
255 xlBrine /= sumx;
256 xlCO2 /= sumx;
257
258 LhsEval result = liquidDensity_(temperature,
259 pressure,
260 xlBrine,
261 xlCO2);
262
263 Valgrind::CheckDefined(result);
264 return result;
265 }
266
267 assert(phaseIdx == gasPhaseIdx);
268
269 // use normalized composition for to calculate the density
270 // (the relations don't seem to take non-normalized
271 // compositions too well...)
272 LhsEval xgBrine = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, BrineIdx))));
273 LhsEval xgCO2 = min(1.0, max(0.0, decay<LhsEval>(fluidState.moleFraction(gasPhaseIdx, CO2Idx))));
274 LhsEval sumx = xgBrine + xgCO2;
275 xgBrine /= sumx;
276 xgCO2 /= sumx;
277
278 LhsEval result = gasDensity_(temperature,
279 pressure,
280 xgBrine,
281 xgCO2);
282 Valgrind::CheckDefined(result);
283 return result;
284 }
285
289 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
290 static LhsEval viscosity(const FluidState& fluidState,
291 const ParameterCache<ParamCacheEval>& /*paramCache*/,
292 unsigned phaseIdx)
293 {
294 assert(phaseIdx < numPhases);
295
296 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
297 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
298
299 if (phaseIdx == liquidPhaseIdx) {
300 // assume pure brine for the liquid phase. TODO: viscosity
301 // of mixture
302 LhsEval result = Brine::liquidViscosity(temperature, pressure);
303 Valgrind::CheckDefined(result);
304 return result;
305 }
306
307 assert(phaseIdx == gasPhaseIdx);
308 LhsEval result = CO2::gasViscosity(temperature, pressure);
309 Valgrind::CheckDefined(result);
310 return result;
311 }
312
316 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
317 static LhsEval fugacityCoefficient(const FluidState& fluidState,
318 const ParameterCache<ParamCacheEval>& /*paramCache*/,
319 unsigned phaseIdx,
320 unsigned compIdx)
321 {
322 assert(phaseIdx < numPhases);
323 assert(compIdx < numComponents);
324
325 if (phaseIdx == gasPhaseIdx)
326 // use the fugacity coefficients of an ideal gas. the
327 // actual value of the fugacity is not relevant, as long
328 // as the relative fluid compositions are observed,
329 return 1.0;
330
331 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
332 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
333 assert(temperature > 0);
334 assert(pressure > 0);
335
336 // Activity model for salt-out effect in brine-CO2 mutual solubility
337 // 2 = Duan & Sun as modified in Spycher & Pruess, Trans. Porous Media, (2009)
338 const int activityModel = 3;
339
340 // calulate the equilibrium composition for the given
341 // temperature and pressure. TODO: calculateMoleFractions()
342 // could use some cleanup.
343 LhsEval xlH2O, xgH2O;
344 LhsEval xlCO2, xgCO2;
346 pressure,
347 LhsEval(Brine_IAPWS::salinity),
348 /*knownPhaseIdx=*/-1,
349 xlCO2,
350 xgH2O,
351 activityModel);
352
353 // normalize the phase compositions
354 xlCO2 = max(0.0, min(1.0, xlCO2));
355 xgH2O = max(0.0, min(1.0, xgH2O));
356
357 xlH2O = 1.0 - xlCO2;
358 xgCO2 = 1.0 - xgH2O;
359
360 if (compIdx == BrineIdx) {
361 Scalar phigH2O = 1.0;
362 return phigH2O * xgH2O / xlH2O;
363 }
364 else {
365 assert(compIdx == CO2Idx);
366
367 Scalar phigCO2 = 1.0;
368 return phigCO2 * xgCO2 / xlCO2;
369 };
370 }
371
375 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
376 static LhsEval diffusionCoefficient(const FluidState& fluidState,
377 const ParameterCache<ParamCacheEval>& /*paramCache*/,
378 unsigned phaseIdx,
379 unsigned /*compIdx*/)
380 {
381 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
382 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
383 if (phaseIdx == liquidPhaseIdx)
384 return BinaryCoeffBrineCO2::liquidDiffCoeff(temperature, pressure);
385
386 assert(phaseIdx == gasPhaseIdx);
387 return BinaryCoeffBrineCO2::gasDiffCoeff(temperature, pressure);
388 }
389
393 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
394 static LhsEval enthalpy(const FluidState& fluidState,
395 const ParameterCache<ParamCacheEval>& /*paramCache*/,
396 unsigned phaseIdx)
397 {
398 assert(phaseIdx < numPhases);
399
400 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
401 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
402
403 if (phaseIdx == liquidPhaseIdx) {
404 const LhsEval& XlCO2 = decay<LhsEval>(fluidState.massFraction(phaseIdx, CO2Idx));
405 const LhsEval& result = liquidEnthalpyBrineCO2_(temperature,
406 pressure,
408 XlCO2);
409 Valgrind::CheckDefined(result);
410 return result;
411 }
412 else {
413 const LhsEval& XCO2 = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, CO2Idx));
414 const LhsEval& XBrine = decay<LhsEval>(fluidState.massFraction(gasPhaseIdx, BrineIdx));
415
416 LhsEval result = 0;
417 result += XBrine * Brine::gasEnthalpy(temperature, pressure);
418 result += XCO2 * CO2::gasEnthalpy(temperature, pressure);
419 Valgrind::CheckDefined(result);
420 return result;
421 }
422 }
423
427 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
428 static LhsEval thermalConductivity(const FluidState& /*fluidState*/,
429 const ParameterCache<ParamCacheEval>& /*paramCache*/,
430 unsigned phaseIdx)
431 {
432 // TODO way too simple!
433 if (phaseIdx == liquidPhaseIdx)
434 return 0.6; // conductivity of water[W / (m K ) ]
435
436 // gas phase
437 return 0.025; // conductivity of air [W / (m K ) ]
438 }
439
452 template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
453 static LhsEval heatCapacity(const FluidState& fluidState,
454 const ParameterCache<ParamCacheEval>& /*paramCache*/,
455 unsigned phaseIdx)
456 {
457 assert(phaseIdx < numPhases);
458
459 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(phaseIdx));
460 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(phaseIdx));
461
462 if(phaseIdx == liquidPhaseIdx)
463 return H2O::liquidHeatCapacity(temperature, pressure);
464 else
465 return CO2::gasHeatCapacity(temperature, pressure);
466 }
467
468private:
469 template <class LhsEval>
470 static LhsEval gasDensity_(const LhsEval& T,
471 const LhsEval& pg,
472 const LhsEval& xgH2O,
473 const LhsEval& xgCO2)
474 {
475 Valgrind::CheckDefined(T);
476 Valgrind::CheckDefined(pg);
477 Valgrind::CheckDefined(xgH2O);
478 Valgrind::CheckDefined(xgCO2);
479
480 return CO2::gasDensity(T, pg);
481 }
482
483 /***********************************************************************/
484 /* */
485 /* Total brine density with dissolved CO2 */
486 /* rho_{b,CO2} = rho_w + contribution(salt) + contribution(CO2) */
487 /* */
488 /***********************************************************************/
489 template <class LhsEval>
490 static LhsEval liquidDensity_(const LhsEval& T,
491 const LhsEval& pl,
492 const LhsEval& xlH2O,
493 const LhsEval& xlCO2)
494 {
495 Valgrind::CheckDefined(T);
496 Valgrind::CheckDefined(pl);
497 Valgrind::CheckDefined(xlH2O);
498 Valgrind::CheckDefined(xlCO2);
499
500 auto cast = [](const auto d)
501 {
502#if HAVE_QUAD
503 if constexpr (std::is_same_v<decltype(d), const quad>)
504 return static_cast<double>(d);
505 else
506#endif
507 return d;
508 };
509
510 auto tostring = [cast](const auto& val) -> std::string
511 {
512 if constexpr (DenseAd::is_evaluation<LhsEval>::value) {
513 return std::to_string(cast(getValue(val.value())));
514 } else {
515 return std::to_string(cast(val));
516 }
517 };
518
519 if(T < 273.15) {
520 const std::string msg =
521 "Liquid density for Brine and CO2 is only "
522 "defined above 273.15K (is +"
523 + tostring(T)+ "K)";
524 throw NumericalProblem(msg);
525 }
526 if(pl >= 2.5e8) {
527 const std::string msg =
528 "Liquid density for Brine and CO2 is only "
529 "defined below 250MPa (is "
530 + tostring(pl) + "Pa)";
531 throw NumericalProblem(msg);
532 }
533
534 const LhsEval& rho_brine = Brine::liquidDensity(T, pl);
535 const LhsEval& rho_pure = H2O::liquidDensity(T, pl);
536 const LhsEval& rho_lCO2 = liquidDensityWaterCO2_(T, pl, xlH2O, xlCO2);
537 const LhsEval& contribCO2 = rho_lCO2 - rho_pure;
538
539 return rho_brine + contribCO2;
540 }
541
542 template <class LhsEval>
543 static LhsEval liquidDensityWaterCO2_(const LhsEval& temperature,
544 const LhsEval& pl,
545 const LhsEval& /*xlH2O*/,
546 const LhsEval& xlCO2)
547 {
548 Scalar M_CO2 = CO2::molarMass();
549 Scalar M_H2O = H2O::molarMass();
550
551 const LhsEval& tempC = temperature - 273.15; /* tempC : temperature in °C */
552 const LhsEval& rho_pure = H2O::liquidDensity(temperature, pl);
553 // calculate the mole fraction of CO2 in the liquid. note that xlH2O is available
554 // as a function parameter, but in the case of a pure gas phase the value of M_T
555 // for the virtual liquid phase can become very large
556 const LhsEval xlH2O = 1.0 - xlCO2;
557 const LhsEval& M_T = M_H2O * xlH2O + M_CO2 * xlCO2;
558 const LhsEval& V_phi =
559 (37.51 +
560 tempC*(-9.585e-2 +
561 tempC*(8.74e-4 -
562 tempC*5.044e-7))) / 1.0e6;
563 return 1/ (xlCO2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
564 }
565
566 template <class LhsEval>
567 static LhsEval liquidEnthalpyBrineCO2_(const LhsEval& T,
568 const LhsEval& p,
569 Scalar S, // salinity
570 const LhsEval& X_CO2_w)
571 {
572 /* X_CO2_w : mass fraction of CO2 in brine */
573
574 /* same function as enthalpy_brine, only extended by CO2 content */
575
576 /*Numerical coefficents from PALLISER*/
577 static Scalar f[] = {
578 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
579 };
580
581 /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
582 static Scalar a[4][3] = {
583 { 9633.6, -4080.0, +286.49 },
584 { +166.58, +68.577, -4.6856 },
585 { -0.90963, -0.36524, +0.249667E-1 },
586 { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
587 };
588
589 LhsEval theta, h_NaCl;
590 LhsEval h_ls1, d_h;
591 LhsEval delta_h;
592 LhsEval delta_hCO2, hg, hw;
593
594 theta = T - 273.15;
595
596 // Regularization
597 Scalar scalarTheta = scalarValue(theta);
598 Scalar S_lSAT = f[0] + scalarTheta*(f[1] + scalarTheta*(f[2] + scalarTheta*f[3]));
599 if (S > S_lSAT)
600 S = S_lSAT;
601
602 hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
603
604 /*DAUBERT and DANNER*/
605 /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
606 +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
607
608 Scalar m = 1E3/58.44 * S/(1-S);
609 int i = 0;
610 int j = 0;
611 d_h = 0;
612
613 for (i = 0; i<=3; i++) {
614 for (j=0; j<=2; j++) {
615 d_h = d_h + a[i][j] * pow(theta, static_cast<Scalar>(i)) * std::pow(m, j);
616 }
617 }
618 /* heat of dissolution for halite according to Michaelides 1971 */
619 delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
620
621 /* Enthalpy of brine without CO2 */
622 h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
623
624 /* heat of dissolution for CO2 according to Fig. 6 in Duan and Sun 2003. (kJ/kg)
625 In the relevant temperature ranges CO2 dissolution is
626 exothermal */
627 delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
628
629 /* enthalpy contribution of CO2 (kJ/kg) */
630 hg = CO2::gasEnthalpy(T, p)/1E3 + delta_hCO2;
631
632 /* Enthalpy of brine with dissolved CO2 */
633 return (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1E3; /*J/kg*/
634 }
635};
636
637} // namespace Opm
638
639#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.
Provides the OPM specific exception classes.
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:47
Scalar Scalar
The type used for scalar quantities.
Definition BaseFluidSystem.hpp:52
Binary coefficients for brine and CO2.
Definition Brine_CO2.hpp:45
static void calculateMoleFractions(const Evaluation &temperature, const Evaluation &pg, const Evaluation &salinity, const int knownPhaseIdx, Evaluation &xlCO2, Evaluation &ygH2O, const int &activityModel, 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:100
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:59
static Evaluation liquidDiffCoeff(const Evaluation &, const Evaluation &)
Binary diffusion coefficent [m^2/s] of CO2 in the brine phase.
Definition Brine_CO2.hpp:76
A two-phase fluid system with water and CO2.
Definition BrineCO2FluidSystem.hpp:69
static const int CO2Idx
The index of the CO2 component.
Definition BrineCO2FluidSystem.hpp:167
static const int BrineIdx
The index of the brine component.
Definition BrineCO2FluidSystem.hpp:165
Brine_Tabulated Brine
The type of the component for brine used by the fluid system.
Definition BrineCO2FluidSystem.hpp:83
static std::string_view componentName(unsigned compIdx)
Return the human readable name of a component.
Definition BrineCO2FluidSystem.hpp:172
static void init()
Initialize the fluid system's static parameters.
Definition BrineCO2FluidSystem.hpp:201
static LhsEval heatCapacity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Specific isobaric heat capacity of a fluid phase [J/kg].
Definition BrineCO2FluidSystem.hpp:453
static const int gasPhaseIdx
The index of the gas phase.
Definition BrineCO2FluidSystem.hpp:100
static LhsEval thermalConductivity(const FluidState &, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Thermal conductivity of a fluid phase [W/(m K)].
Definition BrineCO2FluidSystem.hpp:428
static const int numPhases
The number of phases considered by the fluid system.
Definition BrineCO2FluidSystem.hpp:95
BinaryCoeff::Brine_CO2< Scalar, H2O, CO2 > BinaryCoeffBrineCO2
The binary coefficients for brine and CO2 used by this fluid system.
Definition BrineCO2FluidSystem.hpp:88
::Opm::CO2< Scalar > CO2
The type of the component for pure CO2 used by the fluid system.
Definition BrineCO2FluidSystem.hpp:85
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition BrineCO2FluidSystem.hpp:186
static bool isIdealMixture(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition BrineCO2FluidSystem.hpp:141
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition BrineCO2FluidSystem.hpp:151
static const int numComponents
Number of chemical species in the fluid system.
Definition BrineCO2FluidSystem.hpp:162
static const int liquidPhaseIdx
The index of the liquid phase.
Definition BrineCO2FluidSystem.hpp:98
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:394
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:218
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition BrineCO2FluidSystem.hpp:290
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition BrineCO2FluidSystem.hpp:129
static std::string_view phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition BrineCO2FluidSystem.hpp:105
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition BrineCO2FluidSystem.hpp:239
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:376
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition BrineCO2FluidSystem.hpp:119
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:317
A class for the brine fluid properties.
Definition Brine.hpp:48
static Scalar salinity
The mass fraction of salt assumed to be in the brine.
Definition Brine.hpp:51
static Evaluation gasHeatCapacity(const Evaluation &temperature, const Evaluation &pressure)
Specific isobaric heat capacity of the component [J/kg] as a liquid.
Definition CO2.hpp:261
static Evaluation gasViscosity(Evaluation temperature, const Evaluation &pressure, bool extrapolate=false)
The dynamic viscosity [Pa s] of CO2.
Definition CO2.hpp:208
static std::string_view name()
A human readable name for the CO2.
Definition CO2.hpp:65
static Evaluation gasEnthalpy(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
Specific enthalpy of gaseous CO2 [J/kg].
Definition CO2.hpp:169
static bool gasIsIdeal()
Returns true iff the gas phase is assumed to be ideal.
Definition CO2.hpp:162
static Scalar molarMass()
The mass in [kg] of one mole of CO2.
Definition CO2.hpp:71
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:194
Material properties of pure water .
Definition H2O.hpp:65
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 std::string_view name()
A human readable name for the component.
Definition TabulatedComponent.hpp:215
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 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
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition BrineCO2FluidSystem.hpp:80