3.2-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
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& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
113 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
114 // relative permeabilities
115 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
116 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
117 MPAdapter::relativePermeabilities(relativePermeability_, materialParams, fluidState_, wPhaseIdx);
118
119 typename FluidSystem::ParameterCache paramCache;
120 paramCache.updateAll(fluidState_);
121
122 // porosity
123 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
124
125 if constexpr (enableDiffusion)
126 {
127 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
128 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
129
130 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
131 }
132
133 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
134 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
135 EnergyVolVars::updateEffectiveThermalConductivity();
136 }
137
150 template<class ElemSol, class Problem, class Element, class Scv>
151 void completeFluidState(const ElemSol& elemSol,
152 const Problem& problem,
153 const Element& element,
154 const Scv& scv,
155 FluidState& fluidState,
156 SolidState& solidState)
157 {
158
160 // set the fluid phase temperatures
162 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
163
165 // set the phase saturations
167 auto&& priVars = elemSol[scv.localDofIndex()];
168 Scalar sumSat = 0;
169 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
170 sumSat += priVars[Indices::s0Idx + phaseIdx];
171 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
172 }
174 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
175
177 // set the phase pressures
179 // capillary pressure parameters
180 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
181 fluidState.setWettingPhase(wPhaseIdx);
182 const auto& materialParams =
183 problem.spatialParams().materialLawParams(element, scv, elemSol);
184 // capillary pressures
185 std::vector<Scalar> capPress(numFluidPhases());
186 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
187 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
188 MPAdapter::capillaryPressures(capPress, materialParams, fluidState, wPhaseIdx);
189 // add to the pressure of the first fluid phase
190
191 // depending on which pressure is stored in the primary variables
192 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
193 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
194 // For two phases this means that there is one pressure as primary variable: pw
195 const Scalar pw = priVars[Indices::p0Idx];
196 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
197 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
198 }
199 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
200 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
201 // For two phases this means that there is one pressure as primary variable: pn
202 const Scalar pn = priVars[Indices::p0Idx];
203 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
204 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
205 }
206 else
207 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
208
210 // set the fluid compositions
212 typename FluidSystem::ParameterCache paramCache;
213 paramCache.updateAll(fluidState);
214
215 ComponentVector fug;
216 // retrieve component fugacities
217 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
218 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
219
220 // calculate phase compositions
221 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
222 // initial guess
223 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
224 Scalar x_ij = 1.0/numFluidComps;
225
226 // set initial guess of the component's mole fraction
227 fluidState.setMoleFraction(phaseIdx,
228 compIdx,
229 x_ij);
230 }
231 // calculate the phase composition from the component
232 // fugacities
233 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
234 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
235 }
236
237 // dynamic viscosities
238 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
239 // viscosities
240 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
241 fluidState.setViscosity(phaseIdx, mu);
242
243 // compute and set the enthalpy
244 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
245 fluidState.setEnthalpy(phaseIdx, h);
246 }
247 }
248
253 const FluidState &fluidState() const
254 { return fluidState_; }
255
259 const SolidState &solidState() const
260 { return solidState_; }
261
268 Scalar saturation(int phaseIdx) const
269 { return fluidState_.saturation(phaseIdx); }
270
278 Scalar massFraction(const int phaseIdx, const int compIdx) const
279 { return fluidState_.massFraction(phaseIdx, compIdx); }
280
288 Scalar moleFraction(const int phaseIdx, const int compIdx) const
289 { return fluidState_.moleFraction(phaseIdx, compIdx); }
290
297 Scalar molarity(const int phaseIdx, int compIdx) const
298 { return fluidState_.molarity(phaseIdx, compIdx); }
299
305 Scalar molarDensity(const int phaseIdx) const
306 { return fluidState_.molarDensity(phaseIdx);}
307
314 Scalar pressure(const int phaseIdx) const
315 { return fluidState_.pressure(phaseIdx); }
316
322 Scalar density(const int phaseIdx) const
323 { return fluidState_.density(phaseIdx); }
324
332 Scalar temperature() const
333 { return fluidState_.temperature(0/* phaseIdx*/); }
334
335 Scalar temperature(const int phaseIdx) const
336 { return fluidState_.temperature(phaseIdx); }
337
341 Scalar enthalpy(const int phaseIdx) const
342 { return fluidState_.enthalpy(phaseIdx); }
343
347 Scalar internalEnergy(const int phaseIdx) const
348 { return fluidState_.internalEnergy(phaseIdx); }
349
354 Scalar fluidThermalConductivity(const int phaseIdx) const
355 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
356
360 Scalar fugacity(const int compIdx) const
361 { return fluidState_.fugacity(compIdx); }
362
366 Scalar averageMolarMass(const int phaseIdx) const
367 { return fluidState_.averageMolarMass(phaseIdx); }
368
375 Scalar mobility(const unsigned int phaseIdx) const
376 {
377 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
378 }
379
384 Scalar viscosity(const unsigned int phaseIdx) const
385 { return fluidState_.viscosity(phaseIdx); }
386
393 Scalar relativePermeability(const unsigned int phaseIdx) const
394 { return relativePermeability_[phaseIdx]; }
395
399 Scalar porosity() const
400 { return solidState_.porosity(); }
401
405 const PermeabilityType& permeability() const
406 { return permeability_; }
407
414 bool isPhaseActive(const unsigned int phaseIdx) const
415 {
416 return
417 phasePresentIneq(fluidState(), phaseIdx) -
418 phaseNotPresentIneq(fluidState(), phaseIdx)
419 >= 0;
420 }
421
425 [[deprecated("Signature deprecated. Use diffusionCoefficient(phaseIdx, compIIdx, compJIdx)!")]]
426 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
427 {
428 if (compIdx != phaseIdx)
429 return diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
430 else
431 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficient called for phaseIdx = compIdx");
432 }
433
437 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
438 {
439 typename FluidSystem::ParameterCache paramCache;
440 paramCache.updatePhase(fluidState_, phaseIdx);
441 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
442 }
443
447 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
448 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
449
455 Scalar phaseNcp(const unsigned int phaseIdx) const
456 {
457 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
458 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
459 if (aEval > bEval)
460 return phasePresentIneq(fluidState(), phaseIdx);
461 return phaseNotPresentIneq(fluidState(), phaseIdx);
462 };
463
471 Scalar phasePresentIneq(const FluidState &fluidState,
472 const unsigned int phaseIdx) const
473 { return fluidState.saturation(phaseIdx); }
474
482 Scalar phaseNotPresentIneq(const FluidState &fluidState,
483 const unsigned int phaseIdx) const
484 {
485 // difference of sum of mole fractions in the phase from 100%
486 Scalar a = 1;
487 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
488 a -= fluidState.moleFraction(phaseIdx, compIdx);
489 return a;
490 }
491
492protected:
493 Scalar porosity_;
494 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
495 PermeabilityType permeability_;
496
500
501 // Effective diffusion coefficients for the phases
502 DiffusionCoefficients effectiveDiffCoeff_;
503};
504
505template <class Traits>
507 : public PorousMediumFlowVolumeVariables<Traits>
508 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
509{
512 using Scalar = typename Traits::PrimaryVariables::value_type;
513 using PermeabilityType = typename Traits::PermeabilityType;
514
515 using ModelTraits = typename Traits::ModelTraits;
516 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
517
518 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
519 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
520 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
521
522 using Indices = typename ModelTraits::Indices;
523 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
525 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
526 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
527 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
528
529public:
531 using FluidSystem = typename Traits::FluidSystem;
533 using FluidState = typename Traits::FluidState;
535 using SolidState = typename Traits::SolidState;
537 using SolidSystem = typename Traits::SolidSystem;
538
540 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
542 static constexpr int numFluidComps = ParentType::numFluidComponents();
544
554 template<class ElemSol, class Problem, class Element, class Scv>
555 void update(const ElemSol& elemSol,
556 const Problem& problem,
557 const Element& element,
558 const Scv& scv)
559 {
560 ParentType::update(elemSol, problem, element, scv);
561
562 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
563
564 // calculate the remaining quantities
565 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
566 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
567
568 // relative permeabilities
569 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
570 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
571 MPAdapter::relativePermeabilities(relativePermeability_, materialParams, fluidState_, wPhaseIdx);
572 typename FluidSystem::ParameterCache paramCache;
573 paramCache.updateAll(fluidState_);
574
575 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
576
577 if constexpr (enableDiffusion)
578 {
579 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
580 { return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx); };
581
582 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
583 }
584
585 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
586 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
587 EnergyVolVars::updateEffectiveThermalConductivity();
588 }
589
601 template<class ElemSol, class Problem, class Element, class Scv>
602 void completeFluidState(const ElemSol& elemSol,
603 const Problem& problem,
604 const Element& element,
605 const Scv& scv,
606 FluidState& fluidState,
607 SolidState& solidState)
608 {
610 // set the fluid phase temperatures
612 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
614 // set the phase saturations
616 auto&& priVars = elemSol[scv.localDofIndex()];
617 Scalar sumSat = 0;
618 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
619 sumSat += priVars[Indices::s0Idx + phaseIdx];
620 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
621 }
623 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
624
626 // set the phase pressures
628 // capillary pressure parameters
629 const auto& materialParams =
630 problem.spatialParams().materialLawParams(element, scv, elemSol);
631 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
632 fluidState.setWettingPhase(wPhaseIdx);
633 // capillary pressures
634 std::vector<Scalar> capPress(numFluidPhases());
635 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
636 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
637 MPAdapter::capillaryPressures(capPress, materialParams, fluidState, wPhaseIdx);
638 // add to the pressure of the first fluid phase
639
640 // depending on which pressure is stored in the primary variables
641 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
642 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
643 // For two phases this means that there is one pressure as primary variable: pw
644 const Scalar pw = priVars[Indices::p0Idx];
645 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
646 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
647 }
648 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
649 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
650 // For two phases this means that there is one pressure as primary variable: pn
651 const Scalar pn = priVars[Indices::p0Idx];
652 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
653 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
654 }
655 else
656 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
657
658
660 // set the fluid compositions
662 typename FluidSystem::ParameterCache paramCache;
663 paramCache.updateAll(fluidState);
664
665 ComponentVector fug;
666 // retrieve component fugacities
667 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
668 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
669
670 updateMoleFraction(fluidState,
671 paramCache,
672 priVars);
673
674
675 // dynamic viscosities
676 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
677 // viscosities
678 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
679 fluidState.setViscosity(phaseIdx, mu);
680
681 // compute and set the enthalpy
682 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
683 fluidState.setEnthalpy(phaseIdx, h);
684 }
685 }
686
695 void updateMoleFraction(FluidState & actualFluidState,
696 ParameterCache & paramCache,
697 const typename Traits::PrimaryVariables& priVars)
698 {
699 // setting the mole fractions of the fluid state
700 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
701 {
702 // set the component mole fractions
703 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
704 actualFluidState.setMoleFraction(phaseIdx,
705 compIdx,
706 priVars[Indices::moleFrac00Idx +
707 phaseIdx*numFluidComps +
708 compIdx]);
709 }
710 }
711
712// // For using the ... other way of calculating equilibrium
713// THIS IS ONLY FOR silencing Valgrind but is not used in this model
714 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
715 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
716 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
717 paramCache,
718 phaseIdx,
719 compIdx);
720 actualFluidState.setFugacityCoefficient(phaseIdx,
721 compIdx,
722 phi);
723 }
724
725 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
726 equilFluidState.assign(actualFluidState) ;
727 ConstraintSolver::solve(equilFluidState,
728 paramCache) ;
729
730 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
731 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
732 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
733 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
734 }
735 }
736
737 // compute densities of all phases
738 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
739 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
740 actualFluidState.setDensity(phaseIdx, rho);
741 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
742 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
743 }
744
745 }
746
754 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
755 {
756 return xEquil_[phaseIdx][compIdx] ;
757 }
758
763 const FluidState &fluidState() const
764 { return fluidState_; }
765
769 const SolidState &solidState() const
770 { return solidState_; }
771
778 Scalar saturation(int phaseIdx) const
779 { return fluidState_.saturation(phaseIdx); }
780
788 Scalar massFraction(const int phaseIdx, const int compIdx) const
789 { return fluidState_.massFraction(phaseIdx, compIdx); }
790
798 Scalar moleFraction(const int phaseIdx, const int compIdx) const
799 { return fluidState_.moleFraction(phaseIdx, compIdx); }
800
807 Scalar molarity(const int phaseIdx, int compIdx) const
808 { return fluidState_.molarity(phaseIdx, compIdx); }
809
815 Scalar molarDensity(const int phaseIdx) const
816 { return fluidState_.molarDensity(phaseIdx);}
817
824 Scalar pressure(const int phaseIdx) const
825 { return fluidState_.pressure(phaseIdx); }
826
832 Scalar density(const int phaseIdx) const
833 { return fluidState_.density(phaseIdx); }
834
842 Scalar temperature() const
843 { return fluidState_.temperature(0/* phaseIdx*/); }
844
848 Scalar enthalpy(const int phaseIdx) const
849 { return fluidState_.enthalpy(phaseIdx); }
850
854 Scalar internalEnergy(const int phaseIdx) const
855 { return fluidState_.internalEnergy(phaseIdx); }
856
861 Scalar fluidThermalConductivity(const int phaseIdx) const
862 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
863
867 Scalar fugacity(const int compIdx) const
868 { return fluidState_.fugacity(compIdx); }
869
873 Scalar averageMolarMass(const int phaseIdx) const
874 { return fluidState_.averageMolarMass(phaseIdx); }
875
882 Scalar mobility(const unsigned int phaseIdx) const
883 {
884 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
885 }
886
891 Scalar viscosity(const unsigned int phaseIdx) const
892 { return fluidState_.viscosity(phaseIdx); }
893
900 Scalar relativePermeability(const unsigned int phaseIdx) const
901 { return relativePermeability_[phaseIdx]; }
902
906 Scalar porosity() const
907 { return solidState_.porosity(); }
908
912 const PermeabilityType& permeability() const
913 { return permeability_; }
914
921 bool isPhaseActive(const unsigned int phaseIdx) const
922 {
923 return
924 phasePresentIneq(fluidState(), phaseIdx) -
925 phaseNotPresentIneq(fluidState(), phaseIdx)
926 >= 0;
927 }
928
932 [[deprecated("Will be removed after release 3.2. Use diffusionCoefficient(phaseIdx, compIIdx, compJIdx)!")]]
933 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
934 {
935 if (compIdx != phaseIdx)
936 return diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
937 else
938 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficient called for phaseIdx = compIdx");
939 }
940
944 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
945 {
946 typename FluidSystem::ParameterCache paramCache;
947 paramCache.updatePhase(fluidState_, phaseIdx);
948 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
949 }
950
954 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
955 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
956
962 Scalar phaseNcp(const unsigned int phaseIdx) const
963 {
964 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
965 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
966 if (aEval > bEval)
967 return phasePresentIneq(fluidState(), phaseIdx);
968 return phaseNotPresentIneq(fluidState(), phaseIdx);
969 };
970
978 Scalar phasePresentIneq(const FluidState &fluidState,
979 const unsigned int phaseIdx) const
980 { return fluidState.saturation(phaseIdx); }
981
989 Scalar phaseNotPresentIneq(const FluidState &fluidState,
990 const unsigned int phaseIdx) const
991 {
992 // difference of sum of mole fractions in the phase from 100%
993 Scalar a = 1;
994 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
995 a -= fluidState.moleFraction(phaseIdx, compIdx);
996 return a;
997 }
998
999protected:
1000 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
1001 PermeabilityType permeability_;
1002 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
1003
1007
1008 // Effective diffusion coefficients for the phases
1009 DiffusionCoefficients effectiveDiffCoeff_;
1010};
1011
1012} // end namespace
1013
1014#endif
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...
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.
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
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
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:55
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:71
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:99
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:60
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:42
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:437
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:332
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:366
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:498
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:268
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:314
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:253
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:354
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:399
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:375
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:259
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:393
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:278
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:297
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:335
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:405
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:360
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:341
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:447
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:499
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:322
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:288
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:494
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:81
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:502
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:305
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:455
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:471
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:384
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:151
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:495
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:347
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:482
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:493
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the diffusion coefficient.
Definition: porousmediumflow/mpnc/volumevariables.hh:426
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:414
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:867
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:1001
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:861
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:842
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:882
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:695
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:824
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:769
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:978
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:537
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:906
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:535
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:989
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:533
Scalar enthalpy(const int phaseIdx) const
Returns the enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:848
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:832
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:954
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:1006
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:531
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/mpnc/volumevariables.hh:944
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:912
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:815
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:891
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:854
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:1002
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:798
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:555
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:540
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:778
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:1005
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:962
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:873
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the diffusion coefficient.
Definition: porousmediumflow/mpnc/volumevariables.hh:933
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:1000
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:602
DiffusionCoefficients effectiveDiffCoeff_
Definition: porousmediumflow/mpnc/volumevariables.hh:1009
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:921
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:788
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:754
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:900
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:763
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:807
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.