13#ifndef DUMUX_NCP_FLASH_HH
14#define DUMUX_NCP_FLASH_HH
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
65template <
class Scalar,
class Flu
idSystem>
68 enum { numPhases = FluidSystem::numPhases };
69 enum { numComponents = FluidSystem::numComponents };
71 using ParameterCache =
typename FluidSystem::ParameterCache;
73 static constexpr int numEq = numPhases*(numComponents + 1);
75 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
76 using Vector = Dune::FieldVector<Scalar, numEq>;
86 : wettingPhaseIdx_(wettingPhaseIdx)
95 template <
class Flu
idState>
97 ParameterCache ¶mCache,
102 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
103 sumMoles += globalMolarities[compIdx];
105 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
107 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx)
108 fluidState.setMoleFraction(phaseIdx,
110 globalMolarities[compIdx]/sumMoles);
113 fluidState.setPressure(phaseIdx, 1.0135e5);
116 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
120 paramCache.updateAll(fluidState);
121 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
122 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
123 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
124 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
139 template <
class MaterialLaw,
class Flu
idState>
141 ParameterCache ¶mCache,
142 const MaterialLaw& material,
168 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
170 linearize_(J, b, fluidState, paramCache, material, globalMolarities);
179 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
181 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
192 catch (Dune::FMatrixError &e)
226 Scalar relError =
update_(fluidState, paramCache, material, deltaX);
241 "Flash calculation failed."
242 " {c_alpha^kappa} = {" << globalMolarities <<
"}, T = "
243 << fluidState.temperature(0));
248 template <
class Flu
idState>
251 std::cout <<
"saturations: ";
252 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
253 std::cout << fs.saturation(phaseIdx) <<
" ";
256 std::cout <<
"pressures: ";
257 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
258 std::cout << fs.pressure(phaseIdx) <<
" ";
261 std::cout <<
"densities: ";
262 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
263 std::cout << fs.density(phaseIdx) <<
" ";
266 std::cout <<
"molar densities: ";
267 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
268 std::cout << fs.molarDensity(phaseIdx) <<
" ";
271 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
272 std::cout <<
"composition " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
273 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
274 std::cout << fs.moleFraction(phaseIdx, compIdx) <<
" ";
279 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
280 std::cout <<
"fugacities " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
281 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
282 std::cout << fs.fugacity(phaseIdx, compIdx) <<
" ";
287 std::cout <<
"global component molarities: ";
288 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
290 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
291 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
293 std::cout << sum <<
" ";
298 template <
class MaterialLaw,
class Flu
idState>
301 FluidState &fluidState,
302 ParameterCache ¶mCache,
303 const MaterialLaw& material,
306 FluidState origFluidState(fluidState);
307 ParameterCache origParamCache(paramCache);
317 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
327 setQuantity_(fluidState, paramCache, material, pvIdx, x_i + eps);
335 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
336 J[eqIdx][pvIdx] = tmp[eqIdx];
339 fluidState = origFluidState;
340 paramCache = origParamCache;
347 template <
class Flu
idState>
349 const FluidState &fluidStateEval,
350 const FluidState &fluidState,
356 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
357 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
359 fluidState.fugacity(0, compIdx) -
360 fluidState.fugacity(phaseIdx, compIdx);
365 assert(eqIdx == numComponents*(numPhases - 1));
372 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
374 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
376 fluidState.saturation(phaseIdx)
377 * fluidState.molarity(phaseIdx, compIdx);
380 b[eqIdx] -= globalMolarities[compIdx];
386 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
387 Scalar sumMoleFracEval = 0.0;
388 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
389 sumMoleFracEval += fluidStateEval.moleFraction(phaseIdx, compIdx);
391 if (1.0 - sumMoleFracEval > fluidStateEval.saturation(phaseIdx)) {
392 b[eqIdx] = fluidState.saturation(phaseIdx);
395 Scalar sumMoleFrac = 0.0;
396 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
397 sumMoleFrac += fluidState.moleFraction(phaseIdx, compIdx);
398 b[eqIdx] = 1.0 - sumMoleFrac;
405 template <
class MaterialLaw,
class Flu
idState>
407 ParameterCache ¶mCache,
408 const MaterialLaw& material,
409 const Vector &deltaX)
412 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
414 Scalar delta = deltaX[pvIdx];
418 relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
423 delta = min(0.2, max(-0.2, delta));
428 delta = min(0.15, max(-0.15, delta));
433 delta = min(0.15*fluidState.pressure(0), max(-0.15*fluidState.pressure(0), delta));
441 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
442 Scalar value = fluidState.saturation(phaseIdx);
445 fluidState.setSaturation(phaseIdx, value);
449 value = fluidState.pressure(phaseIdx);
451 fluidState.setPressure(phaseIdx, 0.0);
453 bool allMoleFractionsAreZero =
true;
454 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
455 value = fluidState.moleFraction(phaseIdx, compIdx);
457 allMoleFractionsAreZero =
false;
459 fluidState.setMoleFraction(phaseIdx, compIdx, 0.0);
464 if (allMoleFractionsAreZero)
465 fluidState.setMoleFraction(phaseIdx, 0, 1e-10);
470 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
471 Scalar value = fluidState.saturation(phaseIdx)/(0.95*sumSat);
472 fluidState.setSaturation(phaseIdx, value);
481 template <
class MaterialLaw,
class Flu
idState>
483 ParameterCache ¶mCache,
484 const MaterialLaw& material)
489 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
490 sumSat += fluidState.saturation(phaseIdx);
491 fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
496 const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
497 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
498 fluidState.setPressure(phaseIdx,
499 fluidState.pressure(0)
500 + (pc[phaseIdx] - pc[0]));
503 paramCache.updateAll(fluidState, ParameterCache::Temperature);
506 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
510 fluidState.setDensity(phaseIdx, rho);
511 fluidState.setMolarDensity(phaseIdx, rhoMolar);
513 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
514 Scalar phi = FluidSystem::fugacityCoefficient( fluidState, paramCache, phaseIdx, compIdx);
515 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
521 {
return pvIdx == 0; }
524 {
return 1 <= pvIdx && pvIdx < numPhases; }
527 {
return numPhases <= pvIdx && pvIdx < numPhases + numPhases*numComponents; }
530 template <
class Flu
idState>
533 assert(pvIdx < numEq);
538 return fs.pressure(phaseIdx);
541 else if (pvIdx < numPhases) {
542 int phaseIdx = pvIdx - 1;
543 return fs.saturation(phaseIdx);
548 int phaseIdx = (pvIdx - numPhases)/numComponents;
549 int compIdx = (pvIdx - numPhases)%numComponents;
550 return fs.moleFraction(phaseIdx, compIdx);
555 template <
class MaterialLaw,
class Flu
idState>
557 ParameterCache ¶mCache,
558 const MaterialLaw& material,
562 assert(0 <= pvIdx && pvIdx < numEq);
565 Scalar delta = value - fs.pressure(0);
569 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
570 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
571 paramCache.updateAllPressures(fs);
574 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
577 fs.setDensity(phaseIdx, rho);
578 fs.setMolarDensity(phaseIdx, rhoMolar);
580 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
581 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
582 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
586 else if (pvIdx < numPhases) {
587 Scalar delta = value - fs.saturation(pvIdx - 1);
588 fs.setSaturation(pvIdx - 1, value);
592 fs.setSaturation(numPhases - 1,
593 fs.saturation(numPhases - 1) - delta);
597 const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
598 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
599 fs.setPressure(phaseIdx,
601 + (pc[phaseIdx] - pc[0]));
602 paramCache.updateAllPressures(fs);
605 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
608 fs.setDensity(phaseIdx, rho);
609 fs.setMolarDensity(phaseIdx, rhoMolar);
611 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
612 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
613 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
617 else if (pvIdx < numPhases + numPhases*numComponents)
619 int phaseIdx = (pvIdx - numPhases)/numComponents;
620 int compIdx = (pvIdx - numPhases)%numComponents;
622 fs.setMoleFraction(phaseIdx, compIdx, value);
623 paramCache.updateSingleMoleFraction(fs, phaseIdx, compIdx);
628 fs.setDensity(phaseIdx, rho);
629 fs.setMolarDensity(phaseIdx, rhoMolar);
633 if (!FluidSystem::isIdealMixture(phaseIdx)) {
634 for (
int compIdx2 = 0; compIdx2 < numComponents; ++compIdx2) {
635 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx2);
636 fs.setFugacityCoefficient(phaseIdx, compIdx2, phi);
646 template <
class Flu
idState>
649 assert(pvIdx < numEq);
654 fs.setPressure(phaseIdx, value);
657 else if (pvIdx < numPhases) {
658 int phaseIdx = pvIdx - 1;
659 fs.setSaturation(phaseIdx, value);
664 int phaseIdx = (pvIdx - numPhases)/numComponents;
665 int compIdx = (pvIdx - numPhases)%numComponents;
666 fs.setMoleFraction(phaseIdx, compIdx, value);
670 template <
class Flu
idState>
677 else if (pvIdx < numPhases)
685 int wettingPhaseIdx_ = 0;
Determines the phase compositions, pressures and saturations given the total mass of all components.
Definition: ncpflash.hh:67
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: ncpflash.hh:79
NcpFlash(int wettingPhaseIdx=0)
Construct a new flash.
Definition: ncpflash.hh:85
bool isSaturationIdx_(int pvIdx)
Definition: ncpflash.hh:523
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const Vector &deltaX)
Definition: ncpflash.hh:406
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material)
Definition: ncpflash.hh:482
bool isPressureIdx_(int pvIdx)
Definition: ncpflash.hh:520
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:531
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:299
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: ncpflash.hh:96
void solve(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: ncpflash.hh:140
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: ncpflash.hh:647
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:671
bool isMoleFracIdx_(int pvIdx)
Definition: ncpflash.hh:526
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition: ncpflash.hh:556
void printFluidState_(const FluidState &fs)
Definition: ncpflash.hh:249
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:348
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.
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53