26#ifndef DUMUX_MISCIBLE_MULTIPHASE_COMPOSITION_HH
27#define DUMUX_MISCIBLE_MULTIPHASE_COMPOSITION_HH
29#include <dune/common/fvector.hh>
30#include <dune/common/fmatrix.hh>
57template <
class Scalar,
class Flu
idSystem>
60 static constexpr int numPhases = FluidSystem::numPhases;
61 static constexpr int numComponents = FluidSystem::numComponents;
62 static const int numMajorComponents = FluidSystem::numPhases;
80 template <
class Flu
idState,
class ParameterCache>
81 static void solve(FluidState &fluidState,
82 ParameterCache ¶mCache,
83 int knownPhaseIdx = 0)
90 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
91 assert(FluidSystem::isIdealMixture(phaseIdx));
98 Dune::FieldVector<Scalar, numComponents-numMajorComponents> xKnown(0.0);
99 for (
int knownCompIdx = 0; knownCompIdx < numComponents-numMajorComponents; ++knownCompIdx)
101 xKnown[knownCompIdx] = fluidState.moleFraction(knownPhaseIdx, knownCompIdx + numMajorComponents);
105 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
106 paramCache.updatePhase(fluidState, phaseIdx);
111 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
112 Scalar fugCoeff = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
113 fluidState.setFugacityCoefficient(phaseIdx, compIdx, fugCoeff);
120 Dune::FieldMatrix<Scalar, numComponents*numPhases, numComponents*numPhases> M(0.0);
121 Dune::FieldVector<Scalar, numComponents*numPhases> x(0.0);
122 Dune::FieldVector<Scalar, numComponents*numPhases> b(0.0);
126 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
127 int rowIdx = numComponents*(numPhases - 1) + phaseIdx;
130 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
131 int colIdx = phaseIdx*numComponents + compIdx;
133 M[rowIdx][colIdx] = 1.0;
140 for(
int knownCompIdx = 0; knownCompIdx < numComponents-numMajorComponents; ++knownCompIdx)
142 int rowIdx = numComponents + numPhases + knownCompIdx;
143 int colIdx = knownPhaseIdx*numComponents + knownCompIdx + numMajorComponents;
144 M[rowIdx][colIdx] = 1.0;
145 b[rowIdx] = xKnown[knownCompIdx];
150 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
152 int col1Idx = compIdx;
153 const auto entryPhase0 = fluidState.fugacityCoefficient(0, compIdx)*fluidState.pressure(0);
155 for (
unsigned int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
157 int rowIdx = (phaseIdx - 1)*numComponents + compIdx;
158 int col2Idx = phaseIdx*numComponents + compIdx;
159 M[rowIdx][col1Idx] = entryPhase0;
160 M[rowIdx][col2Idx] = -fluidState.fugacityCoefficient(phaseIdx, compIdx)*fluidState.pressure(phaseIdx);
165 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
168 if (compIdx < numMajorComponents)
177 try { M.solve(x, b); }
178 catch (
const Dune::FMatrixError& e) {
180 "MiscibleMultiPhaseComposition: Failed to solve the linear equation system for the phase composition.");
185 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
186 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
187 int rowIdx = phaseIdx*numComponents + compIdx;
188 fluidState.setMoleFraction(phaseIdx, compIdx, x[rowIdx]);
190 paramCache.updateComposition(fluidState, phaseIdx);
193 fluidState.setDensity(phaseIdx, value);
196 fluidState.setMolarDensity(phaseIdx, value);
Some exceptions thrown in DuMux
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:59
static void solve(FluidState &fluidState, ParameterCache ¶mCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:81