3.1-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 Indices = typename ModelTraits::Indices;
71 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
73
74public:
76 using FluidSystem = typename Traits::FluidSystem;
78 using FluidState = typename Traits::FluidState;
80 using SolidState = typename Traits::SolidState;
82 using SolidSystem = typename Traits::SolidSystem;
83
85 static constexpr int numFluidPhases() { return ModelTraits::numFluidPhases(); }
87 static constexpr int numFluidComps = ParentType::numFluidComponents();
88
98 template<class ElemSol, class Problem, class Element, class Scv>
99 void update(const ElemSol& elemSol,
100 const Problem& problem,
101 const Element& element,
102 const Scv& scv)
103 {
104 ParentType::update(elemSol, problem, element, scv);
105
106 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
107
108 //calculate the remaining quantities
109 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
110 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
111 // relative permeabilities
112 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
113 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
114 MPAdapter::relativePermeabilities(relativePermeability_, materialParams, fluidState_, wPhaseIdx);
115
116 typename FluidSystem::ParameterCache paramCache;
117 paramCache.updateAll(fluidState_);
118
119 if (enableDiffusion)
120 {
121 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
122 {
123 int compIIdx = phaseIdx;
124 for (unsigned int compJIdx = 0; compJIdx < numFluidComps; ++compJIdx)
125 {
126 // binary diffusion coefficients
127 if(compIIdx!= compJIdx)
128 {
129 setDiffusionCoefficient_(phaseIdx, compJIdx,
130 FluidSystem::binaryDiffusionCoefficient(fluidState_,
131 paramCache,
132 phaseIdx,
133 compIIdx,
134 compJIdx));
135 }
136 }
137 }
138 }
139 // porosity
140 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
141 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
142 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
143 }
144
157 template<class ElemSol, class Problem, class Element, class Scv>
158 void completeFluidState(const ElemSol& elemSol,
159 const Problem& problem,
160 const Element& element,
161 const Scv& scv,
162 FluidState& fluidState,
163 SolidState& solidState)
164 {
165
167 // set the fluid phase temperatures
169 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
170
172 // set the phase saturations
174 auto&& priVars = elemSol[scv.localDofIndex()];
175 Scalar sumSat = 0;
176 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
177 sumSat += priVars[Indices::s0Idx + phaseIdx];
178 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
179 }
181 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
182
184 // set the phase pressures
186 // capillary pressure parameters
187 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
188 fluidState.setWettingPhase(wPhaseIdx);
189 const auto& materialParams =
190 problem.spatialParams().materialLawParams(element, scv, elemSol);
191 // capillary pressures
192 std::vector<Scalar> capPress(numFluidPhases());
193 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
194 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
195 MPAdapter::capillaryPressures(capPress, materialParams, fluidState, wPhaseIdx);
196 // add to the pressure of the first fluid phase
197
198 // depending on which pressure is stored in the primary variables
199 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
200 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
201 // For two phases this means that there is one pressure as primary variable: pw
202 const Scalar pw = priVars[Indices::p0Idx];
203 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
204 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
205 }
206 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
207 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
208 // For two phases this means that there is one pressure as primary variable: pn
209 const Scalar pn = priVars[Indices::p0Idx];
210 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
211 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
212 }
213 else
214 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
215
217 // set the fluid compositions
219 typename FluidSystem::ParameterCache paramCache;
220 paramCache.updateAll(fluidState);
221
222 ComponentVector fug;
223 // retrieve component fugacities
224 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
225 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
226
227 // calculate phase compositions
228 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
229 // initial guess
230 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
231 Scalar x_ij = 1.0/numFluidComps;
232
233 // set initial guess of the component's mole fraction
234 fluidState.setMoleFraction(phaseIdx,
235 compIdx,
236 x_ij);
237 }
238 // calculate the phase composition from the component
239 // fugacities
240 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fug);
241 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fug);
242 }
243
244 // dynamic viscosities
245 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
246 // viscosities
247 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
248 fluidState.setViscosity(phaseIdx, mu);
249
250 // compute and set the enthalpy
251 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
252 fluidState.setEnthalpy(phaseIdx, h);
253 }
254 }
255
260 const FluidState &fluidState() const
261 { return fluidState_; }
262
266 const SolidState &solidState() const
267 { return solidState_; }
268
275 Scalar saturation(int phaseIdx) const
276 { return fluidState_.saturation(phaseIdx); }
277
285 Scalar massFraction(const int phaseIdx, const int compIdx) const
286 { return fluidState_.massFraction(phaseIdx, compIdx); }
287
295 Scalar moleFraction(const int phaseIdx, const int compIdx) const
296 { return fluidState_.moleFraction(phaseIdx, compIdx); }
297
304 Scalar molarity(const int phaseIdx, int compIdx) const
305 { return fluidState_.molarity(phaseIdx, compIdx); }
306
312 Scalar molarDensity(const int phaseIdx) const
313 { return fluidState_.molarDensity(phaseIdx);}
314
321 Scalar pressure(const int phaseIdx) const
322 { return fluidState_.pressure(phaseIdx); }
323
329 Scalar density(const int phaseIdx) const
330 { return fluidState_.density(phaseIdx); }
331
339 Scalar temperature() const
340 { return fluidState_.temperature(0/* phaseIdx*/); }
341
342 Scalar temperature(const int phaseIdx) const
343 { return fluidState_.temperature(phaseIdx); }
344
348 Scalar enthalpy(const int phaseIdx) const
349 { return fluidState_.enthalpy(phaseIdx); }
350
354 Scalar internalEnergy(const int phaseIdx) const
355 { return fluidState_.internalEnergy(phaseIdx); }
356
361 Scalar fluidThermalConductivity(const int phaseIdx) const
362 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
363
367 Scalar fugacity(const int compIdx) const
368 { return fluidState_.fugacity(compIdx); }
369
373 Scalar averageMolarMass(const int phaseIdx) const
374 { return fluidState_.averageMolarMass(phaseIdx); }
375
382 Scalar mobility(const unsigned int phaseIdx) const
383 {
384 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
385 }
386
391 Scalar viscosity(const unsigned int phaseIdx) const
392 { return fluidState_.viscosity(phaseIdx); }
393
400 Scalar relativePermeability(const unsigned int phaseIdx) const
401 { return relativePermeability_[phaseIdx]; }
402
406 Scalar porosity() const
407 { return solidState_.porosity(); }
408
412 const PermeabilityType& permeability() const
413 { return permeability_; }
414
421 bool isPhaseActive(const unsigned int phaseIdx) const
422 {
423 return
424 phasePresentIneq(fluidState(), phaseIdx) -
425 phaseNotPresentIneq(fluidState(), phaseIdx)
426 >= 0;
427 }
428
432 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
433 {
434 if (compIdx < phaseIdx)
435 return diffCoefficient_[phaseIdx][compIdx];
436 else if (compIdx > phaseIdx)
437 return diffCoefficient_[phaseIdx][compIdx-1];
438 else
439 DUNE_THROW(Dune::InvalidStateException, "Diffusion coeffiecient called for phaseIdx = compIdx");
440 }
441
447 Scalar phaseNcp(const unsigned int phaseIdx) const
448 {
449 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
450 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
451 if (aEval > bEval)
452 return phasePresentIneq(fluidState(), phaseIdx);
453 return phaseNotPresentIneq(fluidState(), phaseIdx);
454 };
455
463 Scalar phasePresentIneq(const FluidState &fluidState,
464 const unsigned int phaseIdx) const
465 { return fluidState.saturation(phaseIdx); }
466
474 Scalar phaseNotPresentIneq(const FluidState &fluidState,
475 const unsigned int phaseIdx) const
476 {
477 // difference of sum of mole fractions in the phase from 100%
478 Scalar a = 1;
479 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
480 a -= fluidState.moleFraction(phaseIdx, compIdx);
481 return a;
482 }
483
484protected:
485
486 void setDiffusionCoefficient_(int phaseIdx, int compIdx, Scalar d)
487 {
488 if (compIdx < phaseIdx)
489 diffCoefficient_[phaseIdx][compIdx] = std::move(d);
490 else if (compIdx > phaseIdx)
491 diffCoefficient_[phaseIdx][compIdx-1] = std::move(d);
492 else
493 DUNE_THROW(Dune::InvalidStateException, "Diffusion coeffiecient for phaseIdx = compIdx doesn't exist");
494 }
495
496 std::array<std::array<Scalar, numFluidComps-1>, numFluidPhases()> diffCoefficient_;
497 Scalar porosity_;
498 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
499 PermeabilityType permeability_;
500
504
505};
506
507template <class Traits>
509 : public PorousMediumFlowVolumeVariables<Traits>
510 , public EnergyVolumeVariables<Traits, MPNCVolumeVariables<Traits> >
511{
514 using Scalar = typename Traits::PrimaryVariables::value_type;
515 using PermeabilityType = typename Traits::PermeabilityType;
516
517 using ModelTraits = typename Traits::ModelTraits;
518 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
519
520 static constexpr bool enableThermalNonEquilibrium = ModelTraits::enableThermalNonEquilibrium();
521 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
522 static constexpr bool enableDiffusion = ModelTraits::enableMolecularDiffusion();
523
524 using Indices = typename ModelTraits::Indices;
525 using ComponentVector = Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()>;
527
528 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
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 if (enableDiffusion)
575 {
576 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
577 {
578 int compIIdx = phaseIdx;
579 for (unsigned int compJIdx = 0; compJIdx < numFluidComps; ++compJIdx)
580 {
581 // binary diffusion coefficients
582 if(compIIdx!= compJIdx)
583 {
584 setDiffusionCoefficient_(phaseIdx, compJIdx,
585 FluidSystem::binaryDiffusionCoefficient(fluidState_,
586 paramCache,
587 phaseIdx,
588 compIIdx,
589 compJIdx));
590 }
591 }
592 }
593 }
594
595 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
596 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
597 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
598
599 }
600
612 template<class ElemSol, class Problem, class Element, class Scv>
613 void completeFluidState(const ElemSol& elemSol,
614 const Problem& problem,
615 const Element& element,
616 const Scv& scv,
617 FluidState& fluidState,
618 SolidState& solidState)
619 {
621 // set the fluid phase temperatures
623 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
625 // set the phase saturations
627 auto&& priVars = elemSol[scv.localDofIndex()];
628 Scalar sumSat = 0;
629 for (int phaseIdx = 0; phaseIdx < numFluidPhases() - 1; ++phaseIdx) {
630 sumSat += priVars[Indices::s0Idx + phaseIdx];
631 fluidState.setSaturation(phaseIdx, priVars[Indices::s0Idx + phaseIdx]);
632 }
634 fluidState.setSaturation(numFluidPhases() - 1, 1.0 - sumSat);
635
637 // set the phase pressures
639 // capillary pressure parameters
640 const auto& materialParams =
641 problem.spatialParams().materialLawParams(element, scv, elemSol);
642 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
643 fluidState.setWettingPhase(wPhaseIdx);
644 // capillary pressures
645 std::vector<Scalar> capPress(numFluidPhases());
646 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
647 using MPAdapter = MPAdapter<MaterialLaw, numFluidPhases()>;
648 MPAdapter::capillaryPressures(capPress, materialParams, fluidState, wPhaseIdx);
649 // add to the pressure of the first fluid phase
650
651 // depending on which pressure is stored in the primary variables
652 if(pressureFormulation == MpNcPressureFormulation::mostWettingFirst){
653 // This means that the pressures are sorted from the most wetting to the least wetting-1 in the primary variables vector.
654 // For two phases this means that there is one pressure as primary variable: pw
655 const Scalar pw = priVars[Indices::p0Idx];
656 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx)
657 fluidState.setPressure(phaseIdx, pw - capPress[0] + capPress[phaseIdx]);
658 }
659 else if(pressureFormulation == MpNcPressureFormulation::leastWettingFirst){
660 // This means that the pressures are sorted from the least wetting to the most wetting-1 in the primary variables vector.
661 // For two phases this means that there is one pressure as primary variable: pn
662 const Scalar pn = priVars[Indices::p0Idx];
663 for (int phaseIdx = numFluidPhases()-1; phaseIdx >= 0; --phaseIdx)
664 fluidState.setPressure(phaseIdx, pn - capPress[numFluidPhases()-1] + capPress[phaseIdx]);
665 }
666 else
667 DUNE_THROW(Dune::InvalidStateException, "MPNCVolumeVariables do not support the chosen pressure formulation");
668
669
671 // set the fluid compositions
673 typename FluidSystem::ParameterCache paramCache;
674 paramCache.updateAll(fluidState);
675
676 ComponentVector fug;
677 // retrieve component fugacities
678 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
679 fug[compIdx] = priVars[Indices::fug0Idx + compIdx];
680
681 updateMoleFraction(fluidState,
682 paramCache,
683 priVars);
684
685
686 // dynamic viscosities
687 for (int phaseIdx = 0; phaseIdx < numFluidPhases(); ++phaseIdx) {
688 // viscosities
689 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
690 fluidState.setViscosity(phaseIdx, mu);
691
692 // compute and set the enthalpy
693 Scalar h = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
694 fluidState.setEnthalpy(phaseIdx, h);
695 }
696 }
697
706 void updateMoleFraction(FluidState & actualFluidState,
707 ParameterCache & paramCache,
708 const typename Traits::PrimaryVariables& priVars)
709 {
710 // setting the mole fractions of the fluid state
711 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
712 {
713 // set the component mole fractions
714 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
715 actualFluidState.setMoleFraction(phaseIdx,
716 compIdx,
717 priVars[Indices::moleFrac00Idx +
718 phaseIdx*numFluidComps +
719 compIdx]);
720 }
721 }
722
723// // For using the ... other way of calculating equilibrium
724// THIS IS ONLY FOR silencing Valgrind but is not used in this model
725 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx)
726 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx) {
727 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
728 paramCache,
729 phaseIdx,
730 compIdx);
731 actualFluidState.setFugacityCoefficient(phaseIdx,
732 compIdx,
733 phi);
734 }
735
736 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
737 equilFluidState.assign(actualFluidState) ;
738 ConstraintSolver::solve(equilFluidState,
739 paramCache) ;
740
741 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
742 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
743 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
744 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
745 }
746 }
747
748 // compute densities of all phases
749 for(int phaseIdx=0; phaseIdx<numFluidPhases(); ++phaseIdx){
750 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
751 actualFluidState.setDensity(phaseIdx, rho);
752 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
753 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
754 }
755
756 }
757
765 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
766 {
767 return xEquil_[phaseIdx][compIdx] ;
768 }
769
774 const FluidState &fluidState() const
775 { return fluidState_; }
776
780 const SolidState &solidState() const
781 { return solidState_; }
782
789 Scalar saturation(int phaseIdx) const
790 { return fluidState_.saturation(phaseIdx); }
791
799 Scalar massFraction(const int phaseIdx, const int compIdx) const
800 { return fluidState_.massFraction(phaseIdx, compIdx); }
801
809 Scalar moleFraction(const int phaseIdx, const int compIdx) const
810 { return fluidState_.moleFraction(phaseIdx, compIdx); }
811
818 Scalar molarity(const int phaseIdx, int compIdx) const
819 { return fluidState_.molarity(phaseIdx, compIdx); }
820
826 Scalar molarDensity(const int phaseIdx) const
827 { return fluidState_.molarDensity(phaseIdx);}
828
835 Scalar pressure(const int phaseIdx) const
836 { return fluidState_.pressure(phaseIdx); }
837
843 Scalar density(const int phaseIdx) const
844 { return fluidState_.density(phaseIdx); }
845
853 Scalar temperature() const
854 { return fluidState_.temperature(0/* phaseIdx*/); }
855
859 Scalar enthalpy(const int phaseIdx) const
860 { return fluidState_.enthalpy(phaseIdx); }
861
865 Scalar internalEnergy(const int phaseIdx) const
866 { return fluidState_.internalEnergy(phaseIdx); }
867
872 Scalar fluidThermalConductivity(const int phaseIdx) const
873 { return FluidSystem::thermalConductivity(fluidState_, phaseIdx); }
874
878 Scalar fugacity(const int compIdx) const
879 { return fluidState_.fugacity(compIdx); }
880
884 Scalar averageMolarMass(const int phaseIdx) const
885 { return fluidState_.averageMolarMass(phaseIdx); }
886
893 Scalar mobility(const unsigned int phaseIdx) const
894 {
895 return relativePermeability(phaseIdx)/fluidState_.viscosity(phaseIdx);
896 }
897
902 Scalar viscosity(const unsigned int phaseIdx) const
903 { return fluidState_.viscosity(phaseIdx); }
904
911 Scalar relativePermeability(const unsigned int phaseIdx) const
912 { return relativePermeability_[phaseIdx]; }
913
917 Scalar porosity() const
918 { return solidState_.porosity(); }
919
923 const PermeabilityType& permeability() const
924 { return permeability_; }
925
932 bool isPhaseActive(const unsigned int phaseIdx) const
933 {
934 return
935 phasePresentIneq(fluidState(), phaseIdx) -
936 phaseNotPresentIneq(fluidState(), phaseIdx)
937 >= 0;
938 }
939
943 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
944 {
945 if (compIdx < phaseIdx)
946 return diffCoefficient_[phaseIdx][compIdx];
947 else if (compIdx > phaseIdx)
948 return diffCoefficient_[phaseIdx][compIdx-1];
949 else
950 DUNE_THROW(Dune::InvalidStateException, "Diffusion coeffiecient called for phaseIdx = compIdx");
951 }
952
958 Scalar phaseNcp(const unsigned int phaseIdx) const
959 {
960 Scalar aEval = phaseNotPresentIneq(fluidState(), phaseIdx);
961 Scalar bEval = phasePresentIneq(fluidState(), phaseIdx);
962 if (aEval > bEval)
963 return phasePresentIneq(fluidState(), phaseIdx);
964 return phaseNotPresentIneq(fluidState(), phaseIdx);
965 };
966
974 Scalar phasePresentIneq(const FluidState &fluidState,
975 const unsigned int phaseIdx) const
976 { return fluidState.saturation(phaseIdx); }
977
985 Scalar phaseNotPresentIneq(const FluidState &fluidState,
986 const unsigned int phaseIdx) const
987 {
988 // difference of sum of mole fractions in the phase from 100%
989 Scalar a = 1;
990 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
991 a -= fluidState.moleFraction(phaseIdx, compIdx);
992 return a;
993 }
994
995protected:
996
997 void setDiffusionCoefficient_(int phaseIdx, int compIdx, Scalar d)
998 {
999 if (compIdx < phaseIdx)
1000 diffCoefficient_[phaseIdx][compIdx] = std::move(d);
1001 else if (compIdx > phaseIdx)
1002 diffCoefficient_[phaseIdx][compIdx-1] = std::move(d);
1003 else
1004 DUNE_THROW(Dune::InvalidStateException, "Diffusion coeffiecient for phaseIdx = compIdx doesn't exist");
1005 }
1006
1007 std::array<std::array<Scalar, numFluidComps-1>, numFluidPhases()> diffCoefficient_;
1008 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
1009 PermeabilityType permeability_;
1010 std::array<std::array<Scalar, numFluidComps>, numFluidPhases()> xEquil_;
1011
1015};
1016
1017} // end namespace
1018
1019#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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
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
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/mpnc/volumevariables.hh:80
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:339
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:373
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:502
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:275
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:321
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:260
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:361
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:406
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:382
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:266
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:400
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:285
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:99
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:304
Scalar temperature(const int phaseIdx) const
Definition: porousmediumflow/mpnc/volumevariables.hh:342
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:412
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:367
Scalar enthalpy(const int phaseIdx) const
Return enthalpy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:348
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:503
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:329
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:295
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:498
typename Traits::FluidState FluidState
Export the fluid state type.
Definition: porousmediumflow/mpnc/volumevariables.hh:78
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:312
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:447
static constexpr int numFluidPhases()
Return number of phases considered by the model.
Definition: porousmediumflow/mpnc/volumevariables.hh:85
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:463
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:391
void setDiffusionCoefficient_(int phaseIdx, int compIdx, Scalar d)
Definition: porousmediumflow/mpnc/volumevariables.hh:486
std::array< std::array< Scalar, numFluidComps-1 >, numFluidPhases()> diffCoefficient_
Definition: porousmediumflow/mpnc/volumevariables.hh:496
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:158
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:499
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:354
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:82
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:474
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:76
Scalar porosity_
Effective porosity within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:497
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the diffusion coefficient.
Definition: porousmediumflow/mpnc/volumevariables.hh:432
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:421
Scalar fugacity(const int compIdx) const
Returns the fugacity the of the component.
Definition: porousmediumflow/mpnc/volumevariables.hh:878
void setDiffusionCoefficient_(int phaseIdx, int compIdx, Scalar d)
Definition: porousmediumflow/mpnc/volumevariables.hh:997
PermeabilityType permeability_
Definition: porousmediumflow/mpnc/volumevariables.hh:1009
Scalar fluidThermalConductivity(const int phaseIdx) const
Returns the thermal conductivity of a fluid phase in the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:872
Scalar temperature() const
Returns the temperature inside the sub-control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:853
Scalar mobility(const unsigned int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:893
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:706
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:835
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:780
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:974
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:917
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:985
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:859
Scalar density(const int phaseIdx) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:843
SolidState solidState_
Definition: porousmediumflow/mpnc/volumevariables.hh:1014
typename Traits::FluidSystem FluidSystem
Export the underlying fluid system.
Definition: porousmediumflow/mpnc/volumevariables.hh:531
const PermeabilityType & permeability() const
Returns the permeability within the control volume in .
Definition: porousmediumflow/mpnc/volumevariables.hh:923
Scalar molarDensity(const int phaseIdx) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:826
Scalar viscosity(const unsigned int phaseIdx) const
Returns the viscosity of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:902
Scalar internalEnergy(const int phaseIdx) const
Returns the internal energy the of the fluid phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:865
std::array< std::array< Scalar, numFluidComps >, numFluidPhases()> xEquil_
Definition: porousmediumflow/mpnc/volumevariables.hh:1010
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:809
std::array< std::array< Scalar, numFluidComps-1 >, numFluidPhases()> diffCoefficient_
Definition: porousmediumflow/mpnc/volumevariables.hh:1007
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:789
FluidState fluidState_
Mass fractions of each component within each phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:1013
Scalar phaseNcp(const unsigned int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:958
Scalar averageMolarMass(const int phaseIdx) const
Returns the average molar mass the of the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:884
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the diffusion coefficient.
Definition: porousmediumflow/mpnc/volumevariables.hh:943
std::array< Scalar, ModelTraits::numFluidPhases()> relativePermeability_
Effective relative permeability within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:1008
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:613
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:932
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:799
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:765
Scalar relativePermeability(const unsigned int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: porousmediumflow/mpnc/volumevariables.hh:911
const FluidState & fluidState() const
Returns the fluid configuration at the given primary variables.
Definition: porousmediumflow/mpnc/volumevariables.hh:774
Scalar molarity(const int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/mpnc/volumevariables.hh:818
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.