version 3.8
porousmediumflow/mpnc/volumevariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_MPNC_VOLUME_VARIABLES_HH
15#define DUMUX_MPNC_VOLUME_VARIABLES_HH
16
19
23
25
26namespace Dumux {
27
28// forward declaration
29template <class Traits, bool enableChemicalNonEquilibrium>
37template <class Traits>
38using MPNCVolumeVariables = MPNCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium()>;
39
40template <class Traits>
43, public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
44{
47
48 using Scalar = typename Traits::PrimaryVariables::value_type;
49 using PermeabilityType = typename Traits::PermeabilityType;
50
51 using ModelTraits = typename Traits::ModelTraits;
52 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
53
54 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
55 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
56 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
57
58 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
60 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
61 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
62
63public:
65 using Indices = typename Traits::ModelTraits::Indices;
67 using FluidSystem = typename Traits::FluidSystem;
69 using FluidState = typename Traits::FluidState;
71 using SolidState = typename Traits::SolidState;
73 using SolidSystem = typename Traits::SolidSystem;
74
76 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
78 static constexpr int numFluidComps = ParentType::numFluidComponents();
79
89 template<class ElemSol, class Problem, class Element, class Scv>
90 void update(const ElemSol& elemSol,
91 const Problem& problem,
92 const Element& element,
93 const Scv& scv)
94 {
95 ParentType::update(elemSol, problem, element, scv);
96
97 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
98
99 // calculate the remaining quantities
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());
105
106 typename FluidSystem::ParameterCache paramCache;
107 paramCache.updateAll(fluidState_);
108
109 // porosity
110 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
111
112 if constexpr (enableDiffusion)
113 {
114 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
115 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
116
117 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
118 }
119
120 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
121 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
122 EnergyVolVars::updateEffectiveThermalConductivity();
123 }
124
137 template<class ElemSol, class Problem, class Element, class Scv>
138 void completeFluidState(const ElemSol& elemSol,
139 const Problem& problem,
140 const Element& element,
141 const Scv& scv,
142 FluidState& fluidState,
143 SolidState& solidState)
144 {
145
147 // set the fluid phase temperatures
149 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
150
152 // set the phase saturations
154 auto&& priVars = elemSol[scv.localDofIndex()];
155 Scalar sumSat = 0;
156 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
157 sumSat += priVars[Indices::s0Idx + phaseIdx];
158 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
159 }
160 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
161
163 // set the phase pressures
165 // capillary pressure parameters
166 const auto& spatialParams = problem.spatialParams();
167 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
168 fluidState.setWettingPhase(wPhaseIdx);
169 // capillary pressures
170
171 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
172 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
173
174 // add to the pressure of the first fluid phase
175 // depending on which pressure is stored in the primary variables
176 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
177 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
178 // For two phases this means that there is one pressure as primary variable: pw
179 const Scalar pw = priVars[Indices::p0Idx];
180 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
181 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
182 }
183 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
184 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
185 // For two phases this means that there is one pressure as primary variable: pn
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]);
189 }
190 else
191 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
192
194 // set the fluid compositions
196 typename FluidSystem::ParameterCache paramCache;
197 paramCache.updateAll(fluidState);
198
199 ComponentVector fug;
200 // retrieve component fugacities
201 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
202 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
203
204 // calculate phase compositions
205 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
206 // initial guess
207 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
208 // sanity check to make sure we have non-zero fugacity coefficients
209 if (FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx) == 0.0)
210 DUNE_THROW(NumericalProblem, "MPNCVolumeVariables do not support fluidsystems with fugacity coefficients of 0");
211
212 // set initial guess of the component's mole fraction
213 fluidState.setMoleFraction(phaseIdx, compIdx, 1.0/numFluidComps);
214 }
215
216 // calculate the phase composition from the component
217 // fugacities
218 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
219 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
220 }
221
222 // dynamic viscosities
223 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
224 // viscosities
225 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
226 fluidState.setViscosity(phaseIdx, mu);
227
228 // compute and set the enthalpy
229 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
230 fluidState.setEnthalpy(phaseIdx, h);
231 }
232 }
233
238 const FluidState &fluidState() const
239 { return fluidState_; }
240
244 const SolidState &solidState() const
245 { return solidState_; }
246
253 Scalar saturation(int phaseIdx) const
254 { return fluidState_.saturation(phaseIdx); }
255
263 Scalar massFraction(const int phaseIdx, const int compIdx) const
264 { return fluidState_.massFraction(phaseIdx, compIdx); }
265
273 Scalar moleFraction(const int phaseIdx, const int compIdx) const
274 { return fluidState_.moleFraction(phaseIdx, compIdx); }
275
282 Scalar molarity(const int phaseIdx, int compIdx) const
283 { return fluidState_.molarity(phaseIdx, compIdx); }
284
290 Scalar molarDensity(const int phaseIdx) const
291 { return fluidState_.molarDensity(phaseIdx);}
292
299 Scalar pressure(const int phaseIdx) const
300 { return fluidState_.pressure(phaseIdx); }
301
307 Scalar density(const int phaseIdx) const
308 { return fluidState_.density(phaseIdx); }
309
317 Scalar temperature() const
318 { return fluidState_.temperature(0/* phaseIdx*/); }
319
320 Scalar temperature(const int phaseIdx) const
321 { return fluidState_.temperature(phaseIdx); }
322
326 Scalar enthalpy(const int phaseIdx) const
327 { return fluidState_.enthalpy(phaseIdx); }
328
332 Scalar internalEnergy(const int phaseIdx) const
333 { return fluidState_.internalEnergy(phaseIdx); }
334
339 Scalar fluidThermalConductivity(const int phaseIdx) const
340 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
341
345 Scalar fugacity(const int compIdx) const
346 { return fluidState_.fugacity(compIdx); }
347
351 Scalar averageMolarMass(const int phaseIdx) const
352 { return fluidState_.averageMolarMass(phaseIdx); }
353
360 Scalar mobility(const unsigned int phaseIdx) const
361 {
362 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
363 }
364
369 Scalar viscosity(const unsigned int phaseIdx) const
370 { return fluidState_.viscosity(phaseIdx); }
371
378 Scalar relativePermeability(const unsigned int phaseIdx) const
379 { return relativePermeability_[phaseIdx]; }
380
384 Scalar porosity() const
385 { return solidState_.porosity(); }
386
390 const PermeabilityType& permeability() const
391 { return permeability_; }
392
399 bool isPhaseActive(const unsigned int phaseIdx) const
400 {
401 return
402 phasePresentIneq(fluidState(), phaseIdx) -
403 phaseNotPresentIneq(fluidState(), phaseIdx)
404 >= 0;
405 }
406
410 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
411 {
412 typename FluidSystem::ParameterCache paramCache;
413 paramCache.updatePhase(fluidState_, phaseIdx);
414 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
415 }
416
420 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
421 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
422
428 Scalar phaseNcp(const unsigned int phaseIdx) const
429 {
430 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
431 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
432 if (aEval > bEval)
433 return phasePresentIneq(fluidState(), phaseIdx);
434 return phaseNotPresentIneq(fluidState(), phaseIdx);
435 };
436
444 Scalar phasePresentIneq(const FluidState &fluidState,
445 const unsigned int phaseIdx) const
446 { return fluidState.saturation(phaseIdx); }
447
455 Scalar phaseNotPresentIneq(const FluidState &fluidState,
456 const unsigned int phaseIdx) const
457 {
458 // difference of sum of mole fractions in the phase from 100%
459 Scalar a = 1;
460 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
461 a -= fluidState.moleFraction(phaseIdx, compIdx);
462 return a;
463 }
464
465protected:
466 Scalar porosity_;
467 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
468 PermeabilityType permeability_;
469
473
474 // Effective diffusion coefficients for the phases
475 DiffusionCoefficients effectiveDiffCoeff_;
476};
477
478template <class Traits>
480 : public PorousMediumFlowVolumeVariables<Traits>
481 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
482{
485 using Scalar = typename Traits::PrimaryVariables::value_type;
486 using PermeabilityType = typename Traits::PermeabilityType;
487
488 using ModelTraits = typename Traits::ModelTraits;
489 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
490
491 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
492 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
493 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
494
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;
501
502public:
504 using FluidSystem = typename Traits::FluidSystem;
506 using FluidState = typename Traits::FluidState;
508 using SolidState = typename Traits::SolidState;
510 using SolidSystem = typename Traits::SolidSystem;
511
513 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
515 static constexpr int numFluidComps = ParentType::numFluidComponents();
517
527 template<class ElemSol, class Problem, class Element, class Scv>
528 void update(const ElemSol& elemSol,
529 const Problem& problem,
530 const Element& element,
531 const Scv& scv)
532 {
533 ParentType::update(elemSol, problem, element, scv);
534
535 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
536
537 // calculate the remaining quantities
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());
543
544 typename FluidSystem::ParameterCache paramCache;
545 paramCache.updateAll(fluidState_);
546
547 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
548
549 if constexpr (enableDiffusion)
550 {
551 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
552 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
553
554 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
555 }
556
557 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
558 permeability_ = spatialParams.permeability(element, scv, elemSol);
559 EnergyVolVars::updateEffectiveThermalConductivity();
560 }
561
573 template<class ElemSol, class Problem, class Element, class Scv>
574 void completeFluidState(const ElemSol& elemSol,
575 const Problem& problem,
576 const Element& element,
577 const Scv& scv,
578 FluidState& fluidState,
579 SolidState& solidState)
580 {
582 // set the fluid phase temperatures
584 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
586 // set the phase saturations
588 auto&& priVars = elemSol[scv.localDofIndex()];
589 Scalar sumSat = 0;
590 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
591 sumSat += priVars[Indices::s0Idx + phaseIdx];
592 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
593 }
594 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
595
597 // set the phase pressures
599 // capillary pressure parameters
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);
605
606 // add to the pressure of the first fluid phase
607 // depending on which pressure is stored in the primary variables
608 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
609 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
610 // For two phases this means that there is one pressure as primary variable: pw
611 const Scalar pw = priVars[Indices::p0Idx];
612 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
613 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
614 }
615 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
616 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
617 // For two phases this means that there is one pressure as primary variable: pn
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]);
621 }
622 else
623 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
624
625
627 // set the fluid compositions
629 typename FluidSystem::ParameterCache paramCache;
630 paramCache.updateAll(fluidState);
631
632 updateMoleFraction(fluidState, paramCache, priVars);
633
634 // dynamic viscosities
635 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
636 // viscosities
637 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
638 fluidState.setViscosity(phaseIdx, mu);
639
640 // compute and set the enthalpy
641 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
642 fluidState.setEnthalpy(phaseIdx, h);
643 }
644 }
645
654 void updateMoleFraction(FluidState & actualFluidState,
655 ParameterCache & paramCache,
656 const typename Traits::PrimaryVariables& priVars)
657 {
658 // setting the mole fractions of the fluid state
659 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
660 {
661 // set the component mole fractions
662 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
663 actualFluidState.setMoleFraction(phaseIdx,
664 compIdx,
665 priVars[Indices::moleFrac00Idx +
666 phaseIdx*numFluidComps +
667 compIdx]);
668 }
669 }
670
671 // TODO: Can this be removed, are fugacity coefficients used in the constraint solver?
672 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
673 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
674 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
675 paramCache,
676 phaseIdx,
677 compIdx);
678 actualFluidState.setFugacityCoefficient(phaseIdx,
679 compIdx,
680 phi);
681 }
682
683 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
684 equilFluidState.assign(actualFluidState) ;
685 ConstraintSolver::solve(equilFluidState,
686 paramCache) ;
687
688 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
689 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
690 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
691 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
692 }
693 }
694
695 // compute densities of all phases
696 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
697 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
698 actualFluidState.setDensity(phaseIdx, rho);
699 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
700 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
701 }
702
703 }
704
712 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
713 {
714 return xEquil_[phaseIdx][compIdx] ;
715 }
716
721 const FluidState &fluidState() const
722 { return fluidState_; }
723
727 const SolidState &solidState() const
728 { return solidState_; }
729
736 Scalar saturation(int phaseIdx) const
737 { return fluidState_.saturation(phaseIdx); }
738
746 Scalar massFraction(const int phaseIdx, const int compIdx) const
747 { return fluidState_.massFraction(phaseIdx, compIdx); }
748
756 Scalar moleFraction(const int phaseIdx, const int compIdx) const
757 { return fluidState_.moleFraction(phaseIdx, compIdx); }
758
765 Scalar molarity(const int phaseIdx, int compIdx) const
766 { return fluidState_.molarity(phaseIdx, compIdx); }
767
773 Scalar molarDensity(const int phaseIdx) const
774 { return fluidState_.molarDensity(phaseIdx);}
775
782 Scalar pressure(const int phaseIdx) const
783 { return fluidState_.pressure(phaseIdx); }
784
790 Scalar density(const int phaseIdx) const
791 { return fluidState_.density(phaseIdx); }
792
800 Scalar temperature() const
801 { return fluidState_.temperature(0/* phaseIdx*/); }
802
806 Scalar enthalpy(const int phaseIdx) const
807 { return fluidState_.enthalpy(phaseIdx); }
808
812 Scalar internalEnergy(const int phaseIdx) const
813 { return fluidState_.internalEnergy(phaseIdx); }
814
819 Scalar fluidThermalConductivity(const int phaseIdx) const
820 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
821
825 Scalar fugacity(const int compIdx) const
826 { return fluidState_.fugacity(compIdx); }
827
831 Scalar averageMolarMass(const int phaseIdx) const
832 { return fluidState_.averageMolarMass(phaseIdx); }
833
840 Scalar mobility(const unsigned int phaseIdx) const
841 {
842 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
843 }
844
849 Scalar viscosity(const unsigned int phaseIdx) const
850 { return fluidState_.viscosity(phaseIdx); }
851
858 Scalar relativePermeability(const unsigned int phaseIdx) const
859 { return relativePermeability_[phaseIdx]; }
860
864 Scalar porosity() const
865 { return solidState_.porosity(); }
866
870 const PermeabilityType& permeability() const
871 { return permeability_; }
872
879 bool isPhaseActive(const unsigned int phaseIdx) const
880 {
881 return
882 phasePresentIneq(fluidState(), phaseIdx) -
883 phaseNotPresentIneq(fluidState(), phaseIdx)
884 >= 0;
885 }
886
890 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
891 {
892 typename FluidSystem::ParameterCache paramCache;
893 paramCache.updatePhase(fluidState_, phaseIdx);
894 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
895 }
896
900 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
901 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
902
908 Scalar phaseNcp(const unsigned int phaseIdx) const
909 {
910 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
911 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
912 if (aEval > bEval)
913 return phasePresentIneq(fluidState(), phaseIdx);
914 return phaseNotPresentIneq(fluidState(), phaseIdx);
915 };
916
924 Scalar phasePresentIneq(const FluidState &fluidState,
925 const unsigned int phaseIdx) const
926 { return fluidState.saturation(phaseIdx); }
927
935 Scalar phaseNotPresentIneq(const FluidState &fluidState,
936 const unsigned int phaseIdx) const
937 {
938 // difference of sum of mole fractions in the phase from 100%
939 Scalar a = 1;
940 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
941 a -= fluidState.moleFraction(phaseIdx, compIdx);
942 return a;
943 }
944
945protected:
946 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
947 PermeabilityType permeability_;
948 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
949
953
954 // Effective diffusion coefficients for the phases
955 DiffusionCoefficients effectiveDiffCoeff_;
956};
957
958} // end namespace
959
960#endif
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:41
static void guessInitial(FluidState &fluidState, ParameterCache &paramCache, 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 &paramCache, 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 &paramCache, 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
Definition: adapt.hh:17
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.
Enumeration of the formulations accepted by the MpNc model.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.