3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
brineco2.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_CO2_FLUID_SYSTEM_HH
25#define DUMUX_BRINE_CO2_FLUID_SYSTEM_HH
26
27#include <type_traits>
28
29#include <dune/common/exceptions.hh>
30
36
40
42
43#include <dumux/io/name.hh>
44
46
54 template<bool useConstantSalinity>
56
62 template<>
63 struct BrineCO2Indices<true>
64 {
65 static constexpr int BrineIdx = 0;
66 };
67
73 template<>
74 struct BrineCO2Indices<false>
75 {
76 static constexpr int H2OIdx = 0;
77 static constexpr int NaClIdx = 2;
78 static constexpr int comp2Idx = 2;
79 };
80
81} // end namespace Dumux::FluidSystems::Detail
82
83
84namespace Dumux::FluidSystems {
85
90template<bool salinityIsConstant, bool fastButSimplifiedRelations = false>
92{
93 static constexpr bool useConstantSalinity() { return salinityIsConstant; }
94 static constexpr bool useCO2GasDensityAsGasMixtureDensity() { return fastButSimplifiedRelations; }
95};
96
107template< class Scalar,
108 class CO2Table,
110 class Policy = BrineCO2DefaultPolicy</*constantSalinity?*/true> >
112: public Base<Scalar, BrineCO2<Scalar, CO2Table, H2OType, Policy>>
113, public Detail::BrineCO2Indices<Policy::useConstantSalinity()>
114{
117
118 // binary coefficients
120
121 // use constant salinity brine?
122 static constexpr bool useConstantSalinity = Policy::useConstantSalinity();
123
124 // The possible brine types
127 using BrineType = typename std::conditional_t< useConstantSalinity,
130
139 static constexpr int BrineOrH2OIdx = 0;
141 static constexpr int NaClIdx = 2;
142
143public:
145
146 using H2O = H2OType;
147 using Brine = BrineType;
149
150 static constexpr int numComponents = useConstantSalinity ? 2 : 3;
151 static constexpr int numPhases = 2;
152
153 static constexpr int liquidPhaseIdx = 0;
154 static constexpr int gasPhaseIdx = 1;
155 static constexpr int phase0Idx = liquidPhaseIdx;
156 static constexpr int phase1Idx = gasPhaseIdx;
157
158 static constexpr int comp0Idx = 0;
159 static constexpr int comp1Idx = 1;
160
161 // CO2 is always the second component
162 static constexpr int CO2Idx = comp1Idx;
163
164private:
165
166 // Adapter policy for the fluid state corresponding to the brine fluid system
167 struct BrineAdapterPolicy
168 {
169 using FluidSystem = VariableSalinityBrine;
170
171 static constexpr int phaseIdx(int brinePhaseIdx) { return liquidPhaseIdx; }
172 static constexpr int compIdx(int brineCompIdx)
173 {
174 assert(brineCompIdx == VariableSalinityBrine::H2OIdx || brineCompIdx == VariableSalinityBrine::NaClIdx);
175 switch (brineCompIdx)
176 {
177 case VariableSalinityBrine::H2OIdx: return BrineOrH2OIdx;
178 case VariableSalinityBrine::NaClIdx: return NaClIdx;
179 default: return 0; // this will never be reached, only needed to suppress compiler warning
180 }
181 }
182 };
183
184 template<class FluidState>
186
187public:
188
193 static std::string phaseName(int phaseIdx)
194 {
195 switch (phaseIdx)
196 {
197 case liquidPhaseIdx: return IOName::liquidPhase();
198 case gasPhaseIdx: return IOName::gaseousPhase();
199 }
200 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
201 }
202
206 static constexpr bool isMiscible()
207 { return true; }
208
214 static constexpr bool isGas(int phaseIdx)
215 {
216 assert(0 <= phaseIdx && phaseIdx < numPhases);
217
218 return phaseIdx == gasPhaseIdx;
219 }
220
226 static constexpr bool isIdealGas(int phaseIdx)
227 {
228 assert(0 <= phaseIdx && phaseIdx < numPhases);
229 // let the fluids decide
230 if (phaseIdx == gasPhaseIdx)
231 return useConstantSalinity ? (ConstantSalinityBrine::gasIsIdeal() && CO2::gasIsIdeal())
232 : (H2O::gasIsIdeal() && CO2::gasIsIdeal());
233 return false; // not a gas
234 }
235
250 static bool isIdealMixture(int phaseIdx)
251 {
252 assert(0 <= phaseIdx && phaseIdx < numPhases);
253 if (phaseIdx == liquidPhaseIdx)
254 return false;
255 return true;
256 }
257
267 static constexpr bool isCompressible(int phaseIdx)
268 {
269 assert(0 <= phaseIdx && phaseIdx < numPhases);
270 if (phaseIdx == liquidPhaseIdx)
271 return useConstantSalinity ? ConstantSalinityBrine::liquidIsCompressible()
273 return true;
274 }
275
280 static std::string componentName(int compIdx)
281 {
282 assert(0 <= compIdx && compIdx < numComponents);
283 if (useConstantSalinity)
284 {
285 static std::string name[] = { ConstantSalinityBrine::name(), CO2::name() };
286 return name[compIdx];
287 }
288 else
289 {
291 CO2::name(),
293 return name[compIdx];
294 }
295 }
296
301 static Scalar molarMass(int compIdx)
302 {
303 assert(0 <= compIdx && compIdx < numComponents);
304 if (useConstantSalinity)
305 {
306 static const Scalar M[] = { ConstantSalinityBrine::molarMass(), CO2::molarMass() };
307 return M[compIdx];
308 }
309 else
310 {
314 return M[compIdx];
315 }
316 }
317
318 /****************************************
319 * thermodynamic relations
320 ****************************************/
321
322 // Initializing with salinity and default tables
323 static void init()
324 {
325 init(/*startTemp=*/273.15, /*endTemp=*/623.15, /*tempSteps=*/100,
326 /*startPressure=*/1e4, /*endPressure=*/40e6, /*pressureSteps=*/200);
327 }
328
329 // Initializing with custom tables
330 static void init(Scalar startTemp, Scalar endTemp, int tempSteps,
331 Scalar startPressure, Scalar endPressure, int pressureSteps)
332 {
333 std::cout << "The Brine-CO2 fluid system was configured with the following policy:\n";
334 std::cout << " - use constant salinity: " << std::boolalpha << Policy::useConstantSalinity() << "\n";
335 std::cout << " - use CO2 gas density as gas mixture density: " << std::boolalpha << Policy::useCO2GasDensityAsGasMixtureDensity() << std::endl;
336
337 if (H2O::isTabulated)
338 H2O::init(startTemp, endTemp, tempSteps, startPressure, endPressure, pressureSteps);
339 }
340
341 using Base::density;
350 template <class FluidState>
351 static Scalar density(const FluidState& fluidState, int phaseIdx)
352 {
353 Scalar T = fluidState.temperature(phaseIdx);
354 if (phaseIdx == liquidPhaseIdx)
355 return liquidDensityMixture_(fluidState);
356
357 else if (phaseIdx == gasPhaseIdx)
358 {
359 if (Policy::useCO2GasDensityAsGasMixtureDensity())
360 // use the CO2 gas density only and neglect compositional effects
361 return CO2::gasDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
362 else
363 {
364 // assume ideal mixture: steam and CO2 don't "see" each other
365 Scalar rho_gH2O = H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, BrineOrH2OIdx));
366 Scalar rho_gCO2 = CO2::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, CO2Idx));
367 return (rho_gH2O + rho_gCO2);
368 }
369 }
370
371 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
372 }
373
374 using Base::molarDensity;
384 template <class FluidState>
385 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx)
386 {
387 Scalar T = fluidState.temperature(phaseIdx);
388 if (phaseIdx == liquidPhaseIdx)
389 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
390 else if (phaseIdx == gasPhaseIdx)
391 {
392 if (Policy::useCO2GasDensityAsGasMixtureDensity())
393 return CO2::gasMolarDensity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
394 else
395 {
396 // assume ideal mixture: steam and CO2 don't "see" each other
397 Scalar rhoMolar_gH2O = H2O::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, BrineOrH2OIdx));
398 Scalar rhoMolar_gCO2 = CO2::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, CO2Idx));
399 return rhoMolar_gH2O + rhoMolar_gCO2;
400 }
401 }
402 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
403 }
404
405 using Base::viscosity;
416 template <class FluidState>
417 static Scalar viscosity(const FluidState& fluidState, int phaseIdx)
418 {
419 Scalar T = fluidState.temperature(phaseIdx);
420 Scalar p = fluidState.pressure(phaseIdx);
421
422 if (phaseIdx == liquidPhaseIdx)
423 return useConstantSalinity ? ConstantSalinityBrine::liquidViscosity(T, p)
426 else if (phaseIdx == gasPhaseIdx)
427 return CO2::gasViscosity(T, p);
428
429 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
430 }
431
459 template <class FluidState>
460 static Scalar fugacityCoefficient(const FluidState& fluidState,
461 int phaseIdx,
462 int compIdx)
463 {
464 assert(0 <= compIdx && compIdx < numComponents);
465
466 if (phaseIdx == gasPhaseIdx)
467 // use the fugacity coefficients of an ideal gas. the
468 // actual value of the fugacity is not relevant, as long
469 // as the relative fluid compositions are observed,
470 return 1.0;
471
472 else if (phaseIdx == liquidPhaseIdx)
473 {
474 Scalar T = fluidState.temperature(phaseIdx);
475 Scalar pl = fluidState.pressure(liquidPhaseIdx);
476 Scalar pg = fluidState.pressure(gasPhaseIdx);
477
478 assert(T > 0);
479 assert(pl > 0 && pg > 0);
480
481 // calculate the equilibrium composition for given T & p
482 Scalar xlH2O, xgH2O;
483 Scalar xlCO2, xgCO2;
484 const Scalar salinity = useConstantSalinity ? ConstantSalinityBrine::salinity()
485 : fluidState.massFraction(liquidPhaseIdx, NaClIdx);
486 Brine_CO2::calculateMoleFractions(T, pl, salinity, /*knownGasPhaseIdx=*/-1, xlCO2, xgH2O);
487
488 // normalize the phase compositions
489 using std::min;
490 using std::max;
491 xlCO2 = max(0.0, min(1.0, xlCO2));
492 xgH2O = max(0.0, min(1.0, xgH2O));
493 xlH2O = 1.0 - xlCO2;
494 xgCO2 = 1.0 - xgH2O;
495
496 if (compIdx == BrineOrH2OIdx)
497 return (xgH2O/xlH2O)*(pg/pl);
498
499 else if (compIdx == CO2Idx)
500 return (xgCO2/xlCO2)*(pg/pl);
501
502 // NaCl is assumed to stay in the liquid!
503 else if (!useConstantSalinity && compIdx == NaClIdx)
504 return 0.0;
505
506 DUNE_THROW(Dune::InvalidStateException, "Invalid component index.");
507 }
508
509 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
510 }
511
519 template <class FluidState>
520 static Scalar equilibriumMoleFraction(const FluidState& fluidState,
521 const ParameterCache& paramCache,
522 int phaseIdx)
523 {
524 Scalar T = fluidState.temperature(phaseIdx);
525 Scalar p = fluidState.pressure(phaseIdx);
526
527 assert(T > 0);
528 assert(p > 0);
529
530 Scalar xgH2O;
531 Scalar xlCO2;
532
533 // calculate the equilibrium composition for given T & p
534 const Scalar salinity = useConstantSalinity ? ConstantSalinityBrine::salinity()
535 : fluidState.massFraction(liquidPhaseIdx, NaClIdx);
536 Brine_CO2::calculateMoleFractions(T, p, salinity, /*knowgasPhaseIdx=*/-1, xlCO2, xgH2O);
537
538 if (phaseIdx == gasPhaseIdx)
539 return xgH2O;
540 else if (phaseIdx == liquidPhaseIdx)
541 return xlCO2;
542
543 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
544 }
545
546
573 template <class FluidState>
574 static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
575 { DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients"); }
576
586 template <class FluidState>
587 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
588 int phaseIdx,
589 int compIIdx,
590 int compJIdx)
591 {
592 assert(0 <= compIIdx && compIIdx < numComponents);
593 assert(0 <= compJIdx && compJIdx < numComponents);
594
595 if (compIIdx > compJIdx)
596 {
597 using std::swap;
598 swap(compIIdx, compJIdx);
599 }
600
601 Scalar T = fluidState.temperature(phaseIdx);
602 Scalar p = fluidState.pressure(phaseIdx);
603 if (phaseIdx == liquidPhaseIdx)
604 {
605 if (compIIdx == BrineOrH2OIdx && compJIdx == CO2Idx)
606 return Brine_CO2::liquidDiffCoeff(T, p);
607 if (!useConstantSalinity && compIIdx == BrineOrH2OIdx && compJIdx == NaClIdx)
612
613 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components " <<
614 compIIdx << " and " << compJIdx << " in phase " << phaseIdx);
615 }
616 else if (phaseIdx == gasPhaseIdx)
617 {
618 if (compIIdx == BrineOrH2OIdx && compJIdx == CO2Idx)
619 return Brine_CO2::gasDiffCoeff(T, p);
620
621 // NaCl is expected to never be present in the gas phase. we need to
622 // return a diffusion coefficient that does not case numerical problems.
623 // We choose a very small value here.
624 else if (!useConstantSalinity && compIIdx == CO2Idx && compJIdx == NaClIdx)
625 return 1e-12;
626
627 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components " <<
628 compIIdx << " and " << compJIdx << " in phase " << phaseIdx);
629 }
630
631 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
632 }
633
634 using Base::enthalpy;
641 template <class FluidState>
642 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
643 {
644 Scalar T = fluidState.temperature(phaseIdx);
645 Scalar p = fluidState.pressure(phaseIdx);
646
647 if (phaseIdx == liquidPhaseIdx)
648 {
649 // Convert J/kg to kJ/kg
650 const Scalar h_ls1 = useConstantSalinity ? ConstantSalinityBrine::liquidEnthalpy(T, p)/1e3
653
654 // mass fraction of CO2 in Brine
655 const Scalar X_CO2_w = fluidState.massFraction(liquidPhaseIdx, CO2Idx);
656
657 // heat of dissolution for CO2 according to Fig. 6 in Duan and Sun 2003. (kJ/kg)
658 // In the relevant temperature ranges CO2 dissolution is exothermal
659 const Scalar delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
660
661 // enthalpy contribution of water and CO2 (kJ/kg)
662 const Scalar hw = H2O::liquidEnthalpy(T, p)/1e3;
663 const Scalar hg = CO2::liquidEnthalpy(T, p)/1e3 + delta_hCO2;
664
665 // Enthalpy of brine with dissolved CO2 (kJ/kg)
666 return (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1e3;
667 }
668 else if (phaseIdx == gasPhaseIdx)
669 {
670 // we assume NaCl to not enter the gas phase, only consider H2O and CO2
671 return H2O::gasEnthalpy(T, p)*fluidState.massFraction(gasPhaseIdx, BrineOrH2OIdx)
672 + CO2::gasEnthalpy(T, p) *fluidState.massFraction(gasPhaseIdx, CO2Idx);
673 }
674
675 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
676 }
677
684 template <class FluidState>
685 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
686 {
687 const Scalar T = fluidState.temperature(phaseIdx);
688 const Scalar p = fluidState.pressure(phaseIdx);
689
690 if (phaseIdx == liquidPhaseIdx)
691 {
692 if (componentIdx == BrineOrH2OIdx)
693 return H2O::liquidEnthalpy(T, p);
694 else if (componentIdx == CO2Idx)
695 return CO2::liquidEnthalpy(T, p);
696 else if (componentIdx == NaClIdx)
697 DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NaCl is not implemented.");
698 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
699 }
700 else if (phaseIdx == gasPhaseIdx)
701 {
702 if (componentIdx == BrineOrH2OIdx)
703 return H2O::gasEnthalpy(T, p);
704 else if (componentIdx == CO2Idx)
705 return CO2::gasEnthalpy(T, p);
706 else if (componentIdx == NaClIdx)
707 DUNE_THROW(Dune::InvalidStateException, "Implementation assumes NaCl not to be present in gas phase");
708 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
709 }
710 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
711 }
712
723 template <class FluidState>
724 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
725 {
726 if (phaseIdx == liquidPhaseIdx)
727 return useConstantSalinity ? ConstantSalinityBrine::liquidThermalConductivity( fluidState.temperature(phaseIdx),
728 fluidState.pressure(phaseIdx) )
731 else if (phaseIdx == gasPhaseIdx)
732 return CO2::gasThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
733
734 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
735 }
736
737 using Base::heatCapacity;
748 template <class FluidState>
749 static Scalar heatCapacity(const FluidState &fluidState,
750 int phaseIdx)
751 {
752 if(phaseIdx == liquidPhaseIdx)
753 return useConstantSalinity ? ConstantSalinityBrine::liquidHeatCapacity( fluidState.temperature(phaseIdx),
754 fluidState.pressure(phaseIdx) )
757 else if (phaseIdx == gasPhaseIdx)
758 return CO2::liquidHeatCapacity(fluidState.temperature(phaseIdx),
759 fluidState.pressure(phaseIdx));
760
761 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index.");
762 }
763
764private:
765
772 template<class FluidState>
773 static Scalar liquidDensityMixture_(const FluidState& fluidState)
774 {
775 const auto T = fluidState.temperature(liquidPhaseIdx);
776 const auto p = fluidState.pressure(liquidPhaseIdx);
777
778 if (T < 273.15)
779 DUNE_THROW(NumericalProblem, "Liquid density for Brine and CO2 is only "
780 "defined above 273.15K (T = " << T << ")");
781
782 if (p >= 2.5e8)
783 DUNE_THROW(NumericalProblem, "Liquid density for Brine and CO2 is only "
784 "defined below 250MPa (p = " << p << ")");
785
786 // density of pure water
787 Scalar rho_pure = H2O::liquidDensity(T, p);
788
789 // density of water with dissolved CO2 (neglect NaCl)
790 Scalar rho_lCO2;
791 if (useConstantSalinity)
792 {
793 // use normalized composition for to calculate the density
794 // (the relations don't seem to take non-normalized compositions too well...)
795 // TODO Do we really need this normalization???
796 using std::min;
797 using std::max;
798 Scalar xlBrine = min(1.0, max(0.0, fluidState.moleFraction(liquidPhaseIdx, BrineOrH2OIdx)));
799 Scalar xlCO2 = min(1.0, max(0.0, fluidState.moleFraction(liquidPhaseIdx, CO2Idx)));
800 Scalar sumx = xlBrine + xlCO2;
801 xlBrine /= sumx;
802 xlCO2 /= sumx;
803
804 rho_lCO2 = liquidDensityWaterCO2_(T, p, xlBrine, xlCO2);
805 }
806 else
807 {
808 // rescale mole fractions
809 auto xlH2O = fluidState.moleFraction(liquidPhaseIdx, BrineOrH2OIdx);
810 auto xlCO2 = fluidState.moleFraction(liquidPhaseIdx, NaClIdx);
811 const auto sumMoleFrac = xlH2O + xlCO2;
812 xlH2O = xlH2O/sumMoleFrac;
813 xlCO2 = xlCO2/sumMoleFrac;
814 rho_lCO2 = liquidDensityWaterCO2_(T, p, xlH2O, xlCO2);
815 }
816
817 // density of brine (water with nacl)
818 Scalar rho_brine = useConstantSalinity ? ConstantSalinityBrine::liquidDensity(T, p)
819 : VariableSalinityBrine::density( BrineAdapter<FluidState>(fluidState),
820 VariableSalinityBrine::liquidPhaseIdx );
821
822 // contribution of co2 to the density
823 Scalar contribCO2 = rho_lCO2 - rho_pure;
824
825 // Total brine density with dissolved CO2
826 // rho_{b, CO2} = rho_pure + contribution(salt) + contribution(CO2)
827 return rho_brine + contribCO2;
828 }
829
830
841 static Scalar liquidDensityWaterCO2_(Scalar temperature,
842 Scalar pl,
843 Scalar xlH2O,
844 Scalar xlCO2)
845 {
846 const Scalar M_CO2 = CO2::molarMass();
847 const Scalar M_H2O = H2O::molarMass();
848
849 const Scalar tempC = temperature - 273.15; /* tempC : temperature in °C */
850 const Scalar rho_pure = H2O::liquidDensity(temperature, pl);
851
852 // xlH2O is available, but in case of a pure gas phase
853 // the value of M_T for the virtual liquid phase can become very large
854 xlH2O = 1.0 - xlCO2;
855 const Scalar M_T = M_H2O * xlH2O + M_CO2 * xlCO2;
856 const Scalar V_phi = (37.51 +
857 tempC*(-9.585e-2 +
858 tempC*(8.74e-4 -
859 tempC*5.044e-7))) / 1.0e6;
860 return 1/(xlCO2 * V_phi/M_T + M_H2O * xlH2O / (rho_pure * M_T));
861 }
862};
863
864} // end namespace Dumux::FluidSystems
865
866#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A collection of input/output field names for common physical quantities.
Relations valid for an ideal gas.
Adapter class for fluid states with different indices.
A class for the CO2 fluid properties.
Tabulates all thermodynamic properties of a given untabulated chemical species.
Binary coefficients for CO2 and brine.
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:119
Definition: h2o.hh:914
Definition: h2o.hh:914
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Binary coefficients for brine and CO2.
Definition: brine_co2.hh:42
static void calculateMoleFractions(const Scalar temperature, const Scalar pg, const Scalar salinity, const int knownPhaseIdx, Scalar &xlCO2, Scalar &ygH2O)
Returns the mol (!) fraction of CO2 in the liquid phase and the mol (!) fraction of H2O in the gas ph...
Definition: brine_co2.hh:113
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient of CO2 in the brine phase.
Definition: brine_co2.hh:83
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient of water in the CO2 phase.
Definition: brine_co2.hh:57
A class for the brine fluid properties.
Definition: components/brine.hh:56
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure brine.
Definition: components/brine.hh:419
static const Scalar liquidEnthalpy(Scalar T, Scalar p)
Specific enthalpy of liquid brine .
Definition: components/brine.hh:161
static Scalar molarMass()
The molar mass in of brine. This assumes that the salt is pure NaCl.
Definition: components/brine.hh:83
static Scalar salinity()
Return the constant salinity.
Definition: components/brine.hh:73
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of a brine .
Definition: components/brine.hh:443
static const Scalar liquidHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of brine .
Definition: components/brine.hh:220
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: components/brine.hh:303
static std::string name()
A human readable name for the brine.
Definition: components/brine.hh:67
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: components/brine.hh:291
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure brine at a given pressure and temperature .
Definition: components/brine.hh:318
A class for the CO2 fluid properties.
Definition: co2.hh:57
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of CO2.
Definition: co2.hh:399
static Scalar liquidEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of liquid CO2 .
Definition: co2.hh:180
static constexpr Scalar molarMass()
The mass in of one mole of CO2.
Definition: co2.hh:72
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of CO2 gas in at a given pressure and temperature.
Definition: co2.hh:244
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of CO2 at a given pressure and temperature .
Definition: co2.hh:226
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: co2.hh:126
static std::string name()
A human readable name for the CO2.
Definition: co2.hh:66
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of CO2. Equations given in: - Vesovic et al., 1990.
Definition: co2.hh:327
static Scalar liquidHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of the component as a liquid. USE WITH CAUTION! Exploits enthalpy fu...
Definition: co2.hh:304
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of gaseous CO2 .
Definition: co2.hh:162
static std::string name()
A human readable name for the NaCl.
Definition: nacl.hh:52
Tabulates all thermodynamic properties of a given component.
Definition: tabulatedcomponent.hh:684
Adapter class for fluid states with different indices.
Definition: adapter.hh:44
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 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 constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:63
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 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
Class that exports some indices that should be provided by the BrineAir fluid system....
Definition: brineco2.hh:55
Default policy for the Brine-CO2 fluid system.
Definition: brineco2.hh:92
static constexpr bool useCO2GasDensityAsGasMixtureDensity()
Definition: brineco2.hh:94
static constexpr bool useConstantSalinity()
Definition: brineco2.hh:93
A compositional fluid with brine (H2O & NaCl) and carbon dioxide as components in both the liquid and...
Definition: brineco2.hh:114
static constexpr int CO2Idx
Definition: brineco2.hh:162
static constexpr int numComponents
Definition: brineco2.hh:150
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given the phase composition, return the specific phase enthalpy .
Definition: brineco2.hh:642
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: brineco2.hh:267
static constexpr int phase0Idx
index of the first phase
Definition: brineco2.hh:155
static constexpr int comp0Idx
Definition: brineco2.hh:158
H2OType H2O
Definition: brineco2.hh:146
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: brineco2.hh:301
static constexpr int numPhases
Definition: brineco2.hh:151
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: brineco2.hh:749
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: brineco2.hh:280
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given the phase compositions, return the binary diffusion coefficient of two components in a phase.
Definition: brineco2.hh:587
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: brineco2.hh:351
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: brineco2.hh:417
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the molecular diffusion coefficient for a component in a fluid phase .
Definition: brineco2.hh:574
static void init()
Definition: brineco2.hh:323
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: brineco2.hh:214
static constexpr int phase1Idx
index of the second phase
Definition: brineco2.hh:156
static constexpr int gasPhaseIdx
index of the gas phase
Definition: brineco2.hh:154
static constexpr bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: brineco2.hh:226
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: brineco2.hh:193
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: brineco2.hh:724
static constexpr int liquidPhaseIdx
index of the liquid phase
Definition: brineco2.hh:153
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: brineco2.hh:460
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: brineco2.hh:206
BrineType Brine
Definition: brineco2.hh:147
static void init(Scalar startTemp, Scalar endTemp, int tempSteps, Scalar startPressure, Scalar endPressure, int pressureSteps)
Definition: brineco2.hh:330
static constexpr int comp1Idx
Definition: brineco2.hh:159
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: brineco2.hh:685
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: brineco2.hh:250
static Scalar equilibriumMoleFraction(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
Returns the equilibrium concentration of the dissolved component in a phase.
Definition: brineco2.hh:520
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: brineco2.hh:385
The a parameter cache which does nothing.
Definition: nullparametercache.hh:34
Fluid system base class.
A class for the brine fluid properties,.
A fluid system for brine, i.e. H2O with dissolved NaCl.