version 3.8
misciblemultiphasecomposition.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_MISCIBLE_MULTIPHASE_COMPOSITION_HH
15#define DUMUX_MISCIBLE_MULTIPHASE_COMPOSITION_HH
16
17#include <dune/common/fvector.hh>
18#include <dune/common/fmatrix.hh>
19
21
22namespace Dumux {
45template <class Scalar, class FluidSystem>
47{
48 static constexpr int numPhases = FluidSystem::numPhases;
49 static constexpr int numComponents = FluidSystem::numComponents;
50 static const int numMajorComponents = FluidSystem::numPhases;
51
52public:
68 template <class FluidState, class ParameterCache>
69 static void solve(FluidState &fluidState,
70 ParameterCache &paramCache,
71 int knownPhaseIdx = 0)
72 {
73#ifndef NDEBUG
74 // currently this solver can only handle fluid systems which
75 // assume ideal mixtures of all fluids. TODO: relax this
76 // (requires solving a non-linear system of equations, i.e. using
77 // Newton method.)
78 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
79 assert(FluidSystem::isIdealMixture(phaseIdx));
80
81 }
82#endif
83
84 //get the known mole fractions from the fluidState
85 //in a 2pnc system the n>2 mole fractions are primary variables and are already set in the fluidstate
86 Dune::FieldVector<Scalar, numComponents-numMajorComponents> xKnown(0.0);
87 for (int knownCompIdx = 0; knownCompIdx < numComponents-numMajorComponents; ++knownCompIdx)
88 {
89 xKnown[knownCompIdx] = fluidState.moleFraction(knownPhaseIdx, knownCompIdx + numMajorComponents);
90 }
91
92 // compute all fugacity coefficients
93 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
94 paramCache.updatePhase(fluidState, phaseIdx);
95
96 // since we assume ideal mixtures, the fugacity
97 // coefficients of the components cannot depend on
98 // composition, i.e. the parameters in the cache are valid
99 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
100 Scalar fugCoeff = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
101 fluidState.setFugacityCoefficient(phaseIdx, compIdx, fugCoeff);
102 }
103 }
104
105
106 // create the linear system of equations which defines the
107 // mole fractions
108 Dune::FieldMatrix<Scalar, numComponents*numPhases, numComponents*numPhases> M(0.0);
109 Dune::FieldVector<Scalar, numComponents*numPhases> x(0.0);
110 Dune::FieldVector<Scalar, numComponents*numPhases> b(0.0);
111
112 // assemble the equations expressing the assumption that the
113 // sum of all mole fractions in each phase must be 1
114 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
115 int rowIdx = numComponents*(numPhases - 1) + phaseIdx;
116 b[rowIdx] = 1.0;
117
118 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
119 int colIdx = phaseIdx*numComponents + compIdx;
120
121 M[rowIdx][colIdx] = 1.0;
122 }
123 }
124
125 // set the additional equations for the numComponents-numMajorComponents
126 // this is only relevant if numphases != numcomponents e.g. in a 2pnc system
127 // Components, of which the mole fractions are known, set to molefraction(knownCompIdx)=xKnown
128 for(int knownCompIdx = 0; knownCompIdx < numComponents-numMajorComponents; ++knownCompIdx)
129 {
130 int rowIdx = numComponents + numPhases + knownCompIdx;
131 int colIdx = knownPhaseIdx*numComponents + knownCompIdx + numMajorComponents;
132 M[rowIdx][colIdx] = 1.0;
133 b[rowIdx] = xKnown[knownCompIdx];
134 }
135
136 // assemble the equations expressing the fact that the
137 // fugacities of each component are equal in all phases
138 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
139 {
140 int col1Idx = compIdx;
141 const auto entryPhase0 = fluidState.fugacityCoefficient(0, compIdx)*fluidState.pressure(0);
142
143 for (unsigned int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
144 {
145 int rowIdx = (phaseIdx - 1)*numComponents + compIdx;
146 int col2Idx = phaseIdx*numComponents + compIdx;
147 M[rowIdx][col1Idx] = entryPhase0;
148 M[rowIdx][col2Idx] = -fluidState.fugacityCoefficient(phaseIdx, compIdx)*fluidState.pressure(phaseIdx);
149 }
150 }
151
152 // preconditioning of M to reduce condition number
153 for (int compIdx = 0; compIdx < numComponents; compIdx++)
154 {
155 // Multiply row of main component (Raoult's Law) with 10e-5 (order of magn. of pressure)
156 if (compIdx < numMajorComponents)
157 M[compIdx] *= 10e-5;
158
159 // Multiply row of sec. components (Henry's Law) with 10e-9 (order of magn. of Henry constant)
160 else
161 M[compIdx] *= 10e-9;
162 }
163
164 // solve for all mole fractions
165 try { M.solve(x, b); }
166 catch (const Dune::FMatrixError& e) {
167 DUNE_THROW(NumericalProblem,
168 "MiscibleMultiPhaseComposition: Failed to solve the linear equation system for the phase composition.");
169 }
170
171 // set all mole fractions and the the additional quantities in
172 // the fluid state
173 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
174 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
175 int rowIdx = phaseIdx*numComponents + compIdx;
176 fluidState.setMoleFraction(phaseIdx, compIdx, x[rowIdx]);
177 }
178 paramCache.updateComposition(fluidState, phaseIdx);
179
180 Scalar value = FluidSystem::density(fluidState, paramCache, phaseIdx);
181 fluidState.setDensity(phaseIdx, value);
182
183 value = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
184 fluidState.setMolarDensity(phaseIdx, value);
185 }
186 }
187};
188
189} // end namespace Dumux
190
191#endif
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:47
static void solve(FluidState &fluidState, ParameterCache &paramCache, 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:69
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
Some exceptions thrown in DuMux
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17