13#ifndef DUMUX_IMMISCIBLE_FLASH_HH 
   14#define DUMUX_IMMISCIBLE_FLASH_HH 
   16#include <dune/common/fvector.hh> 
   17#include <dune/common/fmatrix.hh> 
   51template <
class Scalar, 
class Flu
idSystem>
 
   54    static constexpr int numPhases = FluidSystem::numPhases;
 
   55    static constexpr int numComponents = FluidSystem::numComponents;
 
   56    static_assert(numPhases == numComponents,
 
   57                  "Immiscibility assumes that the number of phases" 
   58                  " is equal to the number of components");
 
   60    using ParameterCache = 
typename FluidSystem::ParameterCache;
 
   62    static constexpr int numEq = numPhases;
 
   64    using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
 
   65    using Vector = Dune::FieldVector<Scalar, numEq>;
 
   75    : wettingPhaseIdx_(wettingPhaseIdx)
 
 
   84    template <
class Flu
idState>
 
   86                      ParameterCache ¶mCache,
 
   91        for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
 
   92            sumMoles += globalMolarities[compIdx];
 
   94        for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
 
   96            fluidState.setPressure(phaseIdx, 2e5);
 
   99            fluidState.setSaturation(phaseIdx, 1.0/numPhases);
 
 
  113    template <
class MaterialLaw, 
class Flu
idState>
 
  115               ParameterCache& paramCache,
 
  116               const MaterialLaw& material,
 
  122        bool allIncompressible = 
true;
 
  123        for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  124            if (FluidSystem::isCompressible(phaseIdx)) {
 
  125                allIncompressible = 
false;
 
  130        if (allIncompressible) {
 
  131            paramCache.updateAll(fluidState);
 
  132            for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  133                Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
 
  134                Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
 
  135                fluidState.setDensity(phaseIdx, rho);
 
  136                fluidState.setMolarDensity(phaseIdx, rhoMolar);
 
  139                    globalMolarities[phaseIdx]
 
  140                    / fluidState.molarDensity(phaseIdx);
 
  141                fluidState.setSaturation(phaseIdx, saturation);
 
  159        for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
 
  161            linearize_(J, b, fluidState, paramCache, material, globalMolarities);
 
  166            try { J.solve(deltaX, b); }
 
  167            catch (Dune::FMatrixError& e)
 
  199            Scalar relError = 
update_(fluidState, paramCache, material, deltaX);
 
  214                   "Flash calculation failed." 
  215                   " {c_alpha^kappa} = {" << globalMolarities << 
"}, T = " 
  216                   << fluidState.temperature(0));
 
 
  221    template <
class Flu
idState>
 
  224        std::cout << 
"saturations: ";
 
  225        for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  226            std::cout << fs.saturation(phaseIdx) << 
" ";
 
  229        std::cout << 
"pressures: ";
 
  230        for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  231            std::cout << fs.pressure(phaseIdx) << 
" ";
 
  234        std::cout << 
"densities: ";
 
  235        for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  236            std::cout << fs.density(phaseIdx) << 
" ";
 
  239        std::cout << 
"global component molarities: ";
 
  240        for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  242            for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  243                sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
 
  245            std::cout << sum << 
" ";
 
 
  250    template <
class MaterialLaw, 
class Flu
idState>
 
  253                    FluidState &fluidState,
 
  254                    ParameterCache ¶mCache,
 
  255                    const MaterialLaw& material,
 
  258        FluidState origFluidState(fluidState);
 
  259        ParameterCache origParamCache(paramCache);
 
  269        for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
 
  286            for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
 
  287                J[eqIdx][pvIdx] = tmp[eqIdx];
 
  290            fluidState = origFluidState;
 
  291            paramCache = origParamCache;
 
 
  298    template <
class Flu
idState>
 
  300                          const FluidState &fluidStateEval,
 
  301                          const FluidState &fluidState,
 
  305        for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  307                fluidState.saturation(phaseIdx)
 
  308                * fluidState.molarity(phaseIdx, phaseIdx);
 
  310            b[phaseIdx] -= globalMolarities[phaseIdx];
 
 
  314    template <
class MaterialLaw, 
class Flu
idState>
 
  316                   ParameterCache ¶mCache,
 
  317                   const MaterialLaw& material,
 
  318                   const Vector &deltaX)
 
  321        for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
 
  323            Scalar delta = deltaX[pvIdx];
 
  327            relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
 
  332                delta = min(0.2, max(-0.2, delta));
 
  337                delta = min(0.30*fluidState.pressure(0), max(-0.30*fluidState.pressure(0), delta));
 
 
  348    template <
