13#ifndef DUMUX_COMPOSITION_FROM_FUGACITIES_HH 
   14#define DUMUX_COMPOSITION_FROM_FUGACITIES_HH 
   16#include <dune/common/fvector.hh> 
   17#include <dune/common/fmatrix.hh> 
   39template <
class Scalar, 
class Flu
idSystem>
 
   42    static constexpr int numComponents = FluidSystem::numComponents;
 
   43    using ParameterCache = 
typename FluidSystem::ParameterCache;
 
   55    template <
class Flu
idState>
 
   57                             ParameterCache ¶mCache,
 
   61        if (!FluidSystem::isIdealMixture(phaseIdx))
 
   64            for (
unsigned int i = 0; i < numComponents; ++ i)
 
   66                fluidState.setMoleFraction(phaseIdx,
 
 
   83    template <
class Flu
idState>
 
   84    static void solve(FluidState &fluidState,
 
   85                      ParameterCache ¶mCache,
 
   91        if (FluidSystem::isIdealMixture(phaseIdx)) {
 
   97        Dune::FieldVector<Scalar, numComponents> xInit;
 
   98        for (
int i = 0; i < numComponents; ++i) {
 
   99            xInit[i] = fluidState.moleFraction(phaseIdx, i);
 
  107        Dune::FieldMatrix<Scalar, numComponents, numComponents> J;
 
  109        Dune::FieldVector<Scalar, numComponents> x;
 
  111        Dune::FieldVector<Scalar, numComponents> b;
 
  113        paramCache.updatePhase(fluidState, phaseIdx);
 
  117        for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
 
  119            linearize_(J, b, fluidState, paramCache, phaseIdx, targetFug);
 
  123            try { J.solve(x, b); }
 
  124            catch (Dune::FMatrixError &e)
 
  129            Scalar relError = 
update_(fluidState, paramCache, x, b, phaseIdx, targetFug);
 
  131            if (relError < 1e-9) {
 
  132                Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
 
  133                Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
 
  134                fluidState.setDensity(phaseIdx, rho);
 
  135                fluidState.setMolarDensity(phaseIdx, rhoMolar);
 
  143                   "Calculating the " << FluidSystem::phaseName(phaseIdx)
 
  144                   << 
"Phase composition failed. Initial {x} = {" 
  146                   << 
"}, {fug_t} = {" << targetFug << 
"}, p = " << fluidState.pressure(phaseIdx)
 
  147                   << 
", T = " << fluidState.temperature(phaseIdx));
 
 
  155    template <
class Flu
idState>
 
  157                               ParameterCache ¶mCache,
 
  161        for (
int i = 0; i < numComponents; ++ i) {
 
  162            Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
 
  166            Scalar gamma = phi * fluidState.pressure(phaseIdx);
 
  167            fluidState.setFugacityCoefficient(phaseIdx, i, phi);
 
  168            fluidState.setMoleFraction(phaseIdx, i, fugacities[i]/gamma);
 
  171        paramCache.updatePhase(fluidState, phaseIdx);
 
  173        Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
 
  174        fluidState.setDensity(phaseIdx, rho);
 
  175        Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
 
  176        fluidState.setMolarDensity(phaseIdx, rhoMolar);
 
 
  179    template <
class Flu
idState>
 
  180    static void linearize_(Dune::FieldMatrix<Scalar, numComponents, numComponents> &J,
 
  181                           Dune::FieldVector<Scalar, numComponents> &defect,
 
  182                           FluidState &fluidState,
 
  183                           ParameterCache ¶mCache,
 
  192        for (
int i = 0; i < numComponents; ++ i) {
 
  193            Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
 
  197            Scalar f = phi*fluidState.pressure(phaseIdx)*fluidState.moleFraction(phaseIdx, i);
 
  198            fluidState.setFugacityCoefficient(phaseIdx, i, phi);
 
  200            defect[i] = targetFug[i] - f;
 
  204        for (
int i = 0; i < numComponents; ++ i) {
 
  205            const Scalar eps = 1e-11; 
 
  214            Scalar x_i = fluidState.moleFraction(phaseIdx, i);
 
  215            fluidState.setMoleFraction(phaseIdx, i, x_i + eps);
 
  216            paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
 
  220            for (
int j = 0; j < numComponents; ++j) {
 
  222                Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
 
  229                    fluidState.pressure(phaseIdx) *
 
  230                    fluidState.moleFraction(phaseIdx, j);
 
  232                Scalar defJPlusEps = targetFug[j] - f;
 
  236                J[j][i] = (defJPlusEps - defect[j])/eps;
 
  240            fluidState.setMoleFraction(phaseIdx, i, x_i);
 
  241            paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
 
 
  248    template <
class Flu
idState>
 
  250                          ParameterCache ¶mCache,
 
  251                          Dune::FieldVector<Scalar, numComponents> &x,
 
  252                          Dune::FieldVector<Scalar, numComponents> &b,
 
  254                          const Dune::FieldVector<Scalar, numComponents> &targetFug)
 
  257        Dune::FieldVector<Scalar, numComponents> origComp;
 
  263        for (
unsigned int i = 0; i < numComponents; ++i)
 
  265            origComp[i] = fluidState.moleFraction(phaseIdx, i);
 
  266            relError = max(relError, abs(x[i]));
 
  267            sumDelta += abs(x[i]);
 
  271        const Scalar maxDelta = 0.2;
 
  272        if (sumDelta > maxDelta)
 
  273            x /= (sumDelta/maxDelta);
 
  276        for (
unsigned int i = 0; i < numComponents; ++i)
 
  278            Scalar newx = origComp[i] - x[i];
 
  281            if (targetFug[i] > 0)
 
  282                newx = max(0.0, newx);
 
  284            else if (targetFug[i] < 0)
 
  285                newx = min(0.0, newx);
 
  290            fluidState.setMoleFraction(phaseIdx, i, newx);
 
  294        paramCache.updateComposition(fluidState, phaseIdx);
 
 
  299    template <
class Flu
idState>
 
  306        for (
int i = 0; i < numComponents; ++i) {
 
  310                (targetFug[i] - params.fugacity(phaseIdx, i))
 
  312                params.fugacityCoeff(phaseIdx, i) );
 
 
 
Calculates the chemical equilibrium from the component fugacities  in the phase .
Definition compositionfromfugacities.hh:41
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:56
static Scalar calculateDefect_(const FluidState ¶ms, int phaseIdx, const ComponentVector &targetFug)
Definition compositionfromfugacities.hh:300
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:249
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition compositionfromfugacities.hh:46
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:180
static void solveIdealMix_(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugacities)
Definition compositionfromfugacities.hh:156
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:84
Exception thrown if a fixable numerical problem occurs.
Definition exceptions.hh:27
Some exceptions thrown in DuMux