3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_MPNC_VOLUME_VARIABLES_HH
27#define DUMUX_MPNC_VOLUME_VARIABLES_HH
28
31
35
37
38namespace Dumux {
39
40// forward declaration
41template <class Traits, bool enableChemicalNonEquilibrium>
49template <class Traits>
50using MPNCVolumeVariables = MPNCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium()>;
51
52template <class Traits>
55, public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
56{
59
60 using Scalar = typename Traits::PrimaryVariables::value_type;
61 using PermeabilityType = typename Traits::PermeabilityType;
62
63 using ModelTraits = typename Traits::ModelTraits;
64 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
65
66 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
67 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
68 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
69
70 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
72 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
73 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
74
75public:
77 using Indices = typename Traits::ModelTraits::Indices;
79 using FluidSystem = typename Traits::FluidSystem;
81 using FluidState = typename Traits::FluidState;
83 using SolidState = typename Traits::SolidState;
85 using SolidSystem = typename Traits::SolidSystem;
86
88 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
90 static constexpr int numFluidComps = ParentType::numFluidComponents();
91
101 template<class ElemSol, class Problem, class Element, class Scv>
102 void update(const ElemSol& elemSol,
103 const Problem& problem,
104 const Element& element,
105 const Scv& scv)
106 {
107 ParentType::update(elemSol, problem, element, scv);
108
109 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
110
111 // calculate the remaining quantities
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());
117
118 typename FluidSystem::ParameterCache paramCache;
119 paramCache.updateAll(fluidState_);
120
121 // porosity
122 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
123
124 if constexpr (enableDiffusion)
125 {
126 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
127 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
128
129 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
130 }
131
132 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
133 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
134 EnergyVolVars::updateEffectiveThermalConductivity();
135 }
136
149 template<class ElemSol, class Problem, class Element, class Scv>
150 void completeFluidState(const ElemSol& elemSol,
151 const Problem& problem,
152 const Element& element,
153 const Scv& scv,
154 FluidState& fluidState,
155 SolidState& solidState)
156 {
157
159 // set the fluid phase temperatures
161 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
162
164 // set the phase saturations
166 auto&& priVars = elemSol[scv.localDofIndex()];
167 Scalar sumSat = 0;
168 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
169 sumSat += priVars[Indices::s0Idx + phaseIdx];
170 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
171 }
172 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
173
175 // set the phase pressures
177 // capillary pressure parameters
178 const auto& spatialParams = problem.spatialParams();
179 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
180 fluidState.setWettingPhase(wPhaseIdx);
181 // capillary pressures
182
183 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
184 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
185
186 // add to the pressure of the first fluid phase
187 // depending on which pressure is stored in the primary variables
188 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
189 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
190 // For two phases this means that there is one pressure as primary variable: pw
191 const Scalar pw = priVars[Indices::p0Idx];
192 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
193 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
194 }
195 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
196 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
197 // For two phases this means that there is one pressure as primary variable: pn
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]);
201 }
202 else
203 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
204
206 // set the fluid compositions
208 typename FluidSystem::ParameterCache paramCache;
209 paramCache.updateAll(fluidState);
210
211 ComponentVector fug;
212 // retrieve component fugacities
213 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
214 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
215
216 // calculate phase compositions
217 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
218 // initial guess
219 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
220 // sanity check to make sure we have non-zero fugacity coefficients
221 if (FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx) == 0.0)
222 DUNE_THROW(NumericalProblem, "MPNCVolumeVariables do not support fluidsystems with fugacity coefficients of 0");
223
224 // set initial guess of the component's mole fraction
225 fluidState.setMoleFraction(phaseIdx, compIdx, 1.0/numFluidComps);
226 }
227
228 // calculate the phase composition from the component
229 // fugacities
230 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
231 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
232 }
233
234 // dynamic viscosities
235 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
236 // viscosities
237 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
238 fluidState.setViscosity(phaseIdx, mu);
239
240 // compute and set the enthalpy
241 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
242 fluidState.setEnthalpy(phaseIdx, h);
243 }
244 }
245
250 const FluidState &fluidState() const
251 { return fluidState_; }
252
256 const SolidState &solidState() const
257 { return solidState_; }
258
265 Scalar saturation(int phaseIdx) const
266 { return fluidState_.saturation(phaseIdx); }
267
275 Scalar massFraction(const int phaseIdx, const int compIdx) const
276 { return fluidState_.massFraction(phaseIdx, compIdx); }
277
285 Scalar moleFraction(const int phaseIdx, const int compIdx) const
286 { return fluidState_.moleFraction(phaseIdx, compIdx); }
287
294 Scalar molarity(const int phaseIdx, int compIdx) const
295 { return fluidState_.molarity(phaseIdx, compIdx); }
296
302 Scalar molarDensity(const int phaseIdx) const
303 { return fluidState_.molarDensity(phaseIdx);}
304
311 Scalar pressure(const int phaseIdx) const
312 { return fluidState_.pressure(phaseIdx); }
313
319 Scalar density(const int phaseIdx) const
320 { return fluidState_.density(phaseIdx); }
321
329 Scalar temperature() const
330 { return fluidState_.temperature(0/* phaseIdx*/); }
331
332 Scalar temperature(const int phaseIdx) const
333 { return fluidState_.temperature(phaseIdx); }
334
338 Scalar enthalpy(const int phaseIdx) const
339 { return fluidState_.enthalpy(phaseIdx); }
340
344 Scalar internalEnergy(const int phaseIdx) const
345 { return fluidState_.internalEnergy(phaseIdx); }
346
351 Scalar fluidThermalConductivity(const int phaseIdx) const
352 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
353
357 Scalar fugacity(const int compIdx) const
358 { return fluidState_.fugacity(compIdx); }
359
363 Scalar averageMolarMass(const int phaseIdx) const
364 { return fluidState_.averageMolarMass(phaseIdx); }
365
372 Scalar mobility(const unsigned int phaseIdx) const
373 {
374 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
375 }
376
381 Scalar viscosity(const unsigned int phaseIdx) const
382 { return fluidState_.viscosity(phaseIdx); }
383
390 Scalar relativePermeability(const unsigned int phaseIdx) const
391 { return relativePermeability_[phaseIdx]; }
392
396 Scalar porosity() const
397 { return solidState_.porosity(); }
398
402 const PermeabilityType& permeability() const
403 { return permeability_; }
404
411 bool isPhaseActive(const unsigned int phaseIdx) const
412 {
413 return
414 phasePresentIneq(fluidState(), phaseIdx) -
415 phaseNotPresentIneq(fluidState(), phaseIdx)
416 >= 0;
417 }
418
422 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
423 {
424 typename FluidSystem::ParameterCache paramCache;
425 paramCache.updatePhase(fluidState_, phaseIdx);
426 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
427 }
428
432 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
433 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
434
440 Scalar phaseNcp(const unsigned int phaseIdx) const
441 {
442 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
443 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
444 if (aEval > bEval)
445 return phasePresentIneq(fluidState(), phaseIdx);
446 return phaseNotPresentIneq(fluidState(), phaseIdx);
447 };
448
456 Scalar phasePresentIneq(const FluidState &fluidState,
457 const unsigned int phaseIdx) const
458 { return fluidState.saturation(phaseIdx); }
459
467 Scalar phaseNotPresentIneq(const FluidState &fluidState,
468 const unsigned int phaseIdx) const
469 {
470 // difference of sum of mole fractions in the phase from 100%
471 Scalar a = 1;
472 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
473 a -= fluidState.moleFraction(phaseIdx, compIdx);
474 return a;
475 }
476
477protected:
478 Scalar porosity_;
479 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
480 PermeabilityType permeability_;
481
485
486 // Effective diffusion coefficients for the phases
487 DiffusionCoefficients effectiveDiffCoeff_;
488};
489
490template <class Traits>
492 : public PorousMediumFlowVolumeVariables<Traits>
493 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
494{
497 using Scalar = typename Traits::PrimaryVariables::value_type;
498 using PermeabilityType = typename Traits::PermeabilityType;
499
500 using ModelTraits = typename Traits::ModelTraits;
501 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
502
503 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
504 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
505 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
506
507 using Indices = typename ModelTraits::Indices;
508 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
510 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
511 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
512 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
513
514public:
516 using FluidSystem = typename Traits::FluidSystem;
518 using FluidState = typename Traits::FluidState;
520 using SolidState = typename Traits::SolidState;
522 using SolidSystem = typename Traits::SolidSystem;
523
525 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
527 static constexpr int numFluidComps = ParentType::numFluidComponents();
529
539 template<class ElemSol, class Problem, class Element, class Scv>
540 void update(const ElemSol& elemSol,
541 const Problem& problem,
542 const Element& element,
543 const Scv& scv)
544 {
545 ParentType::update(elemSol, problem, element, scv);
546
547 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
548
549 // calculate the remaining quantities
550 const auto& spatialParams = problem.spatialParams();
551 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
552 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
553 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
554 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
555
556 typename FluidSystem::ParameterCache paramCache;
557 paramCache.updateAll(fluidState_);
558
559 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
560
561 if constexpr (enableDiffusion)
562 {
563 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
564 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
565
566 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
567 }
568
569 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
570 permeability_ = spatialParams.permeability(element, scv, elemSol);
571 EnergyVolVars::updateEffectiveThermalConductivity();
572 }
573
585 template<class ElemSol, class Problem, class Element, class Scv>
586 void completeFluidState(const ElemSol& elemSol,
587 const Problem& problem,
588 const Element& element,
589 const Scv& scv,
590 FluidState& fluidState,
591 SolidState& solidState)
592 {
594 // set the fluid phase temperatures
596 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
598 // set the phase saturations
600 auto&& priVars = elemSol[scv.localDofIndex()];
601 Scalar sumSat = 0;
602 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
603 sumSat += priVars[Indices::s0Idx + phaseIdx];
604 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
605 }
606 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
607
609 // set the phase pressures
611 // capillary pressure parameters
612 const auto& spatialParams = problem.spatialParams();
613 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
614 fluidState.setWettingPhase(wPhaseIdx);
615 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
616 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
617
618 // add to the pressure of the first fluid phase
619 // depending on which pressure is stored in the primary variables
620 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
621 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
622 // For two phases this means that there is one pressure as primary variable: pw
623 const Scalar pw = priVars[Indices::p0Idx];
624 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
625 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
626 }
627 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
628 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
629 // For two phases this means that there is one pressure as primary variable: pn
630 const Scalar pn = priVars[Indices::p0Idx];
631 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
632 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
633 }
634 else
635 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
636
637
639 // set the fluid compositions
641 typename FluidSystem::ParameterCache paramCache;
642 paramCache.updateAll(fluidState);
643
644 updateMoleFraction(fluidState, paramCache, priVars);
645
646 // dynamic viscosities
647 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
648 // viscosities
649 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
650 fluidState.setViscosity(phaseIdx, mu);
651
652 // compute and set the enthalpy
653 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
654 fluidState.setEnthalpy(phaseIdx, h);
655 }
656 }
657
666 void updateMoleFraction(FluidState & actualFluidState,
667 ParameterCache & paramCache,
668 const typename Traits::PrimaryVariables& priVars)
669 {
670 // setting the mole fractions of the fluid state
671 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
672 {
673 // set the component mole fractions
674 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
675 actualFluidState.setMoleFraction(phaseIdx,
676 compIdx,
677 priVars[Indices::moleFrac00Idx +
678 phaseIdx*numFluidComps +
679 compIdx]);
680 }
681 }
682
683 // TODO: Can this be removed, are fugacity coefficients used in the constraint solver?
684 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
685 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
686 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
687 paramCache,
688 phaseIdx,
689 compIdx);
690 actualFluidState.setFugacityCoefficient(phaseIdx,
691 compIdx,
692 phi);
693 }
694
695 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
696 equilFluidState.assign(actualFluidState) ;
697 ConstraintSolver::solve(equilFluidState,
698 paramCache) ;
699
700 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
701 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
702 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
703 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
704 }
705 }
706
707 // compute densities of all phases
708 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
709 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
710 actualFluidState.setDensity(phaseIdx, rho);
711 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
712 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
713 }
714
715 }
716
724 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
725 {
726 return xEquil_[phaseIdx][compIdx] ;
727 }
728
733 const FluidState &fluidState() const
734 { return fluidState_; }
735
739 const SolidState &solidState() const
740 { return solidState_; }
741
748 Scalar saturation(int phaseIdx) const
749 { return fluidState_.saturation(phaseIdx); }
750
758 Scalar massFraction(const int phaseIdx, const int compIdx) const
759 { return fluidState_.massFraction(phaseIdx, compIdx); }
760
768 Scalar moleFraction(const int phaseIdx, const int compIdx) const
769 { return fluidState_.moleFraction(phaseIdx, compIdx); }
770
777 Scalar molarity(const int phaseIdx, int compIdx) const
778 { return fluidState_.molarity(phaseIdx, compIdx); }
779
785 Scalar molarDensity(const int phaseIdx) const
786 { return fluidState_.molarDensity(phaseIdx);}
787
794 Scalar pressure(const int phaseIdx) const
795 { return fluidState_.pressure(phaseIdx); }
796
802 Scalar density(const int phaseIdx) const
803 { return fluidState_.density(phaseIdx); }
804
812 Scalar temperature() const
813 { return fluidState_.temperature(0/* phaseIdx*/); }
814
818 Scalar enthalpy(const int phaseIdx) const
819 { return fluidState_.enthalpy(phaseIdx); }
820
824 Scalar internalEnergy(const int phaseIdx) const
825 { return fluidState_.internalEnergy(phaseIdx); }
826
831 Scalar fluidThermalConductivity(const int phaseIdx) const
832 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
833
837 Scalar fugacity(const int compIdx) const
838 { return fluidState_.fugacity(compIdx); }
839
843 Scalar averageMolarMass(const int phaseIdx) const
844 { return fluidState_.averageMolarMass(phaseIdx); }
845
852 Scalar mobility(const unsigned int phaseIdx) const
853 {
854 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
855 }
856
861 Scalar viscosity(const unsigned int phaseIdx) const
862 { return fluidState_.viscosity(phaseIdx); }
863
870 Scalar relativePermeability(const unsigned int phaseIdx) const
871 { return relativePermeability_[phaseIdx]; }
872
876 Scalar porosity() const
877 { return solidState_.porosity(); }
878
882 const PermeabilityType& permeability() const
883 { return permeability_; }
884
891 bool isPhaseActive(const unsigned int phaseIdx) const
892 {
893 return
894 phasePresentIneq(fluidState(), phaseIdx) -
895 phaseNotPresentIneq(fluidState(), phaseIdx)
896 >= 0;
897 }
898
902 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
903 {
904 typename FluidSystem::ParameterCache paramCache;
905 paramCache.updatePhase(fluidState_, phaseIdx);
906 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
907 }
908
912 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
913 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
914
920 Scalar phaseNcp(const unsigned int phaseIdx) const
921 {
922 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
923 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
924 if (aEval > bEval)
925 return phasePresentIneq(fluidState(), phaseIdx);
926 return phaseNotPresentIneq(fluidState(), phaseIdx);
927 };
928
936 Scalar phasePresentIneq(const FluidState &fluidState,
937 const unsigned int phaseIdx) const
938 { return fluidState.saturation(phaseIdx); }
939
947 Scalar phaseNotPresentIneq(const FluidState &fluidState,
948 const unsigned int phaseIdx) const
949 {
950 // difference of sum of mole fractions in the phase from 100%
951 Scalar a = 1;
952 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
953 a -= fluidState.moleFraction(phaseIdx, compIdx);
954 return a;
955 }
956
957protected:
958 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
959 PermeabilityType permeability_;
960 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
961
965
966 // Effective diffusion coefficients for the phases
967 DiffusionCoefficients effectiveDiffCoeff_;
968};
969
970} // end namespace
971
972#endif
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 arbitrary 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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:53
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:68
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:96
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:422
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:329
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:363
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:483
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:265
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:311
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:250
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:351
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:396
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:372
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:256
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:390
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:275
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:294
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:332
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:402
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:357
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:338
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:432
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:484
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:319
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:285
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:479
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:81
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:487
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:302
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:440
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:456
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:381
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:480
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:344
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:467
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:478
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:411
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:837
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:959
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:831
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:812
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:852
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:666
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:794
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:739
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:936
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:522
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:876
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:520
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:947
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:518
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:818
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:802
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:912
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:964
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:516
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:902
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:882
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:785
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:861
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:824
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:960
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:768
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:540
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:525
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:748
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:963
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:920
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:843
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:958
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:586
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:967
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:891
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:758
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:724
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:870
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:733
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:777
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.