25#ifndef DUMUX_NCP_FLASH_HH
26#define DUMUX_NCP_FLASH_HH
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30#include <dune/common/version.hh>
78template <
class Scalar,
class Flu
idSystem>
81 enum { numPhases = FluidSystem::numPhases };
82 enum { numComponents = FluidSystem::numComponents };
84 using ParameterCache =
typename FluidSystem::ParameterCache;
86 static constexpr int numEq = numPhases*(numComponents + 1);
88 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
89 using Vector = Dune::FieldVector<Scalar, numEq>;
99 : wettingPhaseIdx_(wettingPhaseIdx)
108 template <
class Flu
idState>
110 ParameterCache ¶mCache,
115 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
116 sumMoles += globalMolarities[compIdx];
118 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
120 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx)
121 fluidState.setMoleFraction(phaseIdx,
123 globalMolarities[compIdx]/sumMoles);
126 fluidState.setPressure(phaseIdx, 1.0135e5);
129 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
133 paramCache.updateAll(fluidState);
134 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
135 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
136 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
137 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
152 template <
class MaterialLaw,
class Flu
idState>
154 ParameterCache ¶mCache,
155 const typename MaterialLaw::Params &matParams,
158#if DUNE_VERSION_LT(DUNE_COMMON, 2, 7)
159 Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-25);
178 completeFluidState_<MaterialLaw>(fluidState,
191 for (
int nIdx = 0; nIdx < nMax; ++nIdx) {
193 linearize_<MaterialLaw>(J, b, fluidState, paramCache, matParams, globalMolarities);
204 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
206 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
217 catch (Dune::FMatrixError &e)
252 Scalar relError = update_<MaterialLaw>(fluidState, paramCache, matParams, deltaX);
267 "Flash calculation failed."
268 " {c_alpha^kappa} = {" << globalMolarities <<
"}, T = "
269 << fluidState.temperature(0));
274 template <
class Flu
idState>
277 std::cout <<
"saturations: ";
278 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
279 std::cout << fs.saturation(phaseIdx) <<
" ";
282 std::cout <<
"pressures: ";
283 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
284 std::cout << fs.pressure(phaseIdx) <<
" ";
287 std::cout <<
"densities: ";
288 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
289 std::cout << fs.density(phaseIdx) <<
" ";
292 std::cout <<
"molar densities: ";
293 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
294 std::cout << fs.molarDensity(phaseIdx) <<
" ";
297 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
298 std::cout <<
"composition " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
299 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
300 std::cout << fs.moleFraction(phaseIdx, compIdx) <<
" ";
305 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
306 std::cout <<
"fugacities " << FluidSystem::phaseName(phaseIdx) <<
"Phase: ";
307 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
308 std::cout << fs.fugacity(phaseIdx, compIdx) <<
" ";
313 std::cout <<
"global component molarities: ";
314 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
316 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
317 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
319 std::cout << sum <<
" ";
324 template <
class MaterialLaw,
class Flu
idState>
327 FluidState &fluidState,
328 ParameterCache ¶mCache,
329 const typename MaterialLaw::Params &matParams,
332 FluidState origFluidState(fluidState);
333 ParameterCache origParamCache(paramCache);
345 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
355 setQuantity_<MaterialLaw>(fluidState, paramCache, matParams, pvIdx, x_i +
eps);
363 for (
int eqIdx = 0; eqIdx < numEq; ++eqIdx)
364 J[eqIdx][pvIdx] = tmp[eqIdx];
367 fluidState = origFluidState;
368 paramCache = origParamCache;
375 template <
class Flu
idState>
377 const FluidState &fluidStateEval,
378 const FluidState &fluidState,
384 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
385 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
387 fluidState.fugacity(0, compIdx) -
388 fluidState.fugacity(phaseIdx, compIdx);
393 assert(eqIdx == numComponents*(numPhases - 1));
400 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
402 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
404 fluidState.saturation(phaseIdx)
405 * fluidState.molarity(phaseIdx, compIdx);
408 b[eqIdx] -= globalMolarities[compIdx];
414 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
415 Scalar sumMoleFracEval = 0.0;
416 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
417 sumMoleFracEval += fluidStateEval.moleFraction(phaseIdx, compIdx);
419 if (1.0 - sumMoleFracEval > fluidStateEval.saturation(phaseIdx)) {
420 b[eqIdx] = fluidState.saturation(phaseIdx);
423 Scalar sumMoleFrac = 0.0;
424 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
425 sumMoleFrac += fluidState.moleFraction(phaseIdx, compIdx);
426 b[eqIdx] = 1.0 - sumMoleFrac;
433 template <
class MaterialLaw,
class Flu
idState>
435 ParameterCache ¶mCache,
436 const typename MaterialLaw::Params &matParams,
437 const Vector &deltaX)
440 for (
int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
442 Scalar delta = deltaX[pvIdx];
446 relError = max(relError, abs(delta)*
quantityWeight_(fluidState, pvIdx));
451 delta = min(0.2, max(-0.2, delta));
456 delta = min(0.15, max(-0.15, delta));
461 delta = min(0.15*fluidState.pressure(0), max(-0.15*fluidState.pressure(0), delta));
469 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
470 Scalar value = fluidState.saturation(phaseIdx);
473 fluidState.setSaturation(phaseIdx, value);
477 value = fluidState.pressure(phaseIdx);
479 fluidState.setPressure(phaseIdx, 0.0);
481 bool allMoleFractionsAreZero =
true;
482 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
483 value = fluidState.moleFraction(phaseIdx, compIdx);
485 allMoleFractionsAreZero =
false;
487 fluidState.setMoleFraction(phaseIdx, compIdx, 0.0);
492 if (allMoleFractionsAreZero)
493 fluidState.setMoleFraction(phaseIdx, 0, 1e-10);
498 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
499 Scalar value = fluidState.saturation(phaseIdx)/(0.95*sumSat);
500 fluidState.setSaturation(phaseIdx, value);
504 completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
509 template <
class MaterialLaw,
class Flu
idState>
511 ParameterCache ¶mCache,
512 const typename MaterialLaw::Params &matParams)
517 for (
int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
518 sumSat += fluidState.saturation(phaseIdx);
519 fluidState.setSaturation(numPhases - 1, 1.0 - sumSat);
525 MaterialLaw::capillaryPressures(pc, matParams, fluidState, wettingPhaseIdx_);
526 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
527 fluidState.setPressure(phaseIdx,
528 fluidState.pressure(0)
529 + (pc[phaseIdx] - pc[0]));
532 paramCache.updateAll(fluidState, ParameterCache::Temperature);
535 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
539 fluidState.setDensity(phaseIdx, rho);
540 fluidState.setMolarDensity(phaseIdx, rhoMolar);
542 for (
int compIdx = 0; compIdx < numComponents; ++ compIdx) {
543 Scalar phi = FluidSystem::fugacityCoefficient( fluidState, paramCache, phaseIdx, compIdx);
544 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
550 {
return pvIdx == 0; }
553 {
return 1 <= pvIdx && pvIdx < numPhases; }
556 {
return numPhases <= pvIdx && pvIdx < numPhases + numPhases*numComponents; }
559 template <
class Flu
idState>
562 assert(pvIdx < numEq);
567 return fs.pressure(phaseIdx);
570 else if (pvIdx < numPhases) {
571 int phaseIdx = pvIdx - 1;
572 return fs.saturation(phaseIdx);
577 int phaseIdx = (pvIdx - numPhases)/numComponents;
578 int compIdx = (pvIdx - numPhases)%numComponents;
579 return fs.moleFraction(phaseIdx, compIdx);
584 template <
class MaterialLaw,
class Flu
idState>
586 ParameterCache ¶mCache,
587 const typename MaterialLaw::Params &matParams,
591 assert(0 <= pvIdx && pvIdx < numEq);
594 Scalar delta = value - fs.pressure(0);
598 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
599 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
600 paramCache.updateAllPressures(fs);
603 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
606 fs.setDensity(phaseIdx, rho);
607 fs.setMolarDensity(phaseIdx, rhoMolar);
609 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
610 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
611 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
615 else if (pvIdx < numPhases) {
616 Scalar delta = value - fs.saturation(pvIdx - 1);
617 fs.setSaturation(pvIdx - 1, value);
621 fs.setSaturation(numPhases - 1,
622 fs.saturation(numPhases - 1) - delta);
627 MaterialLaw::capillaryPressures(pc, matParams, fs, wettingPhaseIdx_);
628 for (
int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
629 fs.setPressure(phaseIdx,
631 + (pc[phaseIdx] - pc[0]));
632 paramCache.updateAllPressures(fs);
635 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
638 fs.setDensity(phaseIdx, rho);
639 fs.setMolarDensity(phaseIdx, rhoMolar);
641 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
642 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
643 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
647 else if (pvIdx < numPhases + numPhases*numComponents)
649 int phaseIdx = (pvIdx - numPhases)/numComponents;
650 int compIdx = (pvIdx - numPhases)%numComponents;
652 fs.setMoleFraction(phaseIdx, compIdx, value);
653 paramCache.updateSingleMoleFraction(fs, phaseIdx, compIdx);
658 fs.setDensity(phaseIdx, rho);
659 fs.setMolarDensity(phaseIdx, rhoMolar);
663 if (!FluidSystem::isIdealMixture(phaseIdx)) {
664 for (
int compIdx2 = 0; compIdx2 < numComponents; ++compIdx2) {
665 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx2);
666 fs.setFugacityCoefficient(phaseIdx, compIdx2, phi);
676 template <
class Flu
idState>
679 assert(pvIdx < numEq);
684 fs.setPressure(phaseIdx, value);
687 else if (pvIdx < numPhases) {
688 int phaseIdx = pvIdx - 1;
689 fs.setSaturation(phaseIdx, value);
694 int phaseIdx = (pvIdx - numPhases)/numComponents;
695 int compIdx = (pvIdx - numPhases)%numComponents;
696 fs.setMoleFraction(phaseIdx, compIdx, value);
700 template <
class Flu
idState>
707 else if (pvIdx < numPhases)
715 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 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 phase compositions, pressures and saturations given the total mass of all components.
Definition: ncpflash.hh:80
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: ncpflash.hh:92
NcpFlash(int wettingPhaseIdx=0)
Contruct a new flash.
Definition: ncpflash.hh:98
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: ncpflash.hh:153
bool isSaturationIdx_(int pvIdx)
Definition: ncpflash.hh:552
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:325
void setQuantity_(FluidState &fs, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, int pvIdx, Scalar value)
Definition: ncpflash.hh:585
bool isPressureIdx_(int pvIdx)
Definition: ncpflash.hh:549
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:560
void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: ncpflash.hh:109
void completeFluidState_(FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams)
Definition: ncpflash.hh:510
Scalar update_(FluidState &fluidState, ParameterCache ¶mCache, const typename MaterialLaw::Params &matParams, const Vector &deltaX)
Definition: ncpflash.hh:434
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: ncpflash.hh:677
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:701
bool isMoleFracIdx_(int pvIdx)
Definition: ncpflash.hh:555
void printFluidState_(const FluidState &fs)
Definition: ncpflash.hh:275
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:376