class MaterialLaw, 
class Flu
idState>
 
  350                             ParameterCache ¶mCache,
 
  351                             const MaterialLaw& material)
 
  356        for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
 
  357            sumSat += fluidState.saturation(phaseIdx);
 
  362            for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
 
  364                Scalar S = fluidState.saturation(phaseIdx);
 
  365                fluidState.setSaturation(phaseIdx, S/sumSat);
 
  371        fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
 
  376        const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
 
  377        for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
 
  378            fluidState.setPressure(phaseIdx,
 
  379                                   fluidState.pressure(0)
 
  380                                   + (pc[phaseIdx] - pc[0]));
 
  383        paramCache.updateAll(fluidState, ParameterCache::Temperature|ParameterCache::Composition);
 
  386        for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
 
  387            Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
 
  388            Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
 
  389            fluidState.setDensity(phaseIdx, rho);
 
  390            fluidState.setMolarDensity(phaseIdx, rhoMolar);
 
 
  395    { 
return pvIdx == 0; }
 
 
  398    { 
return 1 <= pvIdx && pvIdx < numPhases; }
 
 
  401    template <
class Flu
idState>
 
  404        assert(pvIdx < numEq);
 
  409            return fs.pressure(phaseIdx);
 
  413            assert(pvIdx < numPhases);
 
  414            int phaseIdx = pvIdx - 1;
 
  415            return fs.saturation(phaseIdx);
 
 
  420    template <
class MaterialLaw, 
class Flu
idState>
 
  422                      ParameterCache ¶mCache,
 
  423                      const MaterialLaw& material,
 
  427        assert(0 <= pvIdx && pvIdx < numEq);
 
  431            Scalar delta = value - fs.pressure(0);
 
  435            for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  436                fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
 
  437            paramCache.updateAllPressures(fs);
 
  440            for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  441                Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
 
  442                Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
 
  443                fs.setDensity(phaseIdx, rho);
 
  444                fs.setMolarDensity(phaseIdx, rhoMolar);
 
  448            assert(pvIdx < numPhases);
 
  451            Scalar delta = value - fs.saturation(pvIdx - 1);
 
  452            fs.setSaturation(pvIdx - 1, value);
 
  456            fs.setSaturation(numPhases - 1,
 
  457                             fs.saturation(numPhases - 1) - delta);
 
  461            const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
 
  462            for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
 
  463                fs.setPressure(phaseIdx,
 
  465                               + (pc[phaseIdx] - pc[0]));
 
  466            paramCache.updateAllPressures(fs);
 
  469            for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  470                Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
 
  471                Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
 
  472                fs.setDensity(phaseIdx, rho);
 
  473                fs.setMolarDensity(phaseIdx, rhoMolar);
 
 
  480    template <
class Flu
idState>
 
  483        assert(pvIdx < numEq);
 
  488            fs.setPressure(phaseIdx, value);
 
  492            assert(pvIdx < numPhases);
 
  493            int phaseIdx = pvIdx - 1;
 
  498            value = max(0.0, value);
 
  499            fs.setSaturation(phaseIdx, value);
 
 
  503    template <
class Flu
idState>
 
  511            assert(pvIdx < numPhases);
 
 
  517    int wettingPhaseIdx_ = 0; 
 
 
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition immiscibleflash.hh:85
void printFluidState_(const FluidState &fs)
Definition immiscibleflash.hh:222
void solve(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition immiscibleflash.hh:114
bool isPressureIdx_(int pvIdx)
Definition immiscibleflash.hh:394
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition immiscibleflash.hh:481
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition immiscibleflash.hh:299
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition immiscibleflash.hh:402
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material)
Definition immiscibleflash.hh:349
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition immiscibleflash.hh:251
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const Vector &deltaX)
Definition immiscibleflash.hh:315
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition immiscibleflash.hh:421
bool isSaturationIdx_(int pvIdx)
Definition immiscibleflash.hh:397
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition immiscibleflash.hh:504
ImmiscibleFlash(int wettingPhaseIdx=0)
Construct a new flash.
Definition immiscibleflash.hh:74
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition immiscibleflash.hh:68
Exception thrown if a fixable numerical problem occurs.
Definition exceptions.hh:27
Some exceptions thrown in DuMux
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.