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 enum { numComponents = FluidSystem::numComponents };
56 using ParameterCache =
typename FluidSystem::ParameterCache;
68 template <
class Flu
idState>
70 ParameterCache ¶mCache,
74 if (!FluidSystem::isIdealMixture(phaseIdx))
77 for (
unsigned int i = 0; i < numComponents; ++ i)
79 fluidState.setMoleFraction(phaseIdx,
96 template <
class Flu
idState>
97 static void solve(FluidState &fluidState,
98 ParameterCache ¶mCache,
104 if (FluidSystem::isIdealMixture(phaseIdx)) {
110 Dune::FieldVector<Scalar, numComponents> xInit;
111 for (
int i = 0; i < numComponents; ++i) {
112 xInit[i] = fluidState.moleFraction(phaseIdx, i);
120 Dune::FieldMatrix<Scalar, numComponents, numComponents> J;
122 Dune::FieldVector<Scalar, numComponents> x;
124 Dune::FieldVector<Scalar, numComponents> b;
126 paramCache.updatePhase(fluidState, phaseIdx);
130 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
132 linearize_(J, b, fluidState, paramCache, phaseIdx, targetFug);
136 try { J.solve(x, b); }
137 catch (Dune::FMatrixError &e)
142 Scalar relError =
update_(fluidState, paramCache, x, b, phaseIdx, targetFug);
144 if (relError < 1e-9) {
147 fluidState.setDensity(phaseIdx, rho);
148 fluidState.setMolarDensity(phaseIdx, rhoMolar);
156 "Calculating the " << FluidSystem::phaseName(phaseIdx)
157 <<
"Phase composition failed. Initial {x} = {"
159 <<
"}, {fug_t} = {" << targetFug <<
"}, p = " << fluidState.pressure(phaseIdx)
160 <<
", T = " << fluidState.temperature(phaseIdx));
168 template <
class Flu
idState>
170 ParameterCache ¶mCache,
174 for (
int i = 0; i < numComponents; ++ i) {
175 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
179 Scalar gamma = phi * fluidState.pressure(phaseIdx);
180 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
181 fluidState.setMoleFraction(phaseIdx, i, fugacities[i]/gamma);
184 paramCache.updatePhase(fluidState, phaseIdx);
187 fluidState.setDensity(phaseIdx, rho);
189 fluidState.setMolarDensity(phaseIdx, rhoMolar);
192 template <
class Flu
idState>
193 static void linearize_(Dune::FieldMatrix<Scalar, numComponents, numComponents> &J,
194 Dune::FieldVector<Scalar, numComponents> &defect,
195 FluidState &fluidState,
196 ParameterCache ¶mCache,
205 for (
int i = 0; i < numComponents; ++ i) {
206 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
210 Scalar f = phi*fluidState.pressure(phaseIdx)*fluidState.moleFraction(phaseIdx, i);
211 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
213 defect[i] = targetFug[i] - f;
217 for (
int i = 0; i < numComponents; ++ i) {
218 const Scalar eps = 1e-11;
227 Scalar x_i = fluidState.moleFraction(phaseIdx, i);
228 fluidState.setMoleFraction(phaseIdx, i, x_i + eps);
229 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
233 for (
int j = 0; j < numComponents; ++j) {
235 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
242 fluidState.pressure(phaseIdx) *
243 fluidState.moleFraction(phaseIdx, j);
245 Scalar defJPlusEps = targetFug[j] - f;
249 J[j][i] = (defJPlusEps - defect[j])/eps;
253 fluidState.setMoleFraction(phaseIdx, i, x_i);
254 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
261 template <
class Flu
idState>
263 ParameterCache ¶mCache,
264 Dune::FieldVector<Scalar, numComponents> &x,
265 Dune::FieldVector<Scalar, numComponents> &b,
267 const Dune::FieldVector<Scalar, numComponents> &targetFug)
270 Dune::FieldVector<Scalar, numComponents> origComp;
276 for (
unsigned int i = 0; i < numComponents; ++i)
278 origComp[i] = fluidState.moleFraction(phaseIdx, i);
279 relError = max(relError, abs(x[i]));
280 sumDelta += abs(x[i]);
284 const Scalar maxDelta = 0.2;
285 if (sumDelta > maxDelta)
286 x /= (sumDelta/maxDelta);
289 for (
unsigned int i = 0; i < numComponents; ++i)
291 Scalar newx = origComp[i] - x[i];
294 if (targetFug[i] > 0)
295 newx = max(0.0, newx);
297 else if (targetFug[i] < 0)
298 newx = min(0.0, newx);
303 fluidState.setMoleFraction(phaseIdx, i, newx);
307 paramCache.updateComposition(fluidState, phaseIdx);
312 template <
class Flu
idState>
319 for (
int i = 0; i < numComponents; ++i) {
323 (targetFug[i] - params.fugacity(phaseIdx, i))
325 params.fugacityCoeff(phaseIdx, i) );
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
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:69
static Scalar calculateDefect_(const FluidState ¶ms, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:313
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:262
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: compositionfromfugacities.hh:59
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:193
static void solveIdealMix_(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugacities)
Definition: compositionfromfugacities.hh:169
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:97