3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
compositional.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 *****************************************************************************/
26#ifndef DUMUX_COMPOSITIONAL_FLUID_STATE_HH
27#define DUMUX_COMPOSITIONAL_FLUID_STATE_HH
28
29#include <algorithm>
30#include <cmath>
31#include <type_traits>
32#include <cassert>
33#include <array>
34
35#include <dune/common/exceptions.hh>
36
37namespace Dumux {
38
45template <class ScalarType, class FluidSystem>
47{
48public:
49 static constexpr int numPhases = FluidSystem::numPhases;
50 static constexpr int numComponents = FluidSystem::numComponents;
51
53 using Scalar = ScalarType;
54
57
59 template <class FluidState, typename std::enable_if_t<!std::is_same<FluidState, CompositionalFluidState>::value, int> = 0>
60 CompositionalFluidState(const FluidState &fs)
61 { assign(fs); }
62
63 // copy and move constructor / assignment operator
68
69 /*****************************************************
70 * Generic access to fluid properties (No assumptions
71 * on thermodynamic equilibrium required)
72 *****************************************************/
77 int wettingPhase() const { return wPhaseIdx_; }
78
87 Scalar saturation(int phaseIdx) const
88 { return saturation_[phaseIdx]; }
89
99 Scalar moleFraction(int phaseIdx, int compIdx) const
100 { return moleFraction_[phaseIdx][compIdx]; }
101
117 Scalar massFraction(int phaseIdx, int compIdx) const
118 {
119 // calculate the mass fractions:
120 // for "mass" models this is just a back calculation
121 return sumMoleFractions_[phaseIdx]
122 * moleFraction(phaseIdx, compIdx)
123 * FluidSystem::molarMass(compIdx)
124 / averageMolarMass_[phaseIdx];
125 }
126
131 Scalar phaseMassFraction(int phaseIdx) const
132 {
133 Scalar totalMass = 0.0;
134 for (int pIdx = 0; pIdx < numPhases; ++pIdx)
135 totalMass += saturation(pIdx)*density(pIdx);
136
137 return saturation(phaseIdx)*density(phaseIdx) / totalMass;
138 }
139
148 Scalar averageMolarMass(int phaseIdx) const
149 { return averageMolarMass_[phaseIdx]; }
150
160 Scalar molarity(int phaseIdx, int compIdx) const
161 { return molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
162
180 Scalar fugacity(int phaseIdx, int compIdx) const
181 { return fugacityCoefficient(phaseIdx, compIdx)*moleFraction(phaseIdx, compIdx)*pressure(phaseIdx); }
182
186 Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
187 { return fugacityCoefficient_[phaseIdx][compIdx]; }
188
194 Scalar molarVolume(int phaseIdx) const
195 { return 1.0/molarDensity(phaseIdx); }
196
201 Scalar density(int phaseIdx) const
202 { return density_[phaseIdx]; }
203
208 Scalar molarDensity(int phaseIdx) const
209 { return molarDensity_[phaseIdx]; }
210
214 Scalar temperature(int phaseIdx) const
215 { return temperature_[phaseIdx]; }
216
220 Scalar pressure(int phaseIdx) const
221 { return pressure_[phaseIdx]; }
222
227 Scalar partialPressure(int phaseIdx, int compIdx) const
228 {
229 assert(FluidSystem::isGas(phaseIdx));
230 return moleFraction(phaseIdx, compIdx) * pressure(phaseIdx);
231 }
232
236 Scalar enthalpy(int phaseIdx) const
237 { return enthalpy_[phaseIdx]; }
238
246 Scalar internalEnergy(int phaseIdx) const
247 { return enthalpy_[phaseIdx] - pressure(phaseIdx)/density(phaseIdx); }
248
252 Scalar viscosity(int phaseIdx) const
253 { return viscosity_[phaseIdx]; }
254
255
256 /*****************************************************
257 * Access to fluid properties which only make sense
258 * if assuming thermodynamic equilibrium
259 *****************************************************/
260
265 { return temperature_[0]; }
266
272 Scalar fugacity(int compIdx) const
273 { return fugacity(0, compIdx); }
274
275
276 /*****************************************************
277 * Setter methods. Note that these are not part of the
278 * generic FluidState interface but specific for each
279 * implementation...
280 *****************************************************/
281
290 template <class FluidState>
291 void assign(const FluidState &fs)
292 {
293 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
294 {
295 averageMolarMass_[phaseIdx] = 0;
296 sumMoleFractions_[phaseIdx] = 0;
297 temperature_[phaseIdx] = fs.temperature();
298 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
299 {
300 moleFraction_[phaseIdx][compIdx] = fs.moleFraction(phaseIdx, compIdx);
301 fugacityCoefficient_[phaseIdx][compIdx] = fs.fugacityCoefficient(phaseIdx, compIdx);
302 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
303 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
304 }
305 pressure_[phaseIdx] = fs.pressure(phaseIdx);
306 saturation_[phaseIdx] = fs.saturation(phaseIdx);
307 density_[phaseIdx] = fs.density(phaseIdx);
308 molarDensity_[phaseIdx] = fs.molarDensity(phaseIdx);
309 enthalpy_[phaseIdx] = fs.enthalpy(phaseIdx);
310 viscosity_[phaseIdx] = fs.viscosity(phaseIdx);
311 }
312 wPhaseIdx_ = fs.wettingPhase();
313 }
314
319 {
320 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
321 temperature_[phaseIdx] = value;
322 }
323
328 void setTemperature(const int phaseIdx, const Scalar value)
329 { temperature_[phaseIdx] = value; }
330
334 void setPressure(int phaseIdx, Scalar value)
335 { pressure_[phaseIdx] = value; }
336
340 void setSaturation(int phaseIdx, Scalar value)
341 { saturation_[phaseIdx] = value; }
342
348 void setMoleFraction(int phaseIdx, int compIdx, Scalar value)
349 {
350 moleFraction_[phaseIdx][compIdx] = value;
351
352 // re-calculate the mean molar mass
353 sumMoleFractions_[phaseIdx] = 0.0;
354 averageMolarMass_[phaseIdx] = 0.0;
355 for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx)
356 {
357 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
358 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
359 }
360 }
361
367 void setMassFraction(int phaseIdx, int compIdx, Scalar value)
368 {
369 if (numComponents != 2)
370 DUNE_THROW(Dune::NotImplemented, "This currently only works for 2 components.");
371 else
372 {
373 // calculate average molar mass of the gas phase
374 Scalar M1 = FluidSystem::molarMass(compIdx);
375 Scalar M2 = FluidSystem::molarMass(1-compIdx);
376 Scalar X2 = 1.0-value;
377 Scalar avgMolarMass = M1*M2/(M2 + X2*(M1 - M2));
378
379 moleFraction_[phaseIdx][compIdx] = value * avgMolarMass / M1;
380 moleFraction_[phaseIdx][1-compIdx] = 1.0-moleFraction_[phaseIdx][compIdx];
381
382 // re-calculate the mean molar mass
383 sumMoleFractions_[phaseIdx] = 0.0;
384 averageMolarMass_[phaseIdx] = 0.0;
385 for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
386 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
387 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
388 }
389 }
390 }
391
397 template <class FluidState>
398 void setRelativeHumidity(FluidState &fluidState, int phaseIdx, int compIdx, Scalar value)
399 {
400 // asserts for the assumption under which setting the relative humidity is possible
401 assert(phaseIdx == FluidSystem::nPhaseIdx);
402 assert(compIdx == FluidSystem::wCompIdx);
403 assert(numComponents == 2);
404 assert(FluidSystem::isGas(phaseIdx));
405
406 Scalar moleFraction = value * FluidSystem::vaporPressure(fluidState, FluidSystem::wCompIdx)
407 / fluidState.pressure(phaseIdx);
408 fluidState.setMoleFraction(phaseIdx, FluidSystem::wCompIdx, moleFraction);
409 fluidState.setMoleFraction(phaseIdx, FluidSystem::nCompIdx, 1.0-moleFraction);
410 }
411
416 void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
417 { fugacityCoefficient_[phaseIdx][compIdx] = value; }
418
422 void setDensity(int phaseIdx, Scalar value)
423 { density_[phaseIdx] = value; }
424
428 void setMolarDensity(int phaseIdx, Scalar value)
429 { molarDensity_[phaseIdx] = value; }
430
434 void setEnthalpy(int phaseIdx, Scalar value)
435 { enthalpy_[phaseIdx] = value; }
436
440 void setViscosity(int phaseIdx, Scalar value)
441 { viscosity_[phaseIdx] = value; }
442
446 void setWettingPhase(int phaseIdx)
447 { wPhaseIdx_ = phaseIdx; }
448
449protected:
451 std::array<std::array<Scalar, numComponents>, numPhases> moleFraction_ = {};
452 std::array<std::array<Scalar, numComponents>, numPhases> fugacityCoefficient_ = {};
453 std::array<Scalar, numPhases> averageMolarMass_ = {};
454 std::array<Scalar, numPhases> sumMoleFractions_ = {};
455 std::array<Scalar, numPhases> pressure_ = {};
456 std::array<Scalar, numPhases> saturation_ = {};
457 std::array<Scalar, numPhases> density_ = {};
458 std::array<Scalar, numPhases> molarDensity_ = {};
459 std::array<Scalar, numPhases> enthalpy_ = {};
460 std::array<Scalar, numPhases> viscosity_ = {};
461 std::array<Scalar, numPhases> temperature_ = {};
462
464};
465
466} // end namespace Dumux
467
468#endif
Definition: adapt.hh:29
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Definition: compositional.hh:47
CompositionalFluidState & operator=(CompositionalFluidState &&fs)=default
int wPhaseIdx_
Definition: compositional.hh:463
void setMolarDensity(int phaseIdx, Scalar value)
Set the molar density of a phase .
Definition: compositional.hh:428
CompositionalFluidState(CompositionalFluidState &&fs)=default
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: compositional.hh:220
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: compositional.hh:99
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: compositional.hh:252
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: compositional.hh:214
CompositionalFluidState & operator=(const CompositionalFluidState &fs)=default
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: compositional.hh:246
Scalar molarDensity(int phaseIdx) const
The molar density of the fluid phase in .
Definition: compositional.hh:208
void setRelativeHumidity(FluidState &fluidState, int phaseIdx, int compIdx, Scalar value)
Set the relative humidity of a component in a phase and update the average molar mass according to ...
Definition: compositional.hh:398
void setSaturation(int phaseIdx, Scalar value)
Set the saturation of a phase .
Definition: compositional.hh:340
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: compositional.hh:87
int wettingPhase() const
Returns the index of the most wetting phase in the fluid-solid configuration (for porous medium syste...
Definition: compositional.hh:77
void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
Set the fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:416
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: compositional.hh:291
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: compositional.hh:236
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: compositional.hh:201
void setWettingPhase(int phaseIdx)
Set the index of the most wetting phase.
Definition: compositional.hh:446
Scalar partialPressure(int phaseIdx, int compIdx) const
The partial pressure of a component in a phase .
Definition: compositional.hh:227
std::array< Scalar, numPhases > pressure_
Definition: compositional.hh:455
void setMoleFraction(int phaseIdx, int compIdx, Scalar value)
Set the mole fraction of a component in a phase and update the average molar mass according to the ...
Definition: compositional.hh:348
void setTemperature(Scalar value)
Set the temperature of all phases.
Definition: compositional.hh:318
std::array< std::array< Scalar, numComponents >, numPhases > fugacityCoefficient_
Definition: compositional.hh:452
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: compositional.hh:160
Scalar phaseMassFraction(int phaseIdx) const
Returns the phase mass fraction, i.e. phase mass per total mass .
Definition: compositional.hh:131
static constexpr int numPhases
Definition: compositional.hh:49
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: compositional.hh:148
CompositionalFluidState()=default
default constructor
std::array< Scalar, numPhases > averageMolarMass_
Definition: compositional.hh:453
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: compositional.hh:194
CompositionalFluidState(const FluidState &fs)
copy constructor from arbitrary fluid state
Definition: compositional.hh:60
void setMassFraction(int phaseIdx, int compIdx, Scalar value)
Set the mass fraction of a component in a phase and update the average molar mass according to the ...
Definition: compositional.hh:367
std::array< Scalar, numPhases > saturation_
Definition: compositional.hh:456
Scalar temperature() const
The temperature within the domain .
Definition: compositional.hh:264
std::array< Scalar, numPhases > sumMoleFractions_
Definition: compositional.hh:454
std::array< Scalar, numPhases > density_
Definition: compositional.hh:457
std::array< Scalar, numPhases > enthalpy_
Definition: compositional.hh:459
static constexpr int numComponents
Definition: compositional.hh:50
Scalar fugacity(int compIdx) const
The fugacity of a component .
Definition: compositional.hh:272
std::array< Scalar, numPhases > temperature_
Definition: compositional.hh:461
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: compositional.hh:117
std::array< Scalar, numPhases > molarDensity_
Definition: compositional.hh:458
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:186
CompositionalFluidState(const CompositionalFluidState &fs)=default
std::array< Scalar, numPhases > viscosity_
Definition: compositional.hh:460
ScalarType Scalar
export the scalar type
Definition: compositional.hh:53
void setTemperature(const int phaseIdx, const Scalar value)
Set the temperature of a specific phase. This is not implemented in this fluidstate.
Definition: compositional.hh:328
void setDensity(int phaseIdx, Scalar value)
Set the density of a phase .
Definition: compositional.hh:422
void setViscosity(int phaseIdx, Scalar value)
Set the dynamic viscosity of a phase .
Definition: compositional.hh:440
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: compositional.hh:180
void setEnthalpy(int phaseIdx, Scalar value)
Set the specific enthalpy of a phase .
Definition: compositional.hh:434
void setPressure(int phaseIdx, Scalar value)
Set the fluid pressure of a phase .
Definition: compositional.hh:334
std::array< std::array< Scalar, numComponents >, numPhases > moleFraction_
zero-initialize all data members with braces syntax
Definition: compositional.hh:451