26#ifndef DUMUX_2P2C_VOLUME_VARIABLES_HH
27#define DUMUX_2P2C_VOLUME_VARIABLES_HH
43template <
class Traits,
bool enableChemicalNonEquilibrium,
bool useConstra
intSolver>
52template <
class Traits,
bool useConstra
intSolver = true>
61template <
class Traits,
class Impl>
69 using Scalar =
typename Traits::PrimaryVariables::value_type;
70 using ModelTraits =
typename Traits::ModelTraits;
76 comp0Idx = Traits::FluidSystem::comp0Idx,
77 comp1Idx = Traits::FluidSystem::comp1Idx,
78 phase0Idx = Traits::FluidSystem::phase0Idx,
79 phase1Idx = Traits::FluidSystem::phase1Idx
85 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
86 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
87 bothPhases = ModelTraits::Indices::bothPhases
93 switchIdx = ModelTraits::Indices::switchIdx,
94 pressureIdx = ModelTraits::Indices::pressureIdx
98 static constexpr auto formulation = ModelTraits::priVarFormulation();
100 using PermeabilityType =
typename Traits::PermeabilityType;
103 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
104 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
111 using Indices =
typename ModelTraits::Indices;
120 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
125 static_assert(ModelTraits::numFluidPhases() == 2,
"NumPhases set in the model is not two!");
126 static_assert(ModelTraits::numFluidComponents() == 2,
"NumComponents set in the model is not two!");
138 template<
class ElemSol,
class Problem,
class Element,
class Scv>
139 void update(
const ElemSol& elemSol,
const Problem& problem,
const Element& element,
const Scv& scv)
142 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
144 const auto& spatialParams = problem.spatialParams();
145 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
147 const int wPhaseIdx = fluidState_.wettingPhase();
148 const int nPhaseIdx = 1 - wPhaseIdx;
151 relativePermeability_[wPhaseIdx] = fluidMatrixInteraction.krw(
saturation(wPhaseIdx));
152 relativePermeability_[nPhaseIdx] = fluidMatrixInteraction.krn(
saturation(wPhaseIdx));
156 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
157 permeability_ = spatialParams.permeability(element, scv, elemSol);
159 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
161 return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx);
164 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
166 EnergyVolVars::updateEffectiveThermalConductivity();
182 template<
class ElemSol,
class Problem,
class Element,
class Scv>
184 const Problem& problem,
185 const Element& element,
192 const auto&
priVars = elemSol[scv.localDofIndex()];
193 const auto phasePresence =
priVars.state();
195 const auto& spatialParams = problem.spatialParams();
196 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
201 if (phasePresence == firstPhaseOnly)
206 else if (phasePresence == secondPhaseOnly)
211 else if (phasePresence == bothPhases)
225 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase presence.");
228 pc_ = fluidMatrixInteraction.pc(
fluidState.saturation(wPhaseIdx));
232 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ?
priVars[pressureIdx] + pc_
238 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ?
priVars[pressureIdx] - pc_
247 {
return fluidState_; }
253 {
return solidState_; }
261 {
return fluidState_.averageMolarMass(phaseIdx); }
270 {
return fluidState_.saturation(phaseIdx); }
280 {
return fluidState_.massFraction(phaseIdx, compIdx); }
290 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
299 {
return fluidState_.density(phaseIdx); }
308 {
return fluidState_.viscosity(phaseIdx); }
317 {
return fluidState_.molarDensity(phaseIdx); }
326 {
return fluidState_.pressure(phaseIdx); }
336 {
return fluidState_.temperature(0); }
345 {
return relativePermeability_[phaseIdx]; }
354 {
return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
367 {
return solidState_.porosity(); }
373 {
return permeability_; }
380 typename FluidSystem::ParameterCache paramCache;
381 paramCache.updatePhase(fluidState_, phaseIdx);
382 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
389 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
395 {
return fluidState_.wettingPhase(); }
402 PermeabilityType permeability_;
405 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
408 DiffusionCoefficients effectiveDiffCoeff_;
411 const Impl &
asImp_()
const {
return *
static_cast<const Impl*
>(
this); }
412 Impl &
asImp_() {
return *
static_cast<Impl*
>(
this); }
420template <
class Traits,
bool useConstra
intSolver>
427 using Scalar =
typename Traits::PrimaryVariables::value_type;
428 using ModelTraits =
typename Traits::ModelTraits;
434 comp0Idx = Traits::FluidSystem::comp0Idx,
435 comp1Idx = Traits::FluidSystem::comp1Idx,
436 phase0Idx = Traits::FluidSystem::phase0Idx,
437 phase1Idx = Traits::FluidSystem::phase1Idx
443 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
444 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
445 bothPhases = ModelTraits::Indices::bothPhases
451 switchIdx = ModelTraits::Indices::switchIdx,
452 pressureIdx = ModelTraits::Indices::pressureIdx
456 static constexpr auto formulation = ModelTraits::priVarFormulation();
473 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
478 static_assert(useMoles() || (!useMoles() && useConstraintSolver),
"if !UseMoles, UseConstraintSolver has to be set to true");
482 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
483 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
498 template<
class ElemSol,
class Problem,
class Element,
class Scv>
500 const Problem& problem,
501 const Element& element,
508 const auto&
priVars = elemSol[scv.localDofIndex()];
509 const auto phasePresence =
priVars.state();
512 typename FluidSystem::ParameterCache paramCache;
515 if(!useConstraintSolver)
517 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
519 assert(FluidSystem::isIdealMixture(phaseIdx));
520 for (
int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
521 Scalar phi = FluidSystem::fugacityCoefficient(
fluidState, paramCache, phaseIdx, compIdx);
522 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
528 const Scalar p0 =
fluidState.pressure(phase0Idx);
529 const Scalar p1 =
fluidState.pressure(phase1Idx);
530 if (phasePresence == bothPhases)
535 if(useConstraintSolver)
544 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp0Idx)*p0;
547 const Scalar partPressGas = p1 - partPressLiquid;
550 const Scalar xnn = partPressGas / p1;
551 const Scalar xnw = partPressLiquid / p1;
556 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp1Idx)*p0);
557 const Scalar xww = 1.0 - xwn;
560 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
561 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
562 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
563 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
566 else if (phasePresence == secondPhaseOnly)
581 if (useConstraintSolver)
590 const Scalar xnw =
priVars[switchIdx];
591 const Scalar xnn = 1.0 - xnw;
600 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp0Idx)*p0 );
607 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp1Idx)*p0 );
609 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
610 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
613 else if (phasePresence == firstPhaseOnly)
629 if (useConstraintSolver)
638 const Scalar xwn =
priVars[switchIdx];
643 const Scalar xnw = ( FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp0Idx)*p0 )/p1;
651 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(
fluidState, phase0Idx, comp1Idx)*p0 )/p1;
653 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
654 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
658 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
661 if(!useConstraintSolver)
663 paramCache.updateComposition(
fluidState, phaseIdx);
664 const Scalar rho = FluidSystem::density(
fluidState, paramCache, phaseIdx);
666 Scalar rhoMolar = FluidSystem::molarDensity(
fluidState, paramCache, phaseIdx);
667 fluidState.setMolarDensity(phaseIdx, rhoMolar);
671 const Scalar mu = FluidSystem::viscosity(
fluidState, paramCache, phaseIdx);
673 Scalar h = EnergyVolVars::enthalpy(
fluidState, paramCache, phaseIdx);
686template <
class Traits,
bool useConstra
intSolver>
693 using Scalar =
typename Traits::PrimaryVariables::value_type;
694 using ModelTraits =
typename Traits::ModelTraits;
700 comp0Idx = Traits::FluidSystem::comp0Idx,
701 comp1Idx = Traits::FluidSystem::comp1Idx,
702 phase0Idx = Traits::FluidSystem::phase0Idx,
703 phase1Idx = Traits::FluidSystem::phase1Idx
709 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
710 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
711 bothPhases = ModelTraits::Indices::bothPhases
717 switchIdx = ModelTraits::Indices::switchIdx,
718 pressureIdx = ModelTraits::Indices::pressureIdx
722 static constexpr auto formulation = ModelTraits::priVarFormulation();
724 using PermeabilityType =
typename Traits::PermeabilityType;
741 static constexpr bool useMoles() {
return ModelTraits::useMoles(); }
746 static_assert(useMoles() || (!useMoles() && useConstraintSolver),
"if !UseMoles, UseConstraintSolver has to be set to true");
751 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
752 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
767 template<
class ElemSol,
class Problem,
class Element,
class Scv>
769 const Problem& problem,
770 const Element& element,
777 const auto&
priVars = elemSol[scv.localDofIndex()];
782 typename FluidSystem::ParameterCache paramCache;
790 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
793 const Scalar mu = FluidSystem::viscosity(
fluidState, paramCache, phaseIdx);
795 Scalar h = EnergyVolVars::enthalpy(
fluidState, paramCache, phaseIdx);
808 typename Traits::FluidSystem::ParameterCache & paramCache,
809 const typename Traits::PrimaryVariables&
priVars)
811 const auto phasePresence =
priVars.state();
813 Scalar xwnNonEquil = 0.0;
814 Scalar xwwNonEquil = 0.0;
815 Scalar xnwNonEquil = 0.0;
816 Scalar xnnNonEquil = 0.0;
818 if (phasePresence == bothPhases)
820 xwnNonEquil =
priVars[ModelTraits::numFluidPhases()];
821 xwwNonEquil = 1-xwnNonEquil;
822 xnwNonEquil =
priVars[ModelTraits::numFluidPhases()+comp1Idx];
825 if (actualFluidState.saturation(phase0Idx) < 0.01)
827 const Scalar p1 = actualFluidState.pressure(phase1Idx);
828 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
829 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
831 xnnNonEquil = 1- xnwNonEquil;
833 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
834 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
835 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
836 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
839 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
840 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
842 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
846 actualFluidState.setFugacityCoefficient(phaseIdx,
852 equilFluidState.assign(actualFluidState) ;
854 if(!useConstraintSolver)
856 for (
int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
858 assert(FluidSystem::isIdealMixture(phaseIdx));
859 for (
int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
860 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
861 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
867 const Scalar p0 = equilFluidState.pressure(phase0Idx);
868 const Scalar p1 = equilFluidState.pressure(phase1Idx);
872 if(useConstraintSolver)
881 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
884 const Scalar partPressGas = p1 - partPressLiquid;
887 const Scalar xnn = partPressGas / p1;
888 const Scalar xnw = partPressLiquid / p1;
893 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
894 const Scalar xww = 1.0 - xwn;
897 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
898 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
899 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
900 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
904 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
905 for (
int compIdx=0; compIdx< numFluidComps; ++ compIdx){
906 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
912 DUNE_THROW(Dune::InvalidStateException,
"nonequilibrium is only possible for 2 phases present ");
915 for(
int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
917 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
918 actualFluidState.setDensity(phaseIdx, rho);
919 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
920 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
931 const Scalar
xEquil(
const unsigned int phaseIdx,
const unsigned int compIdx)
const
933 return xEquil_[phaseIdx][compIdx] ;
937 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
The available discretization methods in Dumux.
Computes all quantities of a generic fluid state if a reference phase has been specified.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.
Defines an enumeration for the formulations accepted by the two-phase model.
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition formulation.hh:35
@ p1s0
first phase saturation and second phase pressure as primary variables
Definition formulation.hh:37
@ p0s1
first phase pressure and second phase saturation as primary variables
Definition formulation.hh:36
TwoPTwoCVolumeVariablesImplementation< Traits, Traits::ModelTraits::enableChemicalNonEquilibrium(), useConstraintSolver > TwoPTwoCVolumeVariables
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition porousmediumflow/2p2c/volumevariables.hh:53
EnergyVolumeVariablesImplementation< IsothermalTraits, Impl, IsothermalTraits::ModelTraits::enableEnergyBalance()> EnergyVolumeVariables
Base class for the model specific class which provides access to all volume averaged quantities.
Definition porousmediumflow/nonisothermal/volumevariables.hh:85
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:36
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition computefromreferencephase.hh:76
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:111
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:81
Definition porousmediumflow/2p2c/volumevariables.hh:44
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition porousmediumflow/2p2c/volumevariables.hh:65
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition porousmediumflow/2p2c/volumevariables.hh:260
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:183
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:372
typename ModelTraits::Indices Indices
Export the indices.
Definition porousmediumflow/2p2c/volumevariables.hh:111
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:366
const Impl & asImp_() const
Definition porousmediumflow/2p2c/volumevariables.hh:411
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition porousmediumflow/2p2c/volumevariables.hh:109
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:353
Impl & asImp_()
Definition porousmediumflow/2p2c/volumevariables.hh:412
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:344
int wettingPhase() const
Returns the wetting phase index.
Definition porousmediumflow/2p2c/volumevariables.hh:394
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:139
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition porousmediumflow/2p2c/volumevariables.hh:388
TwoPNCPrimaryVariableSwitch PrimaryVariableSwitch
Export the primary variable switch.
Definition porousmediumflow/2p2c/volumevariables.hh:117
Scalar temperature() const
Returns temperature within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:335
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition porousmediumflow/2p2c/volumevariables.hh:107
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:360
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition porousmediumflow/2p2c/volumevariables.hh:378
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:316
const FluidState & fluidState() const
Definition porousmediumflow/2p2c/volumevariables.hh:246
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition porousmediumflow/2p2c/volumevariables.hh:122
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:298
typename Traits::SolidState SolidState
Export type of solid state.
Definition porousmediumflow/2p2c/volumevariables.hh:113
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:269
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition porousmediumflow/2p2c/volumevariables.hh:120
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition porousmediumflow/2p2c/volumevariables.hh:115
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:289
const SolidState & solidState() const
Definition porousmediumflow/2p2c/volumevariables.hh:252
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:325
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition porousmediumflow/2p2c/volumevariables.hh:307
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:279
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition porousmediumflow/2p2c/volumevariables.hh:462
typename Traits::SolidState SolidState
Export type of solid state.
Definition porousmediumflow/2p2c/volumevariables.hh:466
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:499
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition porousmediumflow/2p2c/volumevariables.hh:464
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition porousmediumflow/2p2c/volumevariables.hh:468
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition porousmediumflow/2p2c/volumevariables.hh:473
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition porousmediumflow/2p2c/volumevariables.hh:475
TwoPNCPrimaryVariableSwitch PrimaryVariableSwitch
Export the primary variable switch.
Definition porousmediumflow/2p2c/volumevariables.hh:470
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:768
TwoPNCPrimaryVariableSwitch PrimaryVariableSwitch
Export the primary variable switch.
Definition porousmediumflow/2p2c/volumevariables.hh:738
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition porousmediumflow/2p2c/volumevariables.hh:732
typename Traits::SolidState SolidState
Export type of solid state.
Definition porousmediumflow/2p2c/volumevariables.hh:734
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition porousmediumflow/2p2c/volumevariables.hh:743
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:931
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition porousmediumflow/2p2c/volumevariables.hh:730
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition porousmediumflow/2p2c/volumevariables.hh:741
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition porousmediumflow/2p2c/volumevariables.hh:736
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:807
The primary variable switch controlling the phase presence state variable.
Definition 2pnc/primaryvariableswitch.hh:41
The isothermal base class.
Definition porousmediumflow/volumevariables.hh:40
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition porousmediumflow/volumevariables.hh:52
const PrimaryVariables & priVars() const
Definition porousmediumflow/volumevariables.hh:76
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:64
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.
The primary variable switch for the 2pnc model.