25#ifndef DUMUX_IMMISCIBLE_FLASH_HH
26#define DUMUX_IMMISCIBLE_FLASH_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30#include <dune/common/version.hh>
64template <
class Scalar,
class Flu
idSystem>
67 static constexpr int numPhases = FluidSystem::numPhases;
68 static constexpr int numComponents = FluidSystem::numComponents;
69 static_assert(numPhases == numComponents,
70 "Immiscibility assumes that the number of phases"
71 " is equal to the number of components");
73 using ParameterCache =
typename FluidSystem::ParameterCache;
75 static constexpr int numEq = numPhases;
77 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
78 using Vector = Dune::FieldVector<Scalar, numEq>;
88 : wettingPhaseIdx_(wettingPhaseIdx)
97 template <
class Flu
idState>
99 ParameterCache ¶mCache,
104 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
105 sumMoles += globalMolarities[compIdx];
107 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
109 fluidState.setPressure(phaseIdx, 2e5);
112 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
126 template <
class MaterialLaw,
class Flu
idState>
128 ParameterCache ¶mCache,
129 const typename MaterialLaw::Params &matParams,
132#if DUNE_VERSION_LT(DUNE_COMMON, 2, 7)
133 Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-25);
139 bool allIncompressible =
true;
140 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
141 if (FluidSystem::isCompressible(phaseIdx)) {
142 allIncompressible =
false;
147 if (allIncompressible) {
148 paramCache.updateAll(fluidState);
149 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
152 fluidState.setDensity(phaseIdx, rho);
153 fluidState.setMolarDensity(phaseIdx, rhoMolar);
156 globalMolarities[phaseIdx]
157 / fluidState.molarDensity(phaseIdx);
158 fluidState.setSaturation(phaseIdx,
saturation);
177 completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
180 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
182 linearize_<MaterialLaw>(J, b, fluidState, paramCache, matParams, globalMolarities);
189 try { J.solve(deltaX, b); }
190 catch (Dune::FMatrixError& e)
223 Scalar relError = update_<MaterialLaw>(fluidState, paramCache, matParams, deltaX);
238 "Flash calculation failed."
239 " {c_alpha^kappa} = {" << globalMolarities <<
"}, T = "
240 << fluidState.temperature(0));
245 template <
class Flu
idState>
248 std::cout <<
"saturations: ";
249 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
250 std::cout << fs.saturation(phaseIdx) <<
" ";
253 std::cout <<
"pressures: ";
254 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
255 std::cout << fs.pressure(phaseIdx) <<
" ";
258 std::cout <<
"densities: ";
259 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
260 std::cout << fs.density(phaseIdx) <<
" ";
263 std::cout <<
"global component molarities: ";
264 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
266 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
267 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
269 std::cout << sum <<
" ";
274 template <
class MaterialLaw,
class Flu
idState>
277 FluidState &fluidState,
278 ParameterCache ¶mCache,
279 const typename MaterialLaw::Params &matParams,
282 FluidState origFluidState(fluidState);
283 ParameterCache origParamCache(paramCache);
295 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
305 setQuantity_<MaterialLaw>(fluidState, paramCache, matParams, pvIdx, x_i +
eps);
312 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
313 J[eqIdx][pvIdx] = tmp[eqIdx];
316 fluidState = origFluidState;
317 paramCache = origParamCache;
324 template <
class Flu
idState>
326 const FluidState &fluidStateEval,
327 const FluidState &fluidState,
331 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
333 fluidState.saturation(phaseIdx)
334 * fluidState.molarity(phaseIdx, phaseIdx);
336 b[phaseIdx] -= globalMolarities[phaseIdx];
340 template <
class MaterialLaw,
class Flu
idState>
342 ParameterCache ¶mCache,
343 const typename MaterialLaw::Params &matParams,
344 const Vector &deltaX)
347 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
349 Scalar delta = deltaX[pvIdx];
353 relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
358 delta = min(0.2, max(-0.2, delta));
363 delta = min(0.30*fluidState.pressure(0), max(-0.30*fluidState.pressure(0), delta));
369 completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
374 template <
class MaterialLaw,
class Flu
idState>
376 ParameterCache ¶mCache,
377 const typename MaterialLaw::Params &matParams)
382 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
383 sumSat += fluidState.saturation(phaseIdx);
388 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
390 Scalar S = fluidState.saturation(phaseIdx);
391 fluidState.setSaturation(phaseIdx, S/sumSat);
397 fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
403 MaterialLaw::capillaryPressures(pc, matParams, fluidState, wettingPhaseIdx_);
404 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
405 fluidState.setPressure(phaseIdx,
406 fluidState.pressure(0)
407 + (pc[phaseIdx] - pc[0]));
410 paramCache.updateAll(fluidState, ParameterCache::Temperature|ParameterCache::Composition);
413 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
416 fluidState.setDensity(phaseIdx, rho);
417 fluidState.setMolarDensity(phaseIdx, rhoMolar);
422 {
return pvIdx == 0; }
425 {
return 1 <= pvIdx && pvIdx < numPhases; }
428 template <
class Flu
idState>
431 assert(pvIdx < numEq);
436 return fs.pressure(phaseIdx);
440 assert(pvIdx < numPhases);
441 int phaseIdx = pvIdx - 1;
442 return fs.saturation(phaseIdx);
447 template <
class MaterialLaw,
class Flu
idState>
449 ParameterCache ¶mCache,
450 const typename MaterialLaw::Params &matParams,
454 assert(0 <= pvIdx && pvIdx < numEq);
458 Scalar delta = value - fs.pressure(0);
462 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
463 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
464 paramCache.updateAllPressures(fs);
467 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
470 fs.setDensity(phaseIdx, rho);
471 fs.setMolarDensity(phaseIdx, rhoMolar);
475 assert(pvIdx < numPhases);
478 Scalar delta = value - fs.saturation(pvIdx - 1);
479 fs.setSaturation(pvIdx - 1, value);
483 fs.setSaturation(numPhases - 1,
484 fs.saturation(numPhases - 1) - delta);
489 MaterialLaw::capillaryPressures(pc, matParams, fs, wettingPhaseIdx_);
490 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
491 fs.setPressure(phaseIdx,
493 + (pc[phaseIdx] - pc[0]));
494 paramCache.updateAllPressures(fs);
497 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
500 fs.setDensity(phaseIdx, rho);
501 fs.setMolarDensity(phaseIdx, rhoMolar);
508 template <
class Flu
idState>
511 assert(pvIdx < numEq);
516 fs.setPressure(phaseIdx, value);
520 assert(pvIdx < numPhases);
521 int phaseIdx = pvIdx - 1;
526 value = max(0.0, value);
527 fs.setSaturation(phaseIdx, value);
531 template <
class Flu
idState>
539 assert(pvIdx < numPhases);
545 int wettingPhaseIdx_ = 0;
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
void SetUndefined(const T &value)
Make the memory on which an object resides undefined.
Definition: valgrind.hh:102
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
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
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Definition: immiscibleflash.hh:66
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: immiscibleflash.hh:98
void printFluidState_(const FluidState &fs)
Definition: immiscibleflash.hh:246
bool isPressureIdx_(int pvIdx)
Definition: immiscibleflash.hh:421
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:509
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:325
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:429
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, const Vector &deltaX)
Definition: immiscibleflash.hh:341
void solve(FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: immiscibleflash.hh:127
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:275
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:448
bool isSaturationIdx_(int pvIdx)
Definition: immiscibleflash.hh:424
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams)
Definition: immiscibleflash.hh:375
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:532
ImmiscibleFlash(int wettingPhaseIdx=0)
Contruct a new flash.
Definition: immiscibleflash.hh:87
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: immiscibleflash.hh:81