14#ifndef DUMUX_MPNC_VOLUME_VARIABLES_HH
15#define DUMUX_MPNC_VOLUME_VARIABLES_HH
29template <
class Traits,
bool enableChemicalNonEquilibrium>
37template <
class Traits>
40template <
class Traits>
48 using Scalar =
typename Traits::PrimaryVariables::value_type;
49 using PermeabilityType =
typename Traits::PermeabilityType;
51 using ModelTraits =
typename Traits::ModelTraits;
52 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
54 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
55 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
56 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
58 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
60 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
61 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
65 using Indices =
typename Traits::ModelTraits::Indices;
76 static constexpr int numFluidPhases() {
return ModelTraits::numFluidPhases(); }
78 static constexpr int numFluidComps = ParentType::numFluidComponents();
89 template<
class ElemSol,
class Problem,
class Element,
class Scv>
91 const Problem& problem,
92 const Element& element,
95 ParentType::update(elemSol, problem, element, scv);
97 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
100 const auto& spatialParams = problem.spatialParams();
101 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
102 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
103 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
104 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
106 typename FluidSystem::ParameterCache paramCache;
107 paramCache.updateAll(fluidState_);
112 if constexpr (enableDiffusion)
114 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
115 {
return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx); };
117 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
120 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
121 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
122 EnergyVolVars::updateEffectiveThermalConductivity();
137 template<
class ElemSol,
class Problem,
class Element,
class Scv>
139 const Problem& problem,
140 const Element& element,
149 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
154 auto&& priVars = elemSol[scv.localDofIndex()];
156 for (
int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
157 sumSat += priVars[Indices::s0Idx + phaseIdx];
158 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
160 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
166 const auto& spatialParams = problem.spatialParams();
167 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
168 fluidState.setWettingPhase(wPhaseIdx);
171 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
172 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
179 const Scalar pw = priVars[Indices::p0Idx];
180 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
181 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
186 const Scalar pn = priVars[Indices::p0Idx];
187 for (
int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
188 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
191 DUNE_THROW(Dune::InvalidStateException,
"MPNCVolumeVariables do not support the chosen pressure formulation");
196 typename FluidSystem::ParameterCache paramCache;
197 paramCache.updateAll(fluidState);
201 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
202 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
205 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
207 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
209 if (FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx) == 0.0)
210 DUNE_THROW(
NumericalProblem,
"MPNCVolumeVariables do not support fluidsystems with fugacity coefficients of 0");
213 fluidState.setMoleFraction(phaseIdx, compIdx, 1.0/numFluidComps);
223 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
226 fluidState.setViscosity(phaseIdx, mu);
229 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
230 fluidState.setEnthalpy(phaseIdx, h);
239 {
return fluidState_; }
245 {
return solidState_; }
254 {
return fluidState_.saturation(phaseIdx); }
264 {
return fluidState_.massFraction(phaseIdx, compIdx); }
274 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
282 Scalar
molarity(
const int phaseIdx,
int compIdx)
const
283 {
return fluidState_.molarity(phaseIdx, compIdx); }
291 {
return fluidState_.molarDensity(phaseIdx);}
300 {
return fluidState_.pressure(phaseIdx); }
308 {
return fluidState_.density(phaseIdx); }
318 {
return fluidState_.temperature(0); }
321 {
return fluidState_.temperature(phaseIdx); }
327 {
return fluidState_.enthalpy(phaseIdx); }
333 {
return fluidState_.internalEnergy(phaseIdx); }
340 {
return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
346 {
return fluidState_.fugacity(compIdx); }
352 {
return fluidState_.averageMolarMass(phaseIdx); }
370 {
return fluidState_.viscosity(phaseIdx); }
379 {
return relativePermeability_[phaseIdx]; }
385 {
return solidState_.porosity(); }
391 {
return permeability_; }
402 phasePresentIneq(fluidState(), phaseIdx) -
403 phaseNotPresentIneq(fluidState(), phaseIdx)
412 typename FluidSystem::ParameterCache paramCache;
413 paramCache.updatePhase(fluidState_, phaseIdx);
414 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
421 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
430 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
431 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
433 return phasePresentIneq(fluidState(), phaseIdx);
434 return phaseNotPresentIneq(fluidState(), phaseIdx);
445 const unsigned int phaseIdx)
const
446 {
return fluidState.saturation(phaseIdx); }
456 const unsigned int phaseIdx)
const
460 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
461 a -= fluidState.moleFraction(phaseIdx, compIdx);
478template <
class Traits>
485 using Scalar =
typename Traits::PrimaryVariables::value_type;
486 using PermeabilityType =
typename Traits::PermeabilityType;
488 using ModelTraits =
typename Traits::ModelTraits;
489 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
491 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
492 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
493 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
495 using Indices =
typename ModelTraits::Indices;
496 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
498 using ParameterCache =
typename Traits::FluidSystem::ParameterCache;
499 using EffDiffModel =
typename Traits::EffectiveDiffusivityModel;
500 using DiffusionCoefficients =
typename Traits::DiffusionType::DiffusionCoefficientsContainer;
515 static constexpr int numFluidComps = ParentType::numFluidComponents();
527 template<
class ElemSol,
class Problem,
class Element,
class Scv>
529 const Problem& problem,
530 const Element& element,
533 ParentType::update(elemSol, problem, element, scv);
535 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
538 const auto& spatialParams = problem.spatialParams();
539 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
540 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
541 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
542 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
544 typename FluidSystem::ParameterCache paramCache;
545 paramCache.updateAll(fluidState_);
549 if constexpr (enableDiffusion)
551 auto getEffectiveDiffusionCoefficient = [&](
int phaseIdx,
int compIIdx,
int compJIdx)
552 {
return EffDiffModel::effectiveDiffusionCoefficient(*
this, phaseIdx, compIIdx, compJIdx); };
554 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
557 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
558 permeability_ = spatialParams.permeability(element, scv, elemSol);
559 EnergyVolVars::updateEffectiveThermalConductivity();
573 template<
class ElemSol,
class Problem,
class Element,
class Scv>
575 const Problem& problem,
576 const Element& element,
584 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
588 auto&& priVars = elemSol[scv.localDofIndex()];
590 for (
int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
591 sumSat += priVars[Indices::s0Idx + phaseIdx];
592 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
594 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
600 const auto& spatialParams = problem.spatialParams();
601 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
602 fluidState.setWettingPhase(wPhaseIdx);
603 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
604 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
611 const Scalar pw = priVars[Indices::p0Idx];
612 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
613 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
618 const Scalar pn = priVars[Indices::p0Idx];
619 for (
int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
620 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
623 DUNE_THROW(Dune::InvalidStateException,
"MPNCVolumeVariables do not support the chosen pressure formulation");
629 typename FluidSystem::ParameterCache paramCache;
630 paramCache.updateAll(fluidState);
632 updateMoleFraction(fluidState, paramCache, priVars);
635 for (
int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
638 fluidState.setViscosity(phaseIdx, mu);
641 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
642 fluidState.setEnthalpy(phaseIdx, h);
655 ParameterCache & paramCache,
656 const typename Traits::PrimaryVariables& priVars)
659 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
662 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
663 actualFluidState.setMoleFraction(phaseIdx,
665 priVars[Indices::moleFrac00Idx +
666 phaseIdx*numFluidComps +
672 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
673 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
674 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
678 actualFluidState.setFugacityCoefficient(phaseIdx,
684 equilFluidState.assign(actualFluidState) ;
685 ConstraintSolver::solve(equilFluidState,
689 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
690 for (
int compIdx=0; compIdx< numFluidComps; ++ compIdx){
691 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
696 for(
int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
698 actualFluidState.setDensity(phaseIdx, rho);
700 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
712 const Scalar
xEquil(
const unsigned int phaseIdx,
const unsigned int compIdx)
const
714 return xEquil_[phaseIdx][compIdx] ;
722 {
return fluidState_; }
728 {
return solidState_; }
737 {
return fluidState_.saturation(phaseIdx); }
747 {
return fluidState_.massFraction(phaseIdx, compIdx); }
757 {
return fluidState_.moleFraction(phaseIdx, compIdx); }
765 Scalar
molarity(
const int phaseIdx,
int compIdx)
const
766 {
return fluidState_.molarity(phaseIdx, compIdx); }
774 {
return fluidState_.molarDensity(phaseIdx);}
783 {
return fluidState_.pressure(phaseIdx); }
791 {
return fluidState_.density(phaseIdx); }
801 {
return fluidState_.temperature(0); }
807 {
return fluidState_.enthalpy(phaseIdx); }
813 {
return fluidState_.internalEnergy(phaseIdx); }
820 {
return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
826 {
return fluidState_.fugacity(compIdx); }
832 {
return fluidState_.averageMolarMass(phaseIdx); }
850 {
return fluidState_.viscosity(phaseIdx); }
859 {
return relativePermeability_[phaseIdx]; }
865 {
return solidState_.porosity(); }
871 {
return permeability_; }
882 phasePresentIneq(fluidState(), phaseIdx) -
883 phaseNotPresentIneq(fluidState(), phaseIdx)
892 typename FluidSystem::ParameterCache paramCache;
893 paramCache.updatePhase(fluidState_, phaseIdx);
894 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
901 {
return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
910 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
911 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
913 return phasePresentIneq(fluidState(), phaseIdx);
914 return phaseNotPresentIneq(fluidState(), phaseIdx);
925 const unsigned int phaseIdx)
const
926 {
return fluidState.saturation(phaseIdx); }
936 const unsigned int phaseIdx)
const
940 for (
int compIdx = 0; compIdx < numFluidComps; ++compIdx)
941 a -= fluidState.moleFraction(phaseIdx, compIdx);
948 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()>
xEquil_;
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:41
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:56
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:84
Definition: porousmediumflow/nonisothermal/volumevariables.hh:63
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:410
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:71
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:317
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:351
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:471
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:253
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:299
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:238
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:339
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:384
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:360
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:244
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:378
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:263
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:90
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:282
typename Traits::ModelTraits::Indices Indices
Export the type encapsulating primary variable indices.
Definition: porousmediumflow/mpnc/volumevariables.hh:65
Scalar temperature(const int phaseIdx) const
Definition: porousmediumflow/mpnc/volumevariables.hh:320
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:390
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:345
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:326
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:420
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:472
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:307
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:273
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:467
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:69
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:475
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:290
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:428
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:76
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:444
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:369
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:138
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:468
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:332
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:73
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:455
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:67
Scalar porosity_
Effective porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:466
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:399
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:825
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:947
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:819
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:800
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:840
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:654
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:782
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:727
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:924
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:510
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:864
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:508
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:935
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:506
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:806
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:790
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:900
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:952
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:504
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:890
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:870
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:773
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:849
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:812
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:948
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:756
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:528
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:513
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:736
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:951
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:908
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:831
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:946
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:574
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:955
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:879
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:746
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:712
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:858
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:721
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:765
Definition: porousmediumflow/mpnc/volumevariables.hh:30
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:47
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:28
Determines the fluid composition given the component fugacities and an arbitrary equation of state.
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
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
std::string relativePermeability(int phaseIdx) noexcept
I/O name of relative permeability for multiphase systems.
Definition: name.hh:80
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.