25#ifndef DUMUX_NCP_FLASH_HH
26#define DUMUX_NCP_FLASH_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
77template <
class Scalar,
class Flu
idSystem>
80 enum { numPhases = FluidSystem::numPhases };
81 enum { numComponents = FluidSystem::numComponents };
83 using ParameterCache =
typename FluidSystem::ParameterCache;
85 static constexpr int numEq = numPhases*(numComponents + 1);
87 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
88 using Vector = Dune::FieldVector<Scalar, numEq>;
98 : wettingPhaseIdx_(wettingPhaseIdx)
107 template <
class Flu
idState>
109 ParameterCache ¶mCache,
114 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
115 sumMoles += globalMolarities[compIdx];
117 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
119 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx)
120 fluidState.setMoleFraction(phaseIdx,
122 globalMolarities[compIdx]/sumMoles);
125 fluidState.setPressure(phaseIdx, 1.0135e5);
128 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
132 paramCache.updateAll(fluidState);
133 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
134 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
135 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
136 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
151 template <
class MaterialLaw,
class Flu
idState>
153 ParameterCache ¶mCache,
154 const MaterialLaw& material,
180 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
182 linearize_(J, b, fluidState, paramCache, material, globalMolarities);
191 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
193 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
204 catch (Dune::FMatrixError &e)
238 Scalar relError =
update_(fluidState, paramCache, material, deltaX);
253 "Flash calculation failed."
254 " {c_alpha^kappa} = {" << globalMolarities <<
"}, T = "
255 << fluidState.temperature(0));
260 template <
class Flu
idState>
263 std::cout <<
"saturations: ";
264 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
265 std::cout << fs.saturation(phaseIdx) <<
" ";
268 std::cout <<
"pressures: ";
269 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
270 std::cout << fs.pressure(phaseIdx) <<
" ";
273 std::cout <<
"densities: ";
274 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
275 std::cout << fs.density(phaseIdx) <<
" ";
278 std::cout <<
"molar densities: ";
279 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
280 std::cout << fs.molarDensity(phaseIdx) <<
" ";
283 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
284 std::cout <<
"composition " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
285 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
286 std::cout << fs.moleFraction(phaseIdx, compIdx) <<
" ";
291 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
292 std::cout <<
"fugacities " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
293 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
294 std::cout << fs.fugacity(phaseIdx, compIdx) <<
" ";
299 std::cout <<
"global component molarities: ";
300 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
302 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
303 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
305 std::cout << sum <<
" ";
310 template <
class MaterialLaw,
class Flu
idState>
313 FluidState &fluidState,
314 ParameterCache ¶mCache,
315 const MaterialLaw& material,
318 FluidState origFluidState(fluidState);
319 ParameterCache origParamCache(paramCache);
329 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
339 setQuantity_(fluidState, paramCache, material, pvIdx, x_i + eps);
347 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
348 J[eqIdx][pvIdx] = tmp[eqIdx];
351 fluidState = origFluidState;
352 paramCache = origParamCache;
359 template <
class Flu
idState>
361 const FluidState &fluidStateEval,
362 const FluidState &fluidState,
368 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
369 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
371 fluidState.fugacity(0, compIdx) -
372 fluidState.fugacity(phaseIdx, compIdx);
377 assert(eqIdx == numComponents*(numPhases - 1));
384 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
386 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
388 fluidState.saturation(phaseIdx)
389 * fluidState.molarity(phaseIdx, compIdx);
392 b[eqIdx] -= globalMolarities[compIdx];
398 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
399 Scalar sumMoleFracEval = 0.0;
400 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
401 sumMoleFracEval += fluidStateEval.moleFraction(phaseIdx, compIdx);
403 if (1.0 - sumMoleFracEval > fluidStateEval.saturation(phaseIdx)) {
404 b[eqIdx] = fluidState.saturation(phaseIdx);
407 Scalar sumMoleFrac = 0.0;
408 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
409 sumMoleFrac += fluidState.moleFraction(phaseIdx, compIdx);
410 b[eqIdx] = 1.0 - sumMoleFrac;
417 template <
class MaterialLaw,
class Flu
idState>
419 ParameterCache ¶mCache,
420 const MaterialLaw& material,
421 const Vector &deltaX)
424 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
426 Scalar delta = deltaX[pvIdx];
430 relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
435 delta = min(0.2, max(-0.2, delta));
440 delta = min(0.15, max(-0.15, delta));
445 delta = min(0.15*fluidState.pressure(0), max(-0.15*fluidState.pressure(0), delta));
453 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
454 Scalar value = fluidState.saturation(phaseIdx);
457 fluidState.setSaturation(phaseIdx, value);
461 value = fluidState.pressure(phaseIdx);
463 fluidState.setPressure(phaseIdx, 0.0);
465 bool allMoleFractionsAreZero =
true;
466 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
467 value = fluidState.moleFraction(phaseIdx, compIdx);
469 allMoleFractionsAreZero =
false;
471 fluidState.setMoleFraction(phaseIdx, compIdx, 0.0);
476 if (allMoleFractionsAreZero)
477 fluidState.setMoleFraction(phaseIdx, 0, 1e-10);
482 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
483 Scalar value = fluidState.saturation(phaseIdx)/(0.95*sumSat);
484 fluidState.setSaturation(phaseIdx, value);
493 template <
class MaterialLaw,
class Flu
idState>
495 ParameterCache ¶mCache,
496 const MaterialLaw& material)
501 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
502 sumSat += fluidState.saturation(phaseIdx);
503 fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
508 const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
509 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
510 fluidState.setPressure(phaseIdx,
511 fluidState.pressure(0)
512 + (pc[phaseIdx] - pc[0]));
515 paramCache.updateAll(fluidState, ParameterCache::Temperature);
518 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
522 fluidState.setDensity(phaseIdx, rho);
523 fluidState.setMolarDensity(phaseIdx, rhoMolar);
525 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
526 Scalar phi = FluidSystem::fugacityCoefficient( fluidState, paramCache, phaseIdx, compIdx);
527 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
533 {
return pvIdx == 0; }
536 {
return 1 <= pvIdx && pvIdx < numPhases; }
539 {
return numPhases <= pvIdx && pvIdx < numPhases + numPhases*numComponents; }
542 template <
class Flu
idState>
545 assert(pvIdx < numEq);
550 return fs.pressure(phaseIdx);
553 else if (pvIdx < numPhases) {
554 int phaseIdx = pvIdx - 1;
555 return fs.saturation(phaseIdx);
560 int phaseIdx = (pvIdx - numPhases)/numComponents;
561 int compIdx = (pvIdx - numPhases)%numComponents;
562 return fs.moleFraction(phaseIdx, compIdx);
567 template <
class MaterialLaw,
class Flu
idState>
569 ParameterCache ¶mCache,
570 const MaterialLaw& material,
574 assert(0 <= pvIdx && pvIdx < numEq);
577 Scalar delta = value - fs.pressure(0);
581 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
582 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
583 paramCache.updateAllPressures(fs);
586 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
589 fs.setDensity(phaseIdx, rho);
590 fs.setMolarDensity(phaseIdx, rhoMolar);
592 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
593 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
594 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
598 else if (pvIdx < numPhases) {
599 Scalar delta = value - fs.saturation(pvIdx - 1);
600 fs.setSaturation(pvIdx - 1, value);
604 fs.setSaturation(numPhases - 1,
605 fs.saturation(numPhases - 1) - delta);
609 const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
610 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
611 fs.setPressure(phaseIdx,
613 + (pc[phaseIdx] - pc[0]));
614 paramCache.updateAllPressures(fs);
617 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
620 fs.setDensity(phaseIdx, rho);
621 fs.setMolarDensity(phaseIdx, rhoMolar);
623 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
624 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
625 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
629 else if (pvIdx < numPhases + numPhases*numComponents)
631 int phaseIdx = (pvIdx - numPhases)/numComponents;
632 int compIdx = (pvIdx - numPhases)%numComponents;
634 fs.setMoleFraction(phaseIdx, compIdx, value);
635 paramCache.updateSingleMoleFraction(fs, phaseIdx, compIdx);
640 fs.setDensity(phaseIdx, rho);
641 fs.setMolarDensity(phaseIdx, rhoMolar);
645 if (!FluidSystem::isIdealMixture(phaseIdx)) {
646 for (
int compIdx2 = 0; compIdx2 < numComponents; ++compIdx2) {
647 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx2);
648 fs.setFugacityCoefficient(phaseIdx, compIdx2, phi);
658 template <
class Flu
idState>
661 assert(pvIdx < numEq);
666 fs.setPressure(phaseIdx, value);
669 else if (pvIdx < numPhases) {
670 int phaseIdx = pvIdx - 1;
671 fs.setSaturation(phaseIdx, value);
676 int phaseIdx = (pvIdx - numPhases)/numComponents;
677 int compIdx = (pvIdx - numPhases)%numComponents;
678 fs.setMoleFraction(phaseIdx, compIdx, value);
682 template <
class Flu
idState>
689 else if (pvIdx < numPhases)
697 int wettingPhaseIdx_ = 0;
Some exceptions thrown in DuMux
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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 phase compositions, pressures and saturations given the total mass of all components.
Definition: ncpflash.hh:79
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: ncpflash.hh:91
NcpFlash(int wettingPhaseIdx=0)
Construct a new flash.
Definition: ncpflash.hh:97
bool isSaturationIdx_(int pvIdx)
Definition: ncpflash.hh:535
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const Vector &deltaX)
Definition: ncpflash.hh:418
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material)
Definition: ncpflash.hh:494
bool isPressureIdx_(int pvIdx)
Definition: ncpflash.hh:532
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:543
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:311
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: ncpflash.hh:108
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:152
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: ncpflash.hh:659
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:683
bool isMoleFracIdx_(int pvIdx)
Definition: ncpflash.hh:538
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition: ncpflash.hh:568
void printFluidState_(const FluidState &fs)
Definition: ncpflash.hh:261
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:360