25#ifndef DUMUX_IMMISCIBLE_FLASH_HH
26#define DUMUX_IMMISCIBLE_FLASH_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
63template <
class Scalar,
class Flu
idSystem>
66 static constexpr int numPhases = FluidSystem::numPhases;
67 static constexpr int numComponents = FluidSystem::numComponents;
68 static_assert(numPhases == numComponents,
69 "Immiscibility assumes that the number of phases"
70 " is equal to the number of components");
72 using ParameterCache =
typename FluidSystem::ParameterCache;
74 static constexpr int numEq = numPhases;
76 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
77 using Vector = Dune::FieldVector<Scalar, numEq>;
87 : wettingPhaseIdx_(wettingPhaseIdx)
96 template <
class Flu
idState>
98 ParameterCache ¶mCache,
103 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
104 sumMoles += globalMolarities[compIdx];
106 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
108 fluidState.setPressure(phaseIdx, 2e5);
111 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
125 template <
class MaterialLaw,
class Flu
idState>
127 ParameterCache& paramCache,
128 const MaterialLaw& material,
134 bool allIncompressible =
true;
135 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 if (FluidSystem::isCompressible(phaseIdx)) {
137 allIncompressible =
false;
142 if (allIncompressible) {
143 paramCache.updateAll(fluidState);
144 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
147 fluidState.setDensity(phaseIdx, rho);
148 fluidState.setMolarDensity(phaseIdx, rhoMolar);
151 globalMolarities[phaseIdx]
152 / fluidState.molarDensity(phaseIdx);
153 fluidState.setSaturation(phaseIdx,
saturation);
171 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
173 linearize_(J, b, fluidState, paramCache, material, globalMolarities);
178 try { J.solve(deltaX, b); }
179 catch (Dune::FMatrixError& e)
211 Scalar relError =
update_(fluidState, paramCache, material, deltaX);
226 "Flash calculation failed."
227 " {c_alpha^kappa} = {" << globalMolarities <<
"}, T = "
228 << fluidState.temperature(0));
233 template <
class Flu
idState>
236 std::cout <<
"saturations: ";
237 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
238 std::cout << fs.saturation(phaseIdx) <<
" ";
241 std::cout <<
"pressures: ";
242 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
243 std::cout << fs.pressure(phaseIdx) <<
" ";
246 std::cout <<
"densities: ";
247 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
248 std::cout << fs.density(phaseIdx) <<
" ";
251 std::cout <<
"global component molarities: ";
252 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
254 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
255 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
257 std::cout << sum <<
" ";
262 template <
class MaterialLaw,
class Flu
idState>
265 FluidState &fluidState,
266 ParameterCache ¶mCache,
267 const MaterialLaw& material,
270 FluidState origFluidState(fluidState);
271 ParameterCache origParamCache(paramCache);
281 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
291 setQuantity_<MaterialLaw>(fluidState, paramCache, material, pvIdx, x_i + eps);
298 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
299 J[eqIdx][pvIdx] = tmp[eqIdx];
302 fluidState = origFluidState;
303 paramCache = origParamCache;
310 template <
class Flu
idState>
312 const FluidState &fluidStateEval,
313 const FluidState &fluidState,
317 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
319 fluidState.saturation(phaseIdx)
320 * fluidState.molarity(phaseIdx, phaseIdx);
322 b[phaseIdx] -= globalMolarities[phaseIdx];
326 template <
class MaterialLaw,
class Flu
idState>
328 ParameterCache ¶mCache,
329 const MaterialLaw& material,
330 const Vector &deltaX)
333 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
335 Scalar delta = deltaX[pvIdx];
339 relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
344 delta = min(0.2, max(-0.2, delta));
349 delta = min(0.30*fluidState.pressure(0), max(-0.30*fluidState.pressure(0), delta));
355 completeFluidState_<MaterialLaw>(fluidState, paramCache, material);
360 template <
class MaterialLaw,
class Flu
idState>
362 ParameterCache ¶mCache,
363 const MaterialLaw& material)
368 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
369 sumSat += fluidState.saturation(phaseIdx);
374 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
376 Scalar S = fluidState.saturation(phaseIdx);
377 fluidState.setSaturation(phaseIdx, S/sumSat);
383 fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
388 const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
389 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
390 fluidState.setPressure(phaseIdx,
391 fluidState.pressure(0)
392 + (pc[phaseIdx] - pc[0]));
395 paramCache.updateAll(fluidState, ParameterCache::Temperature|ParameterCache::Composition);
398 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
401 fluidState.setDensity(phaseIdx, rho);
402 fluidState.setMolarDensity(phaseIdx, rhoMolar);
407 {
return pvIdx == 0; }
410 {
return 1 <= pvIdx && pvIdx < numPhases; }
413 template <
class Flu
idState>
416 assert(pvIdx < numEq);
421 return fs.pressure(phaseIdx);
425 assert(pvIdx < numPhases);
426 int phaseIdx = pvIdx - 1;
427 return fs.saturation(phaseIdx);
432 template <
class MaterialLaw,
class Flu
idState>
434 ParameterCache ¶mCache,
435 const MaterialLaw& material,
439 assert(0 <= pvIdx && pvIdx < numEq);
443 Scalar delta = value - fs.pressure(0);
447 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
448 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
449 paramCache.updateAllPressures(fs);
452 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
455 fs.setDensity(phaseIdx, rho);
456 fs.setMolarDensity(phaseIdx, rhoMolar);
460 assert(pvIdx < numPhases);
463 Scalar delta = value - fs.saturation(pvIdx - 1);
464 fs.setSaturation(pvIdx - 1, value);
468 fs.setSaturation(numPhases - 1,
469 fs.saturation(numPhases - 1) - delta);
473 const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
474 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
475 fs.setPressure(phaseIdx,
477 + (pc[phaseIdx] - pc[0]));
478 paramCache.updateAllPressures(fs);
481 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
484 fs.setDensity(phaseIdx, rho);
485 fs.setMolarDensity(phaseIdx, rhoMolar);
492 template <
class Flu
idState>
495 assert(pvIdx < numEq);
500 fs.setPressure(phaseIdx, value);
504 assert(pvIdx < numPhases);
505 int phaseIdx = pvIdx - 1;
510 value = max(0.0, value);
511 fs.setSaturation(phaseIdx, value);
515 template <
class Flu
idState>
523 assert(pvIdx < numPhases);
529 int wettingPhaseIdx_ = 0;
Some exceptions thrown in DuMux
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
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
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Definition: immiscibleflash.hh:65
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: immiscibleflash.hh:97
void printFluidState_(const FluidState &fs)
Definition: immiscibleflash.hh:234
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:126
bool isPressureIdx_(int pvIdx)
Definition: immiscibleflash.hh:406
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:493
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:311
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:414
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material)
Definition: immiscibleflash.hh:361
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:263
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const Vector &deltaX)
Definition: immiscibleflash.hh:327
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:433
bool isSaturationIdx_(int pvIdx)
Definition: immiscibleflash.hh:409
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:516
ImmiscibleFlash(int wettingPhaseIdx=0)
Contruct a new flash.
Definition: immiscibleflash.hh:86
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: immiscibleflash.hh:80