26#ifndef DUMUX_MPNC_VOLUME_VARIABLES_HH
27#define DUMUX_MPNC_VOLUME_VARIABLES_HH
41template <
class Traits,
bool enableChemicalNonEquilibrium>
49template <
class Traits>
52template <
class Traits>
60 using Scalar =
typename Traits::PrimaryVariables::value_type;
61 using PermeabilityType =
typename Traits::PermeabilityType;
63 using ModelTraits =
typename Traits::ModelTraits;
64 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
66 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
67 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
68 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
70 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
72 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
73 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
77 using Indices =
typename Traits::ModelTraits::Indices;
88 static constexpr int numFluidPhases() {
return ModelTraits::numFluidPhases(); }
90 static constexpr int numFluidComps = ParentType::numFluidComponents();
101 template<
class ElemSol,
class Problem,
class Element,
class Scv>
103 const Problem& problem,
104 const Element& element,
107 ParentType::update(elemSol, problem, element, scv);
109 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
112 const auto& spatialParams = problem.spatialParams();
113 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
114 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
115 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
116 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
118 typename FluidSystem::ParameterCache paramCache;
119 paramCache.updateAll(fluidState_);
124 if constexpr (enableDiffusion)
126 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
127 {
return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx); };
129 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
132 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
133 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
134 EnergyVolVars::updateEffectiveThermalConductivity();
149 template<
class ElemSol,
class Problem,
class Element,
class Scv>
151 const Problem& problem,
152 const Element& element,
161 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
166 auto&& priVars = elemSol[scv.localDofIndex()];
168 for (
int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
169 sumSat += priVars[Indices::s0Idx + phaseIdx];
170 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
172 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
178 const auto& spatialParams = problem.spatialParams();
179 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
180 fluidState.setWettingPhase(wPhaseIdx);
183 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
184 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
191 const Scalar pw = priVars[Indices::p0Idx];
192 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
193 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
198 const Scalar pn = priVars[Indices::p0Idx];
199 for (
int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
200 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
203 DUNE_THROW(Dune::InvalidStateException,
"MPNCVolumeVariables do not support the chosen pressure formulation");
208 typename FluidSystem::ParameterCache paramCache;
209 paramCache.updateAll(fluidState);
213 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
214 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
217 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
219 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
220 Scalar x_ij = 1.0/numFluidComps;
223 fluidState.setMoleFraction(phaseIdx,
234 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
237 fluidState.setViscosity(phaseIdx, mu);
240 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
241 fluidState.setEnthalpy(phaseIdx, h);
250 {
return fluidState_; }
256 {
return solidState_; }
265 {
return fluidState_.saturation(phaseIdx); }
275 {
return fluidState_.massFraction(phaseIdx, compIdx); }
285 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
293 Scalar
molarity(
const int phaseIdx,
int compIdx)
const
294 {
return fluidState_.molarity(phaseIdx, compIdx); }
302 {
return fluidState_.molarDensity(phaseIdx);}
311 {
return fluidState_.pressure(phaseIdx); }
319 {
return fluidState_.density(phaseIdx); }
329 {
return fluidState_.temperature(0); }
332 {
return fluidState_.temperature(phaseIdx); }
338 {
return fluidState_.enthalpy(phaseIdx); }
344 {
return fluidState_.internalEnergy(phaseIdx); }
351 {
return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
357 {
return fluidState_.fugacity(compIdx); }
363 {
return fluidState_.averageMolarMass(phaseIdx); }
381 {
return fluidState_.viscosity(phaseIdx); }
390 {
return relativePermeability_[phaseIdx]; }
396 {
return solidState_.porosity(); }
402 {
return permeability_; }
413 phasePresentIneq(fluidState(), phaseIdx) -
414 phaseNotPresentIneq(fluidState(), phaseIdx)
423 typename FluidSystem::ParameterCache paramCache;
424 paramCache.updatePhase(fluidState_, phaseIdx);
425 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
432 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
441 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
442 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
444 return phasePresentIneq(fluidState(), phaseIdx);
445 return phaseNotPresentIneq(fluidState(), phaseIdx);
456 const unsigned int phaseIdx)
const
457 {
return fluidState.saturation(phaseIdx); }
467 const unsigned int phaseIdx)
const
471 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
472 a -= fluidState.moleFraction(phaseIdx, compIdx);
489template <
class Traits>
496 using Scalar =
typename Traits::PrimaryVariables::value_type;
497 using PermeabilityType =
typename Traits::PermeabilityType;
499 using ModelTraits =
typename Traits::ModelTraits;
500 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
502 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
503 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
504 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
506 using Indices =
typename ModelTraits::Indices;
507 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
509 using ParameterCache =
typename Traits::FluidSystem::ParameterCache;
510 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
511 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
526 static constexpr int numFluidComps = ParentType::numFluidComponents();
538 template<
class ElemSol,
class Problem,
class Element,
class Scv>
540 const Problem& problem,
541 const Element& element,
544 ParentType::update(elemSol, problem, element, scv);
546 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
549 const auto& spatialParams = problem.spatialParams();
550 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
551 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
552 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
553 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
555 typename FluidSystem::ParameterCache paramCache;
556 paramCache.updateAll(fluidState_);
560 if constexpr (enableDiffusion)
562 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
563 {
return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx); };
565 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
568 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
569 permeability_ = spatialParams.permeability(element, scv, elemSol);
570 EnergyVolVars::updateEffectiveThermalConductivity();
584 template<
class ElemSol,
class Problem,
class Element,
class Scv>
586 const Problem& problem,
587 const Element& element,
595 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
599 auto&& priVars = elemSol[scv.localDofIndex()];
601 for (
int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
602 sumSat += priVars[Indices::s0Idx + phaseIdx];
603 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
605 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
611 const auto& spatialParams = problem.spatialParams();
612 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
613 fluidState.setWettingPhase(wPhaseIdx);
614 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
615 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
622 const Scalar pw = priVars[Indices::p0Idx];
623 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
624 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
629 const Scalar pn = priVars[Indices::p0Idx];
630 for (
int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
631 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
634 DUNE_THROW(Dune::InvalidStateException,
"MPNCVolumeVariables do not support the chosen pressure formulation");
640 typename FluidSystem::ParameterCache paramCache;
641 paramCache.updateAll(fluidState);
645 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
646 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
648 updateMoleFraction(fluidState,
654 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
657 fluidState.setViscosity(phaseIdx, mu);
660 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
661 fluidState.setEnthalpy(phaseIdx, h);
674 ParameterCache & paramCache,
675 const typename Traits::PrimaryVariables& priVars)
678 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
681 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
682 actualFluidState.setMoleFraction(phaseIdx,
684 priVars[Indices::moleFrac00Idx +
685 phaseIdx*numFluidComps +
693 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
694 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
695 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
699 actualFluidState.setFugacityCoefficient(phaseIdx,
705 equilFluidState.assign(actualFluidState) ;
706 ConstraintSolver::solve(equilFluidState,
710 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
711 for (
int compIdx=0; compIdx< numFluidComps; ++ compIdx){
712 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
717 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
719 actualFluidState.setDensity(phaseIdx, rho);
721 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
733 const Scalar
xEquil(
const unsigned int phaseIdx,
const unsigned int compIdx)
const
735 return xEquil_[phaseIdx][compIdx] ;
743 {
return fluidState_; }
749 {
return solidState_; }
758 {
return fluidState_.saturation(phaseIdx); }
768 {
return fluidState_.massFraction(phaseIdx, compIdx); }
778 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
786 Scalar
molarity(
const int phaseIdx,
int compIdx)
const
787 {
return fluidState_.molarity(phaseIdx, compIdx); }
795 {
return fluidState_.molarDensity(phaseIdx);}
804 {
return fluidState_.pressure(phaseIdx); }
812 {
return fluidState_.density(phaseIdx); }
822 {
return fluidState_.temperature(0); }
828 {
return fluidState_.enthalpy(phaseIdx); }
834 {
return fluidState_.internalEnergy(phaseIdx); }
841 {
return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
847 {
return fluidState_.fugacity(compIdx); }
853 {
return fluidState_.averageMolarMass(phaseIdx); }
871 {
return fluidState_.viscosity(phaseIdx); }
880 {
return relativePermeability_[phaseIdx]; }
886 {
return solidState_.porosity(); }
892 {
return permeability_; }
903 phasePresentIneq(fluidState(), phaseIdx) -
904 phaseNotPresentIneq(fluidState(), phaseIdx)
913 typename FluidSystem::ParameterCache paramCache;
914 paramCache.updatePhase(fluidState_, phaseIdx);
915 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
922 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
931 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
932 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
934 return phasePresentIneq(fluidState(), phaseIdx);
935 return phaseNotPresentIneq(fluidState(), phaseIdx);
946 const unsigned int phaseIdx)
const
947 {
return fluidState.saturation(phaseIdx); }
957 const unsigned int phaseIdx)
const
961 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
962 a -= fluidState.moleFraction(phaseIdx, compIdx);
969 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()>
xEquil_;
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Determines the fluid composition given the component fugacities and an arbitary equation of state.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.
Enumeration of the formulations accepted by the MpNc model.
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
std::string relativePermeability(int phaseIdx) noexcept
I/O name of relative permeability for multiphase systems.
Definition: name.hh:92
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
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
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:53
static void guessInitial(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &fugVec)
Guess an initial value for the composition of the phase.
Definition: compositionfromfugacities.hh:69
static void solve(FluidState &fluidState, ParameterCache ¶mCache, int phaseIdx, const ComponentVector &targetFug)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: compositionfromfugacities.hh:97
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:59
Definition: porousmediumflow/mpnc/volumevariables.hh:42
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:421
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:83
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:328
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:362
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:482
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:264
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:310
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:249
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:350
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:395
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:371
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:255
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:389
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/mpnc/volumevariables.hh:274
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:102
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:293
typename Traits::ModelTraits::Indices Indices
Export the type encapsulating primary variable indices.
Definition: porousmediumflow/mpnc/volumevariables.hh:77
Scalar temperature(const int phaseIdx) const
Definition: porousmediumflow/mpnc/volumevariables.hh:331
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:401
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:356
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:337
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:431
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:483
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:318
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/mpnc/volumevariables.hh:284
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:478
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:81
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:486
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:301
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:439
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:88
Scalar phasePresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is present.
Definition: porousmediumflow/mpnc/volumevariables.hh:455
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:380
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:150
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:479
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:343
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:85
Scalar phaseNotPresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is not present.
Definition: porousmediumflow/mpnc/volumevariables.hh:466
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:79
Scalar porosity_
Effective porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:477
bool isPhaseActive(const unsigned int phaseIdx) const
Returns true if the fluid state is in the active set for a phase,.
Definition: porousmediumflow/mpnc/volumevariables.hh:410
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:846
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:968
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:840
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:821
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:861
void updateMoleFraction(FluidState &actualFluidState, ParameterCache ¶mCache, const typename Traits::PrimaryVariables &priVars)
Updates composition of all phases in the mutable parameters from the primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:673
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:803
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:748
Scalar phasePresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is present.
Definition: porousmediumflow/mpnc/volumevariables.hh:945
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:521
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:885
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:519
Scalar phaseNotPresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is not present.
Definition: porousmediumflow/mpnc/volumevariables.hh:956
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:517
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:827
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:811
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:921
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:973
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:515
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:911
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:891
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:794
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:870
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:833
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:969
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/mpnc/volumevariables.hh:777
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:539
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:524
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:757
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:972
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:929
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:852
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:967
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:585
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:976
bool isPhaseActive(const unsigned int phaseIdx) const
Returns true if the fluid state is in the active set for a phase,.
Definition: porousmediumflow/mpnc/volumevariables.hh:900
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/mpnc/volumevariables.hh:767
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/mpnc/volumevariables.hh:733
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:879
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:742
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:786
Definition: porousmediumflow/nonisothermal/volumevariables.hh:75
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:40
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.