3.3.0
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
30
33
38
39namespace Dumux {
40
41// forward declaration
42template <class Traits, bool enableChemicalNonEquilibrium>
50template <class Traits>
51using MPNCVolumeVariables = MPNCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium()>;
52
53template <class Traits>
56, public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
57{
60
61 using Scalar = typename Traits::PrimaryVariables::value_type;
62 using PermeabilityType = typename Traits::PermeabilityType;
63
64 using ModelTraits = typename Traits::ModelTraits;
65 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
66
67 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
68 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
69 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
70
71 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
73 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
74 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
75
76 // remove this after deprecation period, i.e. after release 3.3
77 static constexpr auto numPhases = std::integral_constant<int, ModelTraits::numFluidPhases()>{};
78
79public:
81 using Indices = typename Traits::ModelTraits::Indices;
83 using FluidSystem = typename Traits::FluidSystem;
85 using FluidState = typename Traits::FluidState;
87 using SolidState = typename Traits::SolidState;
89 using SolidSystem = typename Traits::SolidSystem;
90
92 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
94 static constexpr int numFluidComps = ParentType::numFluidComponents();
95
105 template<class ElemSol, class Problem, class Element, class Scv>
106 void update(const ElemSol& elemSol,
107 const Problem& problem,
108 const Element& element,
109 const Scv& scv)
110 {
111 ParentType::update(elemSol, problem, element, scv);
112
113 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
114
115 // Old material law interface is deprecated: Replace this by
116 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
117 // after the release of 3.3, when the deprecated interface is no longer supported.
118 // We can safely use the two-p wrapper here without breaking compatibility because the MPAdapter only supports
119 // two phases anyway...
120 const auto fluidMatrixInteraction = Deprecated::makeMPPcKrSw(Scalar{}, problem.spatialParams(), element, scv, elemSol, numPhases);
121
122 //calculate the remaining quantities
123 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
124 // relative permeabilities
125 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
126 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
127
128 typename FluidSystem::ParameterCache paramCache;
129 paramCache.updateAll(fluidState_);
130
131 // porosity
132 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
133
134 if constexpr (enableDiffusion)
135 {
136 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
137 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
138
139 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
140 }
141
142 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
143 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
144 EnergyVolVars::updateEffectiveThermalConductivity();
145 }
146
159 template<class ElemSol, class Problem, class Element, class Scv>
160 void completeFluidState(const ElemSol& elemSol,
161 const Problem& problem,
162 const Element& element,
163 const Scv& scv,
164 FluidState& fluidState,
165 SolidState& solidState)
166 {
167
169 // set the fluid phase temperatures
171 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
172
174 // set the phase saturations
176 auto&& priVars = elemSol[scv.localDofIndex()];
177 Scalar sumSat = 0;
178 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
179 sumSat += priVars[Indices::s0Idx + phaseIdx];
180 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
181 }
182 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
183
185 // set the phase pressures
187 // capillary pressure parameters
188 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
189 fluidState.setWettingPhase(wPhaseIdx);
190 // capillary pressures
191
192 // Old material law interface is deprecated: Replace this by
193 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
194 // after the release of 3.3, when the deprecated interface is no longer supported.
195 // We can safely use the two-p wrapper here without breaking compatibility because the MPAdapter only supports
196 // two phases anyway...
197 const auto fluidMatrixInteraction = Deprecated::makeMPPcKrSw(Scalar{}, problem.spatialParams(), element, scv, elemSol, numPhases);
198 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
199 // add to the pressure of the first fluid phase
200
201 // depending on which pressure is stored in the primary variables
202 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
203 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
204 // For two phases this means that there is one pressure as primary variable: pw
205 const Scalar pw = priVars[Indices::p0Idx];
206 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
207 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
208 }
209 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
210 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
211 // For two phases this means that there is one pressure as primary variable: pn
212 const Scalar pn = priVars[Indices::p0Idx];
213 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
214 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
215 }
216 else
217 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
218
220 // set the fluid compositions
222 typename FluidSystem::ParameterCache paramCache;
223 paramCache.updateAll(fluidState);
224
225 ComponentVector fug;
226 // retrieve component fugacities
227 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
228 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
229
230 // calculate phase compositions
231 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
232 // initial guess
233 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
234 Scalar x_ij = 1.0/numFluidComps;
235
236 // set initial guess of the component's mole fraction
237 fluidState.setMoleFraction(phaseIdx,
238 compIdx,
239 x_ij);
240 }
241 // calculate the phase composition from the component
242 // fugacities
243 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
244 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
245 }
246
247 // dynamic viscosities
248 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
249 // viscosities
250 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
251 fluidState.setViscosity(phaseIdx, mu);
252
253 // compute and set the enthalpy
254 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
255 fluidState.setEnthalpy(phaseIdx, h);
256 }
257 }
258
263 const FluidState &fluidState() const
264 { return fluidState_; }
265
269 const SolidState &solidState() const
270 { return solidState_; }
271
278 Scalar saturation(int phaseIdx) const
279 { return fluidState_.saturation(phaseIdx); }
280
288 Scalar massFraction(const int phaseIdx, const int compIdx) const
289 { return fluidState_.massFraction(phaseIdx, compIdx); }
290
298 Scalar moleFraction(const int phaseIdx, const int compIdx) const
299 { return fluidState_.moleFraction(phaseIdx, compIdx); }
300
307 Scalar molarity(const int phaseIdx, int compIdx) const
308 { return fluidState_.molarity(phaseIdx, compIdx); }
309
315 Scalar molarDensity(const int phaseIdx) const
316 { return fluidState_.molarDensity(phaseIdx);}
317
324 Scalar pressure(const int phaseIdx) const
325 { return fluidState_.pressure(phaseIdx); }
326
332 Scalar density(const int phaseIdx) const
333 { return fluidState_.density(phaseIdx); }
334
342 Scalar temperature() const
343 { return fluidState_.temperature(0/* phaseIdx*/); }
344
345 Scalar temperature(const int phaseIdx) const
346 { return fluidState_.temperature(phaseIdx); }
347
351 Scalar enthalpy(const int phaseIdx) const
352 { return fluidState_.enthalpy(phaseIdx); }
353
357 Scalar internalEnergy(const int phaseIdx) const
358 { return fluidState_.internalEnergy(phaseIdx); }
359
364 Scalar fluidThermalConductivity(const int phaseIdx) const
365 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
366
370 Scalar fugacity(const int compIdx) const
371 { return fluidState_.fugacity(compIdx); }
372
376 Scalar averageMolarMass(const int phaseIdx) const
377 { return fluidState_.averageMolarMass(phaseIdx); }
378
385 Scalar mobility(const unsigned int phaseIdx) const
386 {
387 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
388 }
389
394 Scalar viscosity(const unsigned int phaseIdx) const
395 { return fluidState_.viscosity(phaseIdx); }
396
403 Scalar relativePermeability(const unsigned int phaseIdx) const
404 { return relativePermeability_[phaseIdx]; }
405
409 Scalar porosity() const
410 { return solidState_.porosity(); }
411
415 const PermeabilityType& permeability() const
416 { return permeability_; }
417
424 bool isPhaseActive(const unsigned int phaseIdx) const
425 {
426 return
427 phasePresentIneq(fluidState(), phaseIdx) -
428 phaseNotPresentIneq(fluidState(), phaseIdx)
429 >= 0;
430 }
431
435 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
436 {
437 typename FluidSystem::ParameterCache paramCache;
438 paramCache.updatePhase(fluidState_, phaseIdx);
439 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
440 }
441
445 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
446 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
447
453 Scalar phaseNcp(const unsigned int phaseIdx) const
454 {
455 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
456 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
457 if (aEval > bEval)
458 return phasePresentIneq(fluidState(), phaseIdx);
459 return phaseNotPresentIneq(fluidState(), phaseIdx);
460 };
461
469 Scalar phasePresentIneq(const FluidState &fluidState,
470 const unsigned int phaseIdx) const
471 { return fluidState.saturation(phaseIdx); }
472
480 Scalar phaseNotPresentIneq(const FluidState &fluidState,
481 const unsigned int phaseIdx) const
482 {
483 // difference of sum of mole fractions in the phase from 100%
484 Scalar a = 1;
485 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
486 a -= fluidState.moleFraction(phaseIdx, compIdx);
487 return a;
488 }
489
490protected:
491 Scalar porosity_;
492 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
493 PermeabilityType permeability_;
494
498
499 // Effective diffusion coefficients for the phases
500 DiffusionCoefficients effectiveDiffCoeff_;
501};
502
503template <class Traits>
505 : public PorousMediumFlowVolumeVariables<Traits>
506 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
507{
510 using Scalar = typename Traits::PrimaryVariables::value_type;
511 using PermeabilityType = typename Traits::PermeabilityType;
512
513 using ModelTraits = typename Traits::ModelTraits;
514 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
515
516 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
517 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
518 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
519
520 using Indices = typename ModelTraits::Indices;
521 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
523 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
524 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
525 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
526
527 // remove this after deprecation period, i.e. after release 3.3
528 static constexpr auto numPhases = std::integral_constant<int, ModelTraits::numFluidPhases()>{};
529
530public:
532 using FluidSystem = typename Traits::FluidSystem;
534 using FluidState = typename Traits::FluidState;
536 using SolidState = typename Traits::SolidState;
538 using SolidSystem = typename Traits::SolidSystem;
539
541 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
543 static constexpr int numFluidComps = ParentType::numFluidComponents();
545
555 template<class ElemSol, class Problem, class Element, class Scv>
556 void update(const ElemSol& elemSol,
557 const Problem& problem,
558 const Element& element,
559 const Scv& scv)
560 {
561 ParentType::update(elemSol, problem, element, scv);
562
563 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
564
565 // calculate the remaining quantities
566 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
567
568 // Old material law interface is deprecated: Replace this by
569 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
570 // after the release of 3.3, when the deprecated interface is no longer supported.
571 const auto fluidMatrixInteraction = Deprecated::makeMPPcKrSw(Scalar{}, problem.spatialParams(), element, scv, elemSol, numPhases);
572
573 const auto relPerm = fluidMatrixInteraction.relativePermeabilities(fluidState_, wPhaseIdx);
574 std::copy(relPerm.begin(), relPerm.end(), relativePermeability_.begin());
575
576 typename FluidSystem::ParameterCache paramCache;
577 paramCache.updateAll(fluidState_);
578
579 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
580
581 if constexpr (enableDiffusion)
582 {
583 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
584 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
585
586 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
587 }
588
589 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
590 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
591 EnergyVolVars::updateEffectiveThermalConductivity();
592 }
593
605 template<class ElemSol, class Problem, class Element, class Scv>
606 void completeFluidState(const ElemSol& elemSol,
607 const Problem& problem,
608 const Element& element,
609 const Scv& scv,
610 FluidState& fluidState,
611 SolidState& solidState)
612 {
614 // set the fluid phase temperatures
616 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
618 // set the phase saturations
620 auto&& priVars = elemSol[scv.localDofIndex()];
621 Scalar sumSat = 0;
622 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
623 sumSat += priVars[Indices::s0Idx + phaseIdx];
624 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
625 }
626 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
627
629 // set the phase pressures
631 // capillary pressure parameters
632 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
633 fluidState.setWettingPhase(wPhaseIdx);
634 // capillary pressures
635
636 // Old material law interface is deprecated: Replace this by
637 // const auto& fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
638 // after the release of 3.3, when the deprecated interface is no longer supported.
639 // We can safely use the two-p wrapper here without breaking compatibility because the MPAdapter only supports
640 // two phases anyway...
641 const auto fluidMatrixInteraction = Deprecated::makeMPPcKrSw(Scalar{}, problem.spatialParams(), element, scv, elemSol, numPhases);
642 const auto capPress = fluidMatrixInteraction.capillaryPressures(fluidState, wPhaseIdx);
643 // add to the pressure of the first fluid phase
644
645 // depending on which pressure is stored in the primary variables
646 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
647 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
648 // For two phases this means that there is one pressure as primary variable: pw
649 const Scalar pw = priVars[Indices::p0Idx];
650 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
651 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
652 }
653 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
654 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
655 // For two phases this means that there is one pressure as primary variable: pn
656 const Scalar pn = priVars[Indices::p0Idx];
657 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
658 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
659 }
660 else
661 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
662
663
665 // set the fluid compositions
667 typename FluidSystem::ParameterCache paramCache;
668 paramCache.updateAll(fluidState);
669
670 ComponentVector fug;
671 // retrieve component fugacities
672 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
673 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
674
675 updateMoleFraction(fluidState,
676 paramCache,
677 priVars);
678
679
680 // dynamic viscosities
681 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
682 // viscosities
683 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
684 fluidState.setViscosity(phaseIdx, mu);
685
686 // compute and set the enthalpy
687 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
688 fluidState.setEnthalpy(phaseIdx, h);
689 }
690 }
691
700 void updateMoleFraction(FluidState & actualFluidState,
701 ParameterCache & paramCache,
702 const typename Traits::PrimaryVariables& priVars)
703 {
704 // setting the mole fractions of the fluid state
705 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
706 {
707 // set the component mole fractions
708 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
709 actualFluidState.setMoleFraction(phaseIdx,
710 compIdx,
711 priVars[Indices::moleFrac00Idx +
712 phaseIdx*numFluidComps +
713 compIdx]);
714 }
715 }
716
717// // For using the ... other way of calculating equilibrium
718// THIS IS ONLY FOR silencing Valgrind but is not used in this model
719// TODO: Can this be removed???
720 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
721 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
722 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
723 paramCache,
724 phaseIdx,
725 compIdx);
726 actualFluidState.setFugacityCoefficient(phaseIdx,
727 compIdx,
728 phi);
729 }
730
731 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
732 equilFluidState.assign(actualFluidState) ;
733 ConstraintSolver::solve(equilFluidState,
734 paramCache) ;
735
736 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
737 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
738 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
739 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
740 }
741 }
742
743 // compute densities of all phases
744 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
745 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
746 actualFluidState.setDensity(phaseIdx, rho);
747 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
748 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
749 }
750
751 }
752
760 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
761 {
762 return xEquil_[phaseIdx][compIdx] ;
763 }
764
769 const FluidState &fluidState() const
770 { return fluidState_; }
771
775 const SolidState &solidState() const
776 { return solidState_; }
777
784 Scalar saturation(int phaseIdx) const
785 { return fluidState_.saturation(phaseIdx); }
786
794 Scalar massFraction(const int phaseIdx, const int compIdx) const
795 { return fluidState_.massFraction(phaseIdx, compIdx); }
796
804 Scalar moleFraction(const int phaseIdx, const int compIdx) const
805 { return fluidState_.moleFraction(phaseIdx, compIdx); }
806
813 Scalar molarity(const int phaseIdx, int compIdx) const
814 { return fluidState_.molarity(phaseIdx, compIdx); }
815
821 Scalar molarDensity(const int phaseIdx) const
822 { return fluidState_.molarDensity(phaseIdx);}
823
830 Scalar pressure(const int phaseIdx) const
831 { return fluidState_.pressure(phaseIdx); }
832
838 Scalar density(const int phaseIdx) const
839 { return fluidState_.density(phaseIdx); }
840
848 Scalar temperature() const
849 { return fluidState_.temperature(0/* phaseIdx*/); }
850
854 Scalar enthalpy(const int phaseIdx) const
855 { return fluidState_.enthalpy(phaseIdx); }
856
860 Scalar internalEnergy(const int phaseIdx) const
861 { return fluidState_.internalEnergy(phaseIdx); }
862
867 Scalar fluidThermalConductivity(const int phaseIdx) const
868 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
869
873 Scalar fugacity(const int compIdx) const
874 { return fluidState_.fugacity(compIdx); }
875
879 Scalar averageMolarMass(const int phaseIdx) const
880 { return fluidState_.averageMolarMass(phaseIdx); }
881
888 Scalar mobility(const unsigned int phaseIdx) const
889 {
890 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
891 }
892
897 Scalar viscosity(const unsigned int phaseIdx) const
898 { return fluidState_.viscosity(phaseIdx); }
899
906 Scalar relativePermeability(const unsigned int phaseIdx) const
907 { return relativePermeability_[phaseIdx]; }
908
912 Scalar porosity() const
913 { return solidState_.porosity(); }
914
918 const PermeabilityType& permeability() const
919 { return permeability_; }
920
927 bool isPhaseActive(const unsigned int phaseIdx) const
928 {
929 return
930 phasePresentIneq(fluidState(), phaseIdx) -
931 phaseNotPresentIneq(fluidState(), phaseIdx)
932 >= 0;
933 }
934
938 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
939 {
940 typename FluidSystem::ParameterCache paramCache;
941 paramCache.updatePhase(fluidState_, phaseIdx);
942 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
943 }
944
948 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
949 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
950
956 Scalar phaseNcp(const unsigned int phaseIdx) const
957 {
958 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
959 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
960 if (aEval > bEval)
961 return phasePresentIneq(fluidState(), phaseIdx);
962 return phaseNotPresentIneq(fluidState(), phaseIdx);
963 };
964
972 Scalar phasePresentIneq(const FluidState &fluidState,
973 const unsigned int phaseIdx) const
974 { return fluidState.saturation(phaseIdx); }
975
983 Scalar phaseNotPresentIneq(const FluidState &fluidState,
984 const unsigned int phaseIdx) const
985 {
986 // difference of sum of mole fractions in the phase from 100%
987 Scalar a = 1;
988 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
989 a -= fluidState.moleFraction(phaseIdx, compIdx);
990 return a;
991 }
992
993protected:
994 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
995 PermeabilityType permeability_;
996 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
997
1001
1002 // Effective diffusion coefficients for the phases
1003 DiffusionCoefficients effectiveDiffCoeff_;
1004};
1005
1006} // end namespace
1007
1008#endif
Helpers for deprecation.
Determines the fluid composition given the component fugacities and an arbitary equation of state.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
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:43
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:435
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:87
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:342
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:376
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:496
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:278
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:324
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:263
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:364
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:409
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:385
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:269
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:403
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:288
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:106
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:307
typename Traits::ModelTraits::Indices Indices
Export the type encapsulating primary variable indices.
Definition: porousmediumflow/mpnc/volumevariables.hh:81
Scalar temperature(const int phaseIdx) const
Definition: porousmediumflow/mpnc/volumevariables.hh:345
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:415
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:370
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:351
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:445
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:497
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:332
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:298
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:492
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:85
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:500
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:315
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:453
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:92
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:469
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:394
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:160
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:493
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:357
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:89
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:480
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:83
Scalar porosity_
Effective porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:491
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:424
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:873
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:995
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:867
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:848
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:888
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:700
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:830
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:775
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:972
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:538
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:912
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:536
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:983
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:534
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:854
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:838
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:948
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:1000
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:532
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:938
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:918
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:821
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:897
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:860
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:996
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:804
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:556
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:541
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:784
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:999
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:956
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:879
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:994
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:606
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:1003
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:927
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:794
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:760
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:906
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:769
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:813
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.