Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 Scalar x_ij = 1.0/numFluidComps;
221
222 // set initial guess of the component's mole fraction
223 fluidState.setMoleFraction(phaseIdx,
224 compIdx,
225 x_ij);
226 }
227 // calculate the phase composition from the component
228 // fugacities
229 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
230 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
231 }
232
233 // dynamic viscosities
234 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
235 // viscosities
236 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
237 fluidState.setViscosity(phaseIdx, mu);
238
239 // compute and set the enthalpy
240 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
241 fluidState.setEnthalpy(phaseIdx, h);
242 }
243 }
244
249 const FluidState &fluidState() const
250 { return fluidState_; }
251
255 const SolidState &solidState() const
256 { return solidState_; }
257
264 Scalar saturation(int phaseIdx) const
265 { return fluidState_.saturation(phaseIdx); }
266
274 Scalar massFraction(const int phaseIdx, const int compIdx) const
275 { return fluidState_.massFraction(phaseIdx, compIdx); }
276
284 Scalar moleFraction(const int phaseIdx, const int compIdx) const
285 { return fluidState_.moleFraction(phaseIdx, compIdx); }
286
293 Scalar molarity(const int phaseIdx, int compIdx) const
294 { return fluidState_.molarity(phaseIdx, compIdx); }
295
301 Scalar molarDensity(const int phaseIdx) const
302 { return fluidState_.molarDensity(phaseIdx);}
303
310 Scalar pressure(const int phaseIdx) const
311 { return fluidState_.pressure(phaseIdx); }
312
318 Scalar density(const int phaseIdx) const
319 { return fluidState_.density(phaseIdx); }
320
328 Scalar temperature() const
329 { return fluidState_.temperature(0/* phaseIdx*/); }
330
331 Scalar temperature(const int phaseIdx) const
332 { return fluidState_.temperature(phaseIdx); }
333
337 Scalar enthalpy(const int phaseIdx) const
338 { return fluidState_.enthalpy(phaseIdx); }
339
343 Scalar internalEnergy(const int phaseIdx) const
344 { return fluidState_.internalEnergy(phaseIdx); }
345
350 Scalar fluidThermalConductivity(const int phaseIdx) const
351 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
352
356 Scalar fugacity(const int compIdx) const
357 { return fluidState_.fugacity(compIdx); }
358
362 Scalar averageMolarMass(const int phaseIdx) const
363 { return fluidState_.averageMolarMass(phaseIdx); }
364
371 Scalar mobility(const unsigned int phaseIdx) const
372 {
373 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
374 }
375
380 Scalar viscosity(const unsigned int phaseIdx) const
381 { return fluidState_.viscosity(phaseIdx); }
382
389 Scalar relativePermeability(const unsigned int phaseIdx) const
390 { return relativePermeability_[phaseIdx]; }
391
395 Scalar porosity() const
396 { return solidState_.porosity(); }
397
401 const PermeabilityType& permeability() const
402 { return permeability_; }
403
410 bool isPhaseActive(const unsigned int phaseIdx) const
411 {
412 return
413 phasePresentIneq(fluidState(), phaseIdx) -
414 phaseNotPresentIneq(fluidState(), phaseIdx)
415 >= 0;
416 }
417
421 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
422 {
423 typename FluidSystem::ParameterCache paramCache;
424 paramCache.updatePhase(fluidState_, phaseIdx);
425 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
426 }
427
431 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
432 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
433
439 Scalar phaseNcp(const unsigned int phaseIdx) const
440 {
441 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
442 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
443 if (aEval > bEval)
444 return phasePresentIneq(fluidState(), phaseIdx);
445 return phaseNotPresentIneq(fluidState(), phaseIdx);
446 };
447
455 Scalar phasePresentIneq(const FluidState &fluidState,
456 const unsigned int phaseIdx) const
457 { return fluidState.saturation(phaseIdx); }
458
466 Scalar phaseNotPresentIneq(const FluidState &fluidState,
467 const unsigned int phaseIdx) const
468 {
469 // difference of sum of mole fractions in the phase from 100%
470 Scalar a = 1;
471 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
472 a -= fluidState.moleFraction(phaseIdx, compIdx);
473 return a;
474 }
475
476protected:
477 Scalar porosity_;
478 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
479 PermeabilityType permeability_;
480
484
485 // Effective diffusion coefficients for the phases
486 DiffusionCoefficients effectiveDiffCoeff_;
487};
488
489template <class Traits>
491 : public PorousMediumFlowVolumeVariables<Traits>
492 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
493{
496 using Scalar = typename Traits::PrimaryVariables::value_type;
497 using PermeabilityType = typename Traits::PermeabilityType;
498
499 using ModelTraits = typename Traits::ModelTraits;
500 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
501
502 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
503 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
504 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
505
506 using Indices = typename ModelTraits::Indices;
507 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
509 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
510 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
511 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
512
513public:
515 using FluidSystem = typename Traits::FluidSystem;
517 using FluidState = typename Traits::FluidState;
519 using SolidState = typename Traits::SolidState;
521 using SolidSystem = typename Traits::SolidSystem;
522
524 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
526 static constexpr int numFluidComps = ParentType::numFluidComponents();
528
538 template<class ElemSol, class Problem, class Element, class Scv>
539 void update(const ElemSol& elemSol,
540 const Problem& problem,
541 const Element& element,
542 const Scv& scv)
543 {
544 ParentType::update(elemSol, problem, element, scv);
545
546 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
547
548 // calculate the remaining quantities
549 const auto& spatialParams = problem.spatialParams();
550 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
551 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
552 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
553 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
554
555 typename FluidSystem::ParameterCache paramCache;
556 paramCache.updateAll(fluidState_);
557
558 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
559
560 if constexpr (enableDiffusion)
561 {
562 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
563 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
564
565 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
566 }
567
568 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
569 permeability_ = spatialParams.permeability(element, scv, elemSol);
570 EnergyVolVars::updateEffectiveThermalConductivity();
571 }
572
584 template<class ElemSol, class Problem, class Element, class Scv>
585 void completeFluidState(const ElemSol& elemSol,
586 const Problem& problem,
587 const Element& element,
588 const Scv& scv,
589 FluidState& fluidState,
590 SolidState& solidState)
591 {
593 // set the fluid phase temperatures
595 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
597 // set the phase saturations
599 auto&& priVars = elemSol[scv.localDofIndex()];
600 Scalar sumSat = 0;
601 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
602 sumSat += priVars[Indices::s0Idx + phaseIdx];
603 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
604 }
605 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
606
608 // set the phase pressures
610 // capillary pressure parameters
611 const auto& spatialParams = problem.spatialParams();
612 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
613 fluidState.setWettingPhase(wPhaseIdx);
614 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
615 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
616
617 // add to the pressure of the first fluid phase
618 // depending on which pressure is stored in the primary variables
619 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
620 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
621 // For two phases this means that there is one pressure as primary variable: pw
622 const Scalar pw = priVars[Indices::p0Idx];
623 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
624 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
625 }
626 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
627 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
628 // For two phases this means that there is one pressure as primary variable: pn
629 const Scalar pn = priVars[Indices::p0Idx];
630 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
631 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
632 }
633 else
634 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
635
636
638 // set the fluid compositions
640 typename FluidSystem::ParameterCache paramCache;
641 paramCache.updateAll(fluidState);
642
643 ComponentVector fug;
644 // retrieve component fugacities
645 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
646 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
647
648 updateMoleFraction(fluidState,
649 paramCache,
650 priVars);
651
652
653 // dynamic viscosities
654 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
655 // viscosities
656 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
657 fluidState.setViscosity(phaseIdx, mu);
658
659 // compute and set the enthalpy
660 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
661 fluidState.setEnthalpy(phaseIdx, h);
662 }
663 }
664
673 void updateMoleFraction(FluidState & actualFluidState,
674 ParameterCache & paramCache,
675 const typename Traits::PrimaryVariables& priVars)
676 {
677 // setting the mole fractions of the fluid state
678 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
679 {
680 // set the component mole fractions
681 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
682 actualFluidState.setMoleFraction(phaseIdx,
683 compIdx,
684 priVars[Indices::moleFrac00Idx +
685 phaseIdx*numFluidComps +
686 compIdx]);
687 }
688 }
689
690// // For using the ... other way of calculating equilibrium
691// THIS IS ONLY FOR silencing Valgrind but is not used in this model
692// TODO: Can this be removed???
693 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
694 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
695 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
696 paramCache,
697 phaseIdx,
698 compIdx);
699 actualFluidState.setFugacityCoefficient(phaseIdx,
700 compIdx,
701 phi);
702 }
703
704 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
705 equilFluidState.assign(actualFluidState) ;
706 ConstraintSolver::solve(equilFluidState,
707 paramCache) ;
708
709 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
710 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
711 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
712 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
713 }
714 }
715
716 // compute densities of all phases
717 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
718 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
719 actualFluidState.setDensity(phaseIdx, rho);
720 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
721 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
722 }
723
724 }
725
733 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
734 {
735 return xEquil_[phaseIdx][compIdx] ;
736 }
737
742 const FluidState &fluidState() const
743 { return fluidState_; }
744
748 const SolidState &solidState() const
749 { return solidState_; }
750
757 Scalar saturation(int phaseIdx) const
758 { return fluidState_.saturation(phaseIdx); }
759
767 Scalar massFraction(const int phaseIdx, const int compIdx) const
768 { return fluidState_.massFraction(phaseIdx, compIdx); }
769
777 Scalar moleFraction(const int phaseIdx, const int compIdx) const
778 { return fluidState_.moleFraction(phaseIdx, compIdx); }
779
786 Scalar molarity(const int phaseIdx, int compIdx) const
787 { return fluidState_.molarity(phaseIdx, compIdx); }
788
794 Scalar molarDensity(const int phaseIdx) const
795 { return fluidState_.molarDensity(phaseIdx);}
796
803 Scalar pressure(const int phaseIdx) const
804 { return fluidState_.pressure(phaseIdx); }
805
811 Scalar density(const int phaseIdx) const
812 { return fluidState_.density(phaseIdx); }
813
821 Scalar temperature() const
822 { return fluidState_.temperature(0/* phaseIdx*/); }
823
827 Scalar enthalpy(const int phaseIdx) const
828 { return fluidState_.enthalpy(phaseIdx); }
829
833 Scalar internalEnergy(const int phaseIdx) const
834 { return fluidState_.internalEnergy(phaseIdx); }
835
840 Scalar fluidThermalConductivity(const int phaseIdx) const
841 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
842
846 Scalar fugacity(const int compIdx) const
847 { return fluidState_.fugacity(compIdx); }
848
852 Scalar averageMolarMass(const int phaseIdx) const
853 { return fluidState_.averageMolarMass(phaseIdx); }
854
861 Scalar mobility(const unsigned int phaseIdx) const
862 {
863 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
864 }
865
870 Scalar viscosity(const unsigned int phaseIdx) const
871 { return fluidState_.viscosity(phaseIdx); }
872
879 Scalar relativePermeability(const unsigned int phaseIdx) const
880 { return relativePermeability_[phaseIdx]; }
881
885 Scalar porosity() const
886 { return solidState_.porosity(); }
887
891 const PermeabilityType& permeability() const
892 { return permeability_; }
893
900 bool isPhaseActive(const unsigned int phaseIdx) const
901 {
902 return
903 phasePresentIneq(fluidState(), phaseIdx) -
904 phaseNotPresentIneq(fluidState(), phaseIdx)
905 >= 0;
906 }
907
911 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
912 {
913 typename FluidSystem::ParameterCache paramCache;
914 paramCache.updatePhase(fluidState_, phaseIdx);
915 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
916 }
917
921 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
922 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
923
929 Scalar phaseNcp(const unsigned int phaseIdx) const
930 {
931 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
932 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
933 if (aEval > bEval)
934 return phasePresentIneq(fluidState(), phaseIdx);
935 return phaseNotPresentIneq(fluidState(), phaseIdx);
936 };
937
945 Scalar phasePresentIneq(const FluidState &fluidState,
946 const unsigned int phaseIdx) const
947 { return fluidState.saturation(phaseIdx); }
948
956 Scalar phaseNotPresentIneq(const FluidState &fluidState,
957 const unsigned int phaseIdx) const
958 {
959 // difference of sum of mole fractions in the phase from 100%
960 Scalar a = 1;
961 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
962 a -= fluidState.moleFraction(phaseIdx, compIdx);
963 return a;
964 }
965
966protected:
967 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
968 PermeabilityType permeability_;
969 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
970
974
975 // Effective diffusion coefficients for the phases
976 DiffusionCoefficients effectiveDiffCoeff_;
977};
978
979} // end namespace
980
981#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 arbitary equation of state.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.
Enumeration of the formulations accepted by the MpNc model.
void updateSolidVolumeFractions(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, SolidState &solidState, const int solidVolFracOffset)
update the solid volume fractions (inert and reacitve) and set them in the solidstate
Definition: updatesolidvolumefractions.hh:36
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
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:69
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:97
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:59
Definition: porousmediumflow/mpnc/volumevariables.hh:42
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:421
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:83
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:328
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:362
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:482
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:264
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:310
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:249
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:350
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:395
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:371
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:255
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:389
Scalar massFraction(const int phaseIdx, const int compIdx) const
Returns the mass fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:274
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:102
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:293
typename Traits::ModelTraits::Indices Indices
Export the type encapsulating primary variable indices.
Definition: porousmediumflow/mpnc/volumevariables.hh:77
Scalar temperature(const int phaseIdx) const
Definition: porousmediumflow/mpnc/volumevariables.hh:331
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:401
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:356
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:337
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:431
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:483
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:318
Scalar moleFraction(const int phaseIdx, const int compIdx) const
Returns the mole fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:284
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:478
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:81
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:486
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:301
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:439
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:88
Scalar phasePresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is present.
Definition: porousmediumflow/mpnc/volumevariables.hh:455
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:380
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:150
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:479
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:343
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:85
Scalar phaseNotPresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is not present.
Definition: porousmediumflow/mpnc/volumevariables.hh:466
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:79
Scalar porosity_
Effective porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:477
bool isPhaseActive(const unsigned int phaseIdx) const
Returns true if the fluid state is in the active set for a phase,.
Definition: porousmediumflow/mpnc/volumevariables.hh:410
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:846
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:968
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:840
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:821
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:861
void updateMoleFraction(FluidState &actualFluidState, ParameterCache &paramCache, const typename Traits::PrimaryVariables &priVars)
Updates composition of all phases in the mutable parameters from the primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:673
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:803
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:748
Scalar phasePresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is present.
Definition: porousmediumflow/mpnc/volumevariables.hh:945
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:521
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:885
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:519
Scalar phaseNotPresentIneq(const FluidState &fluidState, const unsigned int phaseIdx) const
Returns the value of the inequality where a phase is not present.
Definition: porousmediumflow/mpnc/volumevariables.hh:956
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:517
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:827
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:811
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:921
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:973
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:515
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:911
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:891
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:794
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:870
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:833
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:969
Scalar moleFraction(const int phaseIdx, const int compIdx) const
Returns the mole fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:777
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:539
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:524
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:757
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:972
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:929
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:852
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:967
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:585
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:976
bool isPhaseActive(const unsigned int phaseIdx) const
Returns true if the fluid state is in the active set for a phase,.
Definition: porousmediumflow/mpnc/volumevariables.hh:900
Scalar massFraction(const int phaseIdx, const int compIdx) const
Returns the mass fraction of a given component in a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:767
const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
The mole fraction we would have in the case of chemical equilibrium / on the interface.
Definition: porousmediumflow/mpnc/volumevariables.hh:733
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:879
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:742
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:786
Definition: porousmediumflow/nonisothermal/volumevariables.hh:75
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:40
Base class for the model specific class which provides access to all volume averaged quantities.
Base class for the model specific class which provides access to all volume averaged quantities.