25#ifndef DUMUX_COMPOSITION_FROM_FUGACITIES_HH
26#define DUMUX_COMPOSITION_FROM_FUGACITIES_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
51template <
class Scalar,
class Flu
idSystem>
54 static constexpr int numComponents = FluidSystem::numComponents;
55 using ParameterCache =
typename FluidSystem::ParameterCache;
67 template <
class Flu
idState>
69 ParameterCache ¶mCache,
73 if (!FluidSystem::isIdealMixture(phaseIdx))
76 for (
unsigned int i = 0; i < numComponents; ++ i)
78 fluidState.setMoleFraction(phaseIdx,
95 template <
class Flu
idState>
96 static void solve(FluidState &fluidState,
97 ParameterCache ¶mCache,
103 if (FluidSystem::isIdealMixture(phaseIdx)) {
109 Dune::FieldVector<Scalar, numComponents> xInit;
110 for (
int i = 0; i < numComponents; ++i) {
111 xInit[i] = fluidState.moleFraction(phaseIdx, i);
119 Dune::FieldMatrix<Scalar, numComponents, numComponents> J;
121 Dune::FieldVector<Scalar, numComponents> x;
123 Dune::FieldVector<Scalar, numComponents> b;
125 paramCache.updatePhase(fluidState, phaseIdx);
129 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
131 linearize_(J, b, fluidState, paramCache, phaseIdx, targetFug);
135 try { J.solve(x, b); }
136 catch (Dune::FMatrixError &e)
141 Scalar relError =
update_(fluidState, paramCache, x, b, phaseIdx, targetFug);
143 if (relError < 1e-9) {
146 fluidState.setDensity(phaseIdx, rho);
147 fluidState.setMolarDensity(phaseIdx, rhoMolar);
155 "Calculating the " << FluidSystem::phaseName(phaseIdx)
156 <<
"Phase composition failed. Initial {x} = {"
158 <<
"}, {fug_t} = {" << targetFug <<
"}, p = " << fluidState.pressure(phaseIdx)
159 <<
", T = " << fluidState.temperature(phaseIdx));
167 template <
class Flu
idState>
169 ParameterCache ¶mCache,
173 for (
int i = 0; i < numComponents; ++ i) {
174 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
178 Scalar gamma = phi * fluidState.pressure(phaseIdx);
179 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
180 fluidState.setMoleFraction(phaseIdx, i, fugacities[i]/gamma);
183 paramCache.updatePhase(fluidState, phaseIdx);
186 fluidState.setDensity(phaseIdx, rho);
188 fluidState.setMolarDensity(phaseIdx, rhoMolar);
191 template <
class Flu
idState>
192 static void linearize_(Dune::FieldMatrix<Scalar, numComponents, numComponents> &J,
193 Dune::FieldVector<Scalar, numComponents> &defect,
194 FluidState &fluidState,
195 ParameterCache ¶mCache,
204 for (
int i = 0; i < numComponents; ++ i) {
205 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
209 Scalar f = phi*fluidState.pressure(phaseIdx)*fluidState.moleFraction(phaseIdx, i);
210 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
212 defect[i] = targetFug[i] - f;
216 for (
int i = 0; i < numComponents; ++ i) {
217 const Scalar eps = 1e-11;
226 Scalar x_i = fluidState.moleFraction(phaseIdx, i);
227 fluidState.setMoleFraction(phaseIdx, i, x_i + eps);
228 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
232 for (
int j = 0; j < numComponents; ++j) {
234 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
241 fluidState.pressure(phaseIdx) *
242 fluidState.moleFraction(phaseIdx, j);
244 Scalar defJPlusEps = targetFug[j] - f;
248 J[j][i] = (defJPlusEps - defect[j])/eps;
252 fluidState.setMoleFraction(phaseIdx, i, x_i);
253 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
260 template <
class Flu
idState>
262 ParameterCache ¶mCache,
263 Dune::FieldVector<Scalar, numComponents> &x,
264 Dune::FieldVector<Scalar, numComponents> &b,
266 const Dune::FieldVector<Scalar, numComponents> &targetFug)
269 Dune::FieldVector<Scalar, numComponents> origComp;
275 for (
unsigned int i = 0; i < numComponents; ++i)
277 origComp[i] = fluidState.moleFraction(phaseIdx, i);
278 relError = max(relError, abs(x[i]));
279 sumDelta += abs(x[i]);
283 const Scalar maxDelta = 0.2;
284 if (sumDelta > maxDelta)
285 x /= (sumDelta/maxDelta);
288 for (
unsigned int i = 0; i < numComponents; ++i)
290 Scalar newx = origComp[i] - x[i];
293 if (targetFug[i] > 0)
294 newx = max(0.0, newx);
296 else if (targetFug[i] < 0)
297 newx = min(0.0, newx);
302 fluidState.setMoleFraction(phaseIdx, i, newx);
306 paramCache.updateComposition(fluidState, phaseIdx);
311 template <
class Flu
idState>
318 for (
int i = 0; i < numComponents; ++i) {
322 (targetFug[i] - params.fugacity(phaseIdx, i))
324 params.fugacityCoeff(phaseIdx, i) );
Some exceptions thrown in DuMux
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:53
static void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugVec)
Guess an initial value for the composition of the phase.
Definition: compositionfromfugacities.hh:68
static Scalar calculateDefect_(const FluidState ¶ms, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:312
static Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, Dune::FieldVector< Scalar, numComponents > &x, Dune::FieldVector< Scalar, numComponents > &b, int phaseIdx, const Dune::FieldVector< Scalar, numComponents > &targetFug)
Definition: compositionfromfugacities.hh:261
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: compositionfromfugacities.hh:58
static void linearize_(Dune::FieldMatrix< Scalar, numComponents, numComponents > &J, Dune::FieldVector< Scalar, numComponents > &defect, FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:192
static void solveIdealMix_(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugacities)
Definition: compositionfromfugacities.hh:168
static void solve(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &targetFug)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: compositionfromfugacities.hh:96