version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_COMPOSITIONAL_FLUID_STATE_HH
15#define DUMUX_COMPOSITIONAL_FLUID_STATE_HH
16
17#include <algorithm>
18#include <cmath>
19#include <type_traits>
20#include <cassert>
21#include <array>
22
23#include <dune/common/exceptions.hh>
24
25namespace Dumux {
26
33template <class ScalarType, class FluidSystem>
35{
36public:
37 static constexpr int numPhases = FluidSystem::numPhases;
38 static constexpr int numComponents = FluidSystem::numComponents;
39
41 using Scalar = ScalarType;
42
45
47 template <class FluidState, typename std::enable_if_t<!std::is_same<FluidState, CompositionalFluidState>::value, int> = 0>
48 CompositionalFluidState(const FluidState &fs)
49 { assign(fs); }
50
51 // copy and move constructor / assignment operator
56
57 /*****************************************************
58 * Generic access to fluid properties (No assumptions
59 * on thermodynamic equilibrium required)
60 *****************************************************/
65 int wettingPhase() const { return wPhaseIdx_; }
66
75 Scalar saturation(int phaseIdx) const
76 { return saturation_[phaseIdx]; }
77
87 Scalar moleFraction(int phaseIdx, int compIdx) const
88 { return moleFraction_[phaseIdx][compIdx]; }
89
105 Scalar massFraction(int phaseIdx, int compIdx) const
106 {
107 // calculate the mass fractions:
108 // for "mass" models this is just a back calculation
109 return sumMoleFractions_[phaseIdx]
110 * moleFraction(phaseIdx, compIdx)
111 * FluidSystem::molarMass(compIdx)
112 / averageMolarMass_[phaseIdx];
113 }
114
119 Scalar phaseMassFraction(int phaseIdx) const
120 {
121 Scalar totalMass = 0.0;
122 for (int pIdx = 0; pIdx < numPhases; ++pIdx)
123 totalMass += saturation(pIdx)*density(pIdx);
124
125 return saturation(phaseIdx)*density(phaseIdx) / totalMass;
126 }
127
136 Scalar averageMolarMass(int phaseIdx) const
137 { return averageMolarMass_[phaseIdx]; }
138
148 Scalar molarity(int phaseIdx, int compIdx) const
149 { return molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); }
150
168 Scalar fugacity(int phaseIdx, int compIdx) const
169 { return fugacityCoefficient(phaseIdx, compIdx)*moleFraction(phaseIdx, compIdx)*pressure(phaseIdx); }
170
174 Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
175 { return fugacityCoefficient_[phaseIdx][compIdx]; }
176
182 Scalar molarVolume(int phaseIdx) const
183 { return 1.0/molarDensity(phaseIdx); }
184
189 Scalar density(int phaseIdx) const
190 { return density_[phaseIdx]; }
191
196 Scalar molarDensity(int phaseIdx) const
197 { return molarDensity_[phaseIdx]; }
198
202 Scalar temperature(int phaseIdx) const
203 { return temperature_[phaseIdx]; }
204
208 Scalar pressure(int phaseIdx) const
209 { return pressure_[phaseIdx]; }
210
215 Scalar partialPressure(int phaseIdx, int compIdx) const
216 {
217 assert(FluidSystem::isGas(phaseIdx));
218 return moleFraction(phaseIdx, compIdx) * pressure(phaseIdx);
219 }
220
224 Scalar enthalpy(int phaseIdx) const
225 { return enthalpy_[phaseIdx]; }
226
234 Scalar internalEnergy(int phaseIdx) const
235 { return enthalpy_[phaseIdx] - pressure(phaseIdx)/density(phaseIdx); }
236
240 Scalar viscosity(int phaseIdx) const
241 { return viscosity_[phaseIdx]; }
242
243
244 /*****************************************************
245 * Access to fluid properties which only make sense
246 * if assuming thermodynamic equilibrium
247 *****************************************************/
248
253 { return temperature_[0]; }
254
260 Scalar fugacity(int compIdx) const
261 { return fugacity(0, compIdx); }
262
263
264 /*****************************************************
265 * Setter methods. Note that these are not part of the
266 * generic FluidState interface but specific for each
267 * implementation...
268 *****************************************************/
269
278 template <class FluidState>
279 void assign(const FluidState &fs)
280 {
281 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
282 {
283 averageMolarMass_[phaseIdx] = 0;
284 sumMoleFractions_[phaseIdx] = 0;
285 temperature_[phaseIdx] = fs.temperature();
286 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
287 {
288 moleFraction_[phaseIdx][compIdx] = fs.moleFraction(phaseIdx, compIdx);
289 fugacityCoefficient_[phaseIdx][compIdx] = fs.fugacityCoefficient(phaseIdx, compIdx);
290 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx);
291 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx];
292 }
293 pressure_[phaseIdx] = fs.pressure(phaseIdx);
294 saturation_[phaseIdx] = fs.saturation(phaseIdx);
295 density_[phaseIdx] = fs.density(phaseIdx);
296 molarDensity_[phaseIdx] = fs.molarDensity(phaseIdx);
297 enthalpy_[phaseIdx] = fs.enthalpy(phaseIdx);
298 viscosity_[phaseIdx] = fs.viscosity(phaseIdx);
299 }
300 wPhaseIdx_ = fs.wettingPhase();
301 }
302
307 {
308 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
309 temperature_[phaseIdx] = value;
310 }
311
316 void setTemperature(const int phaseIdx, const Scalar value)
317 { temperature_[phaseIdx] = value; }
318
322 void setPressure(int phaseIdx, Scalar value)
323 { pressure_[phaseIdx] = value; }
324
328 void setSaturation(int phaseIdx, Scalar value)
329 { saturation_[phaseIdx] = value; }
330
336 void setMoleFraction(int phaseIdx, int compIdx, Scalar value)
337 {
338 moleFraction_[phaseIdx][compIdx] = value;
339
340 // re-calculate the mean molar mass
341 sumMoleFractions_[phaseIdx] = 0.0;
342 averageMolarMass_[phaseIdx] = 0.0;
343 for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx)
344 {
345 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
346 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
347 }
348 }
349
355 void setMassFraction(int phaseIdx, int compIdx, Scalar value)
356 {
357 if (numComponents != 2)
358 DUNE_THROW(Dune::NotImplemented, "This currently only works for 2 components.");
359 else
360 {
361 // calculate average molar mass of the gas phase
362 Scalar M1 = FluidSystem::molarMass(compIdx);
363 Scalar M2 = FluidSystem::molarMass(1-compIdx);
364 Scalar X2 = 1.0-value;
365 Scalar avgMolarMass = M1*M2/(M2 + X2*(M1 - M2));
366
367 moleFraction_[phaseIdx][compIdx] = value * avgMolarMass / M1;
368 moleFraction_[phaseIdx][1-compIdx] = 1.0-moleFraction_[phaseIdx][compIdx];
369
370 // re-calculate the mean molar mass
371 sumMoleFractions_[phaseIdx] = 0.0;
372 averageMolarMass_[phaseIdx] = 0.0;
373 for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx) {
374 sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx];
375 averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx);
376 }
377 }
378 }
379
384 void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
385 { fugacityCoefficient_[phaseIdx][compIdx] = value; }
386
390 void setDensity(int phaseIdx, Scalar value)
391 { density_[phaseIdx] = value; }
392
396 void setMolarDensity(int phaseIdx, Scalar value)
397 { molarDensity_[phaseIdx] = value; }
398
402 void setEnthalpy(int phaseIdx, Scalar value)
403 { enthalpy_[phaseIdx] = value; }
404
408 void setViscosity(int phaseIdx, Scalar value)
409 { viscosity_[phaseIdx] = value; }
410
414 void setWettingPhase(int phaseIdx)
415 { wPhaseIdx_ = phaseIdx; }
416
417protected:
419 std::array<std::array<Scalar, numComponents>, numPhases> moleFraction_ = {};
420 std::array<std::array<Scalar, numComponents>, numPhases> fugacityCoefficient_ = {};
421 std::array<Scalar, numPhases> averageMolarMass_ = {};
422 std::array<Scalar, numPhases> sumMoleFractions_ = {};
423 std::array<Scalar, numPhases> pressure_ = {};
424 std::array<Scalar, numPhases> saturation_ = {};
425 std::array<Scalar, numPhases> density_ = {};
426 std::array<Scalar, numPhases> molarDensity_ = {};
427 std::array<Scalar, numPhases> enthalpy_ = {};
428 std::array<Scalar, numPhases> viscosity_ = {};
429 std::array<Scalar, numPhases> temperature_ = {};
430
432};
433
434} // end namespace Dumux
435
436#endif
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Definition: compositional.hh:35
CompositionalFluidState & operator=(CompositionalFluidState &&fs)=default
int wPhaseIdx_
Definition: compositional.hh:431
void setMolarDensity(int phaseIdx, Scalar value)
Set the molar density of a phase .
Definition: compositional.hh:396
CompositionalFluidState(CompositionalFluidState &&fs)=default
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: compositional.hh:208
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: compositional.hh:87
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: compositional.hh:240
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: compositional.hh:202
CompositionalFluidState & operator=(const CompositionalFluidState &fs)=default
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: compositional.hh:234
Scalar molarDensity(int phaseIdx) const
The molar density of the fluid phase in .
Definition: compositional.hh:196
void setSaturation(int phaseIdx, Scalar value)
Set the saturation of a phase .
Definition: compositional.hh:328
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: compositional.hh:75
int wettingPhase() const
Returns the index of the most wetting phase in the fluid-solid configuration (for porous medium syste...
Definition: compositional.hh:65
void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
Set the fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:384
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: compositional.hh:279
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: compositional.hh:224
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: compositional.hh:189
void setWettingPhase(int phaseIdx)
Set the index of the most wetting phase.
Definition: compositional.hh:414
Scalar partialPressure(int phaseIdx, int compIdx) const
The partial pressure of a component in a phase .
Definition: compositional.hh:215
std::array< Scalar, numPhases > pressure_
Definition: compositional.hh:423
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:336
void setTemperature(Scalar value)
Set the temperature of all phases.
Definition: compositional.hh:306
std::array< std::array< Scalar, numComponents >, numPhases > fugacityCoefficient_
Definition: compositional.hh:420
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: compositional.hh:148
Scalar phaseMassFraction(int phaseIdx) const
Returns the phase mass fraction, i.e. phase mass per total mass .
Definition: compositional.hh:119
static constexpr int numPhases
Definition: compositional.hh:37
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: compositional.hh:136
CompositionalFluidState()=default
default constructor
std::array< Scalar, numPhases > averageMolarMass_
Definition: compositional.hh:421
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: compositional.hh:182
CompositionalFluidState(const FluidState &fs)
copy constructor from arbitrary fluid state
Definition: compositional.hh:48
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:355
std::array< Scalar, numPhases > saturation_
Definition: compositional.hh:424
Scalar temperature() const
The temperature within the domain .
Definition: compositional.hh:252
std::array< Scalar, numPhases > sumMoleFractions_
Definition: compositional.hh:422
std::array< Scalar, numPhases > density_
Definition: compositional.hh:425
std::array< Scalar, numPhases > enthalpy_
Definition: compositional.hh:427
static constexpr int numComponents
Definition: compositional.hh:38
Scalar fugacity(int compIdx) const
The fugacity of a component .
Definition: compositional.hh:260
std::array< Scalar, numPhases > temperature_
Definition: compositional.hh:429
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: compositional.hh:105
std::array< Scalar, numPhases > molarDensity_
Definition: compositional.hh:426
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:174
CompositionalFluidState(const CompositionalFluidState &fs)=default
std::array< Scalar, numPhases > viscosity_
Definition: compositional.hh:428
ScalarType Scalar
export the scalar type
Definition: compositional.hh:41
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:316
void setDensity(int phaseIdx, Scalar value)
Set the density of a phase .
Definition: compositional.hh:390
void setViscosity(int phaseIdx, Scalar value)
Set the dynamic viscosity of a phase .
Definition: compositional.hh:408
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: compositional.hh:168
void setEnthalpy(int phaseIdx, Scalar value)
Set the specific enthalpy of a phase .
Definition: compositional.hh:402
void setPressure(int phaseIdx, Scalar value)
Set the fluid pressure of a phase .
Definition: compositional.hh:322
std::array< std::array< Scalar, numComponents >, numPhases > moleFraction_
zero-initialize all data members with braces syntax
Definition: compositional.hh:419
Definition: adapt.hh:17