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>
30#include <dune/common/version.hh>
53template <
class Scalar,
class Flu
idSystem>
56 enum { numComponents = FluidSystem::numComponents };
58 using ParameterCache =
typename FluidSystem::ParameterCache;
70 template <
class Flu
idState>
72 ParameterCache ¶mCache,
76 if (!FluidSystem::isIdealMixture(phaseIdx))
79 for (
unsigned int i = 0; i < numComponents; ++ i)
81 fluidState.setMoleFraction(phaseIdx,
98 template <
class Flu
idState>
99 static void solve(FluidState &fluidState,
100 ParameterCache ¶mCache,
106 if (FluidSystem::isIdealMixture(phaseIdx)) {
110#if DUNE_VERSION_LT(DUNE_COMMON, 2, 7)
111 Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-25);
114 Dune::FieldVector<Scalar, numComponents> xInit;
115 for (
int i = 0; i < numComponents; ++i) {
116 xInit[i] = fluidState.moleFraction(phaseIdx, i);
124 Dune::FieldMatrix<Scalar, numComponents, numComponents> J;
126 Dune::FieldVector<Scalar, numComponents> x;
128 Dune::FieldVector<Scalar, numComponents> b;
130 paramCache.updatePhase(fluidState, phaseIdx);
134 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
136 linearize_(J, b, fluidState, paramCache, phaseIdx, targetFug);
142 try { J.solve(x, b); }
143 catch (Dune::FMatrixError &e)
152 Scalar relError =
update_(fluidState, paramCache, x, b, phaseIdx, targetFug);
155 if (relError < 1e-9) {
158 fluidState.setDensity(phaseIdx, rho);
159 fluidState.setMolarDensity(phaseIdx, rhoMolar);
167 "Calculating the " << FluidSystem::phaseName(phaseIdx)
168 <<
"Phase composition failed. Initial {x} = {"
170 <<
"}, {fug_t} = {" << targetFug <<
"}, p = " << fluidState.pressure(phaseIdx)
171 <<
", T = " << fluidState.temperature(phaseIdx));
179 template <
class Flu
idState>
181 ParameterCache ¶mCache,
185 for (
int i = 0; i < numComponents; ++ i) {
186 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
190 Scalar gamma = phi * fluidState.pressure(phaseIdx);
191 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
192 fluidState.setMoleFraction(phaseIdx, i, fugacities[i]/gamma);
195 paramCache.updatePhase(fluidState, phaseIdx);
198 fluidState.setDensity(phaseIdx, rho);
200 fluidState.setMolarDensity(phaseIdx, rhoMolar);
203 template <
class Flu
idState>
204 static void linearize_(Dune::FieldMatrix<Scalar, numComponents, numComponents> &J,
205 Dune::FieldVector<Scalar, numComponents> &defect,
206 FluidState &fluidState,
207 ParameterCache ¶mCache,
216 for (
int i = 0; i < numComponents; ++ i) {
217 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
221 Scalar f = phi*fluidState.pressure(phaseIdx)*fluidState.moleFraction(phaseIdx, i);
222 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
224 defect[i] = targetFug[i] - f;
228 for (
int i = 0; i < numComponents; ++ i) {
229 const Scalar
eps = 1e-11;
238 Scalar x_i = fluidState.moleFraction(phaseIdx, i);
239 fluidState.setMoleFraction(phaseIdx, i, x_i +
eps);
240 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
244 for (
int j = 0; j < numComponents; ++j) {
246 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
253 fluidState.pressure(phaseIdx) *
254 fluidState.moleFraction(phaseIdx, j);
256 Scalar defJPlusEps = targetFug[j] - f;
260 J[j][i] = (defJPlusEps - defect[j])/
eps;
264 fluidState.setMoleFraction(phaseIdx, i, x_i);
265 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
272 template <
class Flu
idState>
274 ParameterCache ¶mCache,
275 Dune::FieldVector<Scalar, numComponents> &x,
276 Dune::FieldVector<Scalar, numComponents> &b,
278 const Dune::FieldVector<Scalar, numComponents> &targetFug)
281 Dune::FieldVector<Scalar, numComponents> origComp;
287 for (
unsigned int i = 0; i < numComponents; ++i)
289 origComp[i] = fluidState.moleFraction(phaseIdx, i);
290 relError = max(relError, abs(x[i]));
291 sumDelta += abs(x[i]);
295 const Scalar maxDelta = 0.2;
296 if (sumDelta > maxDelta)
297 x /= (sumDelta/maxDelta);
300 for (
unsigned int i = 0; i < numComponents; ++i)
302 Scalar newx = origComp[i] - x[i];
305 if (targetFug[i] > 0)
306 newx = max(0.0, newx);
308 else if (targetFug[i] < 0)
309 newx = min(0.0, newx);
314 fluidState.setMoleFraction(phaseIdx, i, newx);
318 paramCache.updateComposition(fluidState, phaseIdx);
323 template <
class Flu
idState>
330 for (
int i = 0; i < numComponents; ++i) {
334 (targetFug[i] - params.fugacity(phaseIdx, i))
336 params.fugacityCoeff(phaseIdx, i) );
Some exceptions thrown in DuMux
Some templates to wrap the valgrind macros.
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
make the local view function available whenever we use the grid geometry
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
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
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:55
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:71
static Scalar calculateDefect_(const FluidState ¶ms, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:324
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:273
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: compositionfromfugacities.hh:61
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:204
static void solveIdealMix_(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugacities)
Definition: compositionfromfugacities.hh:180
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:99