14#ifndef DUMUX_2P2C_VOLUME_VARIABLES_HH
15#define DUMUX_2P2C_VOLUME_VARIABLES_HH
31template <
class Traits,
bool enableChemicalNonEquilibrium,
bool useConstra
intSolver>
40template <
class Traits,
bool useConstra
intSolver = true>
49template <
class Traits,
class Impl>
57 using Scalar =
typename Traits::PrimaryVariables::value_type;
58 using ModelTraits =
typename Traits::ModelTraits;
64 comp0Idx = Traits::FluidSystem::comp0Idx,
65 comp1Idx = Traits::FluidSystem::comp1Idx,
66 phase0Idx = Traits::FluidSystem::phase0Idx,
67 phase1Idx = Traits::FluidSystem::phase1Idx
73 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
74 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
75 bothPhases = ModelTraits::Indices::bothPhases
81 switchIdx = ModelTraits::Indices::switchIdx,
82 pressureIdx = ModelTraits::Indices::pressureIdx
86 static constexpr auto formulation = ModelTraits::priVarFormulation();
88 using PermeabilityType =
typename Traits::PermeabilityType;
91 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
92 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
99 using Indices =
typename ModelTraits::Indices;
108 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
113 static_assert(ModelTraits::numFluidPhases() == 2,
"NumPhases set in the model is not two!");
114 static_assert(ModelTraits::numFluidComponents() == 2,
"NumComponents set in the model is not two!");
126 template<
class ElemSol,
class Problem,
class Element,
class Scv>
127 void update(
const ElemSol& elemSol,
const Problem& problem,
const Element& element,
const Scv& scv)
130 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
132 const auto& spatialParams = problem.spatialParams();
133 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
135 const int wPhaseIdx = fluidState_.wettingPhase();
136 const int nPhaseIdx = 1 - wPhaseIdx;
139 relativePermeability_[wPhaseIdx] = fluidMatrixInteraction.krw(
saturation(wPhaseIdx));
140 relativePermeability_[nPhaseIdx] = fluidMatrixInteraction.krn(
saturation(wPhaseIdx));
144 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
145 permeability_ = spatialParams.permeability(element, scv, elemSol);
147 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
149 return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx);
152 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
154 EnergyVolVars::updateEffectiveThermalConductivity();
170 template<
class ElemSol,
class Problem,
class Element,
class Scv>
172 const Problem& problem,
173 const Element& element,
180 const auto&
priVars = elemSol[scv.localDofIndex()];
183 const auto& spatialParams = problem.spatialParams();
184 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
185 const auto wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
213 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase presence.");
216 pc_ = fluidMatrixInteraction.pc(
fluidState.saturation(wPhaseIdx));
220 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ?
priVars[pressureIdx] + pc_
226 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ?
priVars[pressureIdx] - pc_
235 {
return fluidState_; }
241 {
return solidState_; }
249 {
return fluidState_.averageMolarMass(phaseIdx); }
258 {
return fluidState_.saturation(phaseIdx); }
268 {
return fluidState_.massFraction(phaseIdx, compIdx); }
278 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
287 {
return fluidState_.density(phaseIdx); }
296 {
return fluidState_.viscosity(phaseIdx); }
305 {
return fluidState_.molarDensity(phaseIdx); }
314 {
return fluidState_.pressure(phaseIdx); }
324 {
return fluidState_.temperature(0); }
333 {
return relativePermeability_[phaseIdx]; }
342 {
return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
355 {
return solidState_.porosity(); }
361 {
return permeability_; }
368 typename FluidSystem::ParameterCache paramCache;
369 paramCache.updatePhase(fluidState_, phaseIdx);
370 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
377 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
383 {
return fluidState_.wettingPhase(); }
390 PermeabilityType permeability_;
393 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
396 DiffusionCoefficients effectiveDiffCoeff_;
399 const Impl &
asImp_()
const {
return *
static_cast<const Impl*
>(
this); }
400 Impl &
asImp_() {
return *
static_cast<Impl*
>(
this); }
408template <
class Traits,
bool useConstra
intSolver>
415 using Scalar =
typename Traits::PrimaryVariables::value_type;
416 using ModelTraits =
typename Traits::ModelTraits;
418 static constexpr int numFluidComps = ParentType::numFluidComponents();
422 comp0Idx = Traits::FluidSystem::comp0Idx,
423 comp1Idx = Traits::FluidSystem::comp1Idx,
424 phase0Idx = Traits::FluidSystem::phase0Idx,
425 phase1Idx = Traits::FluidSystem::phase1Idx
431 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
432 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
433 bothPhases = ModelTraits::Indices::bothPhases
439 switchIdx = ModelTraits::Indices::switchIdx,
440 pressureIdx = ModelTraits::Indices::pressureIdx
444 static constexpr auto formulation = ModelTraits::priVarFormulation();
461 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
466 static_assert(useMoles() || (!useMoles() && useConstraintSolver),
"if !UseMoles, UseConstraintSolver has to be set to true");
470 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
471 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
486 template<
class ElemSol,
class Problem,
class Element,
class Scv>
488 const Problem& problem,
489 const Element& element,
494 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
496 const auto& priVars = elemSol[scv.localDofIndex()];
500 typename FluidSystem::ParameterCache paramCache;
503 if(!useConstraintSolver)
505 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
507 assert(FluidSystem::isIdealMixture(phaseIdx));
508 for (
int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
509 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
510 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
516 const Scalar p0 = fluidState.pressure(phase0Idx);
517 const Scalar p1 = fluidState.pressure(phase1Idx);
523 if(useConstraintSolver)
532 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0;
535 const Scalar partPressGas = p1 - partPressLiquid;
538 const Scalar xnn = partPressGas / p1;
539 const Scalar xnw = partPressLiquid / p1;
544 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0);
545 const Scalar xww = 1.0 - xwn;
548 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
549 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
550 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
551 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
559 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1 - priVars[switchIdx]);
560 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
564 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
569 if (useConstraintSolver)
578 const Scalar xnw = priVars[switchIdx];
579 const Scalar xnn = 1.0 - xnw;
588 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 );
595 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 );
597 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
598 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
607 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
608 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
612 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
617 if (useConstraintSolver)
626 const Scalar xwn = priVars[switchIdx];
631 const Scalar xnw = ( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 )/p1;
639 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 )/p1;
641 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
642 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
646 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
649 if(!useConstraintSolver)
651 paramCache.updateComposition(fluidState, phaseIdx);
653 fluidState.setDensity(phaseIdx, rho);
655 fluidState.setMolarDensity(phaseIdx, rhoMolar);
660 fluidState.setViscosity(phaseIdx,mu);
661 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
662 fluidState.setEnthalpy(phaseIdx, h);
674template <
class Traits,
bool useConstra
intSolver>
681 using Scalar =
typename Traits::PrimaryVariables::value_type;
682 using ModelTraits =
typename Traits::ModelTraits;
684 static constexpr int numFluidComps = ParentType::numFluidComponents();
688 comp0Idx = Traits::FluidSystem::comp0Idx,
689 comp1Idx = Traits::FluidSystem::comp1Idx,
690 phase0Idx = Traits::FluidSystem::phase0Idx,
691 phase1Idx = Traits::FluidSystem::phase1Idx
697 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
698 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
699 bothPhases = ModelTraits::Indices::bothPhases
705 switchIdx = ModelTraits::Indices::switchIdx,
706 pressureIdx = ModelTraits::Indices::pressureIdx
710 static constexpr auto formulation = ModelTraits::priVarFormulation();
712 using PermeabilityType =
typename Traits::PermeabilityType;
729 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
734 static_assert(useMoles() || (!useMoles() && useConstraintSolver),
"if !UseMoles, UseConstraintSolver has to be set to true");
739 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
740 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
755 template<
class ElemSol,
class Problem,
class Element,
class Scv>
757 const Problem& problem,
758 const Element& element,
763 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
765 const auto& priVars = elemSol[scv.localDofIndex()];
770 typename FluidSystem::ParameterCache paramCache;
771 paramCache.updateAll(fluidState);
773 updateMoleFraction(fluidState,
778 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
782 fluidState.setViscosity(phaseIdx,mu);
783 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
784 fluidState.setEnthalpy(phaseIdx, h);
796 typename Traits::FluidSystem::ParameterCache & paramCache,
797 const typename Traits::PrimaryVariables& priVars)
801 Scalar xwnNonEquil = 0.0;
802 Scalar xwwNonEquil = 0.0;
803 Scalar xnwNonEquil = 0.0;
804 Scalar xnnNonEquil = 0.0;
808 xwnNonEquil = priVars[ModelTraits::numFluidPhases()];
809 xwwNonEquil = 1-xwnNonEquil;
810 xnwNonEquil = priVars[ModelTraits::numFluidPhases()+comp1Idx];
813 if (actualFluidState.saturation(phase0Idx) < 0.01)
815 const Scalar p1 = actualFluidState.pressure(phase1Idx);
816 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
817 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
819 xnnNonEquil = 1- xnwNonEquil;
821 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
822 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
823 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
824 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
827 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
828 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
830 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
834 actualFluidState.setFugacityCoefficient(phaseIdx,
840 equilFluidState.assign(actualFluidState) ;
842 if(!useConstraintSolver)
844 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
846 assert(FluidSystem::isIdealMixture(phaseIdx));
847 for (
int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
848 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
849 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
855 const Scalar p0 = equilFluidState.pressure(phase0Idx);
856 const Scalar p1 = equilFluidState.pressure(phase1Idx);
860 if(useConstraintSolver)
869 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
872 const Scalar partPressGas = p1 - partPressLiquid;
875 const Scalar xnn = partPressGas / p1;
876 const Scalar xnw = partPressLiquid / p1;
881 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
882 const Scalar xww = 1.0 - xwn;
885 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
886 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
887 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
888 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
892 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
893 for (
int compIdx=0; compIdx< numFluidComps; ++ compIdx){
894 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
900 DUNE_THROW(Dune::InvalidStateException,
"nonequilibrium is only possible for 2 phases present ");
903 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
906 actualFluidState.setDensity(phaseIdx, rho);
908 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
919 const Scalar
xEquil(
const unsigned int phaseIdx,
const unsigned int compIdx)
const
921 return xEquil_[phaseIdx][compIdx] ;
925 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
The primary variable switch for the 2pnc model.
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:64
static void solve(FluidState &fluidState, ParameterCache ¶mCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:99
Definition: porousmediumflow/nonisothermal/volumevariables.hh:63
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:47
static void solve(FluidState &fluidState, ParameterCache ¶mCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:69
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:28
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition: porousmediumflow/volumevariables.hh:40
const PrimaryVariables & priVars() const
Returns the vector of primary variables.
Definition: porousmediumflow/volumevariables.hh:64
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/volumevariables.hh:52
The primary variable switch controlling the phase presence state variable.
Definition: 2pnc/primaryvariableswitch.hh:29
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition: porousmediumflow/2p2c/volumevariables.hh:53
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/2p2c/volumevariables.hh:248
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:171
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:360
typename ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/2p2c/volumevariables.hh:99
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:354
const Impl & asImp_() const
Definition: porousmediumflow/2p2c/volumevariables.hh:399
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:97
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:341
Impl & asImp_()
Definition: porousmediumflow/2p2c/volumevariables.hh:400
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:332
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/2p2c/volumevariables.hh:382
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:127
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:376
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:323
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:95
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:348
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:366
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:304
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:234
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:110
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:286
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:101
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:257
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:108
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:103
Scalar moleFraction(const int phaseIdx, const int compIdx) const
Returns the mole fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:277
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:240
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:313
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:295
Scalar massFraction(const int phaseIdx, const int compIdx) const
Returns the mass fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:267
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:450
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:454
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:487
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:452
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:456
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:461
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:463
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:756
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:720
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:722
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:731
const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
The mole fraction we would have in the case of chemical equilibrium / on the interface.
Definition: porousmediumflow/2p2c/volumevariables.hh:919
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:718
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:729
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:724
void updateMoleFraction(FluidState &actualFluidState, typename Traits::FluidSystem::ParameterCache ¶mCache, const typename Traits::PrimaryVariables &priVars)
Updates composition of all phases from the primary variables.
Definition: porousmediumflow/2p2c/volumevariables.hh:795
Definition: porousmediumflow/2p2c/volumevariables.hh:32
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Computes all quantities of a generic fluid state if a reference phase has been specified.
void updateSolidVolumeFractions(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, SolidState &solidState, const int solidVolFracOffset)
update the solid volume fractions (inert and reacitve) and set them in the solidstate
Definition: updatesolidvolumefractions.hh:24
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition: formulation.hh:23
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
The available discretization methods in Dumux.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
std::string phasePresence() noexcept
I/O name of phase presence.
Definition: name.hh:135
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:62
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
Base class for the model specific class which provides access to all volume averaged quantities.
Base class for the model specific class which provides access to all volume averaged quantities.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.