3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/2p2c/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_2P2C_VOLUME_VARIABLES_HH
27#define DUMUX_2P2C_VOLUME_VARIABLES_HH
28
32
39
40namespace Dumux {
41
42// forward declaration
43template <class Traits, bool enableChemicalNonEquilibrium, bool useConstraintSolver>
45
46
52template <class Traits, bool useConstraintSolver = true>
53using TwoPTwoCVolumeVariables = TwoPTwoCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium(), useConstraintSolver>;
54
61template <class Traits, class Impl>
64, public EnergyVolumeVariables<Traits, TwoPTwoCVolumeVariables<Traits> >
65{
68
69 using Scalar = typename Traits::PrimaryVariables::value_type;
70 using ModelTraits = typename Traits::ModelTraits;
71
72 static constexpr int numFluidComps = ParentType::numFluidComponents();
73 // component indices
74 enum
75 {
76 comp0Idx = Traits::FluidSystem::comp0Idx,
77 comp1Idx = Traits::FluidSystem::comp1Idx,
78 phase0Idx = Traits::FluidSystem::phase0Idx,
79 phase1Idx = Traits::FluidSystem::phase1Idx
80 };
81
82 // phase presence indices
83 enum
84 {
85 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
86 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
87 bothPhases = ModelTraits::Indices::bothPhases
88 };
89
90 // primary variable indices
91 enum
92 {
93 switchIdx = ModelTraits::Indices::switchIdx,
94 pressureIdx = ModelTraits::Indices::pressureIdx
95 };
96
97 // formulations
98 static constexpr auto formulation = ModelTraits::priVarFormulation();
99
100 using PermeabilityType = typename Traits::PermeabilityType;
103 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
104 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
105public:
107 using FluidState = typename Traits::FluidState;
109 using FluidSystem = typename Traits::FluidSystem;
111 using SolidState = typename Traits::SolidState;
113 using SolidSystem = typename Traits::SolidSystem;
116
118 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
120 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
121
122 // check for permissive combinations
123 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
124 static_assert(ModelTraits::numFluidComponents() == 2, "NumComponents set in the model is not two!");
125 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
126
136 template<class ElemSol, class Problem, class Element, class Scv>
137 void update(const ElemSol& elemSol, const Problem& problem, const Element& element, const Scv& scv)
138 {
139 ParentType::update(elemSol, problem, element, scv);
140 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
141
142 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
143 const auto& matParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
144
145 const int wPhaseIdx = fluidState_.wettingPhase();
146 const int nPhaseIdx = 1 - wPhaseIdx;
147
148 // relative permeabilities -> require wetting phase saturation as parameter!
149 relativePermeability_[wPhaseIdx] = MaterialLaw::krw(matParams, saturation(wPhaseIdx));
150 relativePermeability_[nPhaseIdx] = MaterialLaw::krn(matParams, saturation(wPhaseIdx));
151
152 // porosity & permeabilty
153 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
154 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
155 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
156
157 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
158 {
159 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
160 };
161
162 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
163
164 EnergyVolVars::updateEffectiveThermalConductivity();
165 }
166
180 template<class ElemSol, class Problem, class Element, class Scv>
181 void completeFluidState(const ElemSol& elemSol,
182 const Problem& problem,
183 const Element& element,
184 const Scv& scv,
187 {
188 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
189
190 const auto& priVars = elemSol[scv.localDofIndex()];
191 const auto phasePresence = priVars.state();
192
193 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
194 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
195 const auto wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
196 fluidState.setWettingPhase(wPhaseIdx);
197
198 // set the saturations
199 if (phasePresence == firstPhaseOnly)
200 {
201 fluidState.setSaturation(phase0Idx, 1.0);
202 fluidState.setSaturation(phase1Idx, 0.0);
203 }
204 else if (phasePresence == secondPhaseOnly)
205 {
206 fluidState.setSaturation(phase0Idx, 0.0);
207 fluidState.setSaturation(phase1Idx, 1.0);
208 }
209 else if (phasePresence == bothPhases)
210 {
211 if (formulation == TwoPFormulation::p0s1)
212 {
213 fluidState.setSaturation(phase1Idx, priVars[switchIdx]);
214 fluidState.setSaturation(phase0Idx, 1 - priVars[switchIdx]);
215 }
216 else
217 {
218 fluidState.setSaturation(phase0Idx, priVars[switchIdx]);
219 fluidState.setSaturation(phase1Idx, 1 - priVars[switchIdx]);
220 }
221 }
222 else
223 DUNE_THROW(Dune::InvalidStateException, "Invalid phase presence.");
224
225 // set pressures of the fluid phases
226 pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
227 if (formulation == TwoPFormulation::p0s1)
228 {
229 fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
230 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] + pc_
231 : priVars[pressureIdx] - pc_);
232 }
233 else
234 {
235 fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
236 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] - pc_
237 : priVars[pressureIdx] + pc_);
238 }
239 }
240
244 const FluidState &fluidState() const
245 { return fluidState_; }
246
250 const SolidState &solidState() const
251 { return solidState_; }
252
258 Scalar averageMolarMass(int phaseIdx) const
259 { return fluidState_.averageMolarMass(phaseIdx); }
260
267 Scalar saturation(const int phaseIdx) const
268 { return fluidState_.saturation(phaseIdx); }
269
277 Scalar massFraction(const int phaseIdx, const int compIdx) const
278 { return fluidState_.massFraction(phaseIdx, compIdx); }
279
287 Scalar moleFraction(const int phaseIdx, const int compIdx) const
288 { return fluidState_.moleFraction(phaseIdx, compIdx); }
289
296 Scalar density(const int phaseIdx) const
297 { return fluidState_.density(phaseIdx); }
298
305 Scalar viscosity(const int phaseIdx) const
306 { return fluidState_.viscosity(phaseIdx); }
307
314 Scalar molarDensity(const int phaseIdx) const
315 { return fluidState_.molarDensity(phaseIdx); }
316
323 Scalar pressure(const int phaseIdx) const
324 { return fluidState_.pressure(phaseIdx); }
325
333 Scalar temperature() const
334 { return fluidState_.temperature(/*phaseIdx=*/0); }
335
342 Scalar relativePermeability(const int phaseIdx) const
343 { return relativePermeability_[phaseIdx]; }
344
351 Scalar mobility(const int phaseIdx) const
352 { return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
353
358 Scalar capillaryPressure() const
359 { return pc_; }
360
364 Scalar porosity() const
365 { return solidState_.porosity(); }
366
370 const PermeabilityType& permeability() const
371 { return permeability_; }
372
376 [[deprecated("Will be removed after release 3.2. Use diffusionCoefficient(phaseIdx, compIIdx, compJIdx)!")]]
377 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
378 { return diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx); }
379
383 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
384 {
385 typename FluidSystem::ParameterCache paramCache;
386 paramCache.updatePhase(fluidState_, phaseIdx);
387 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
388 }
389
393 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
394 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
395
399 int wettingPhase() const
400 { return fluidState_.wettingPhase(); }
401
402private:
403 FluidState fluidState_;
404 SolidState solidState_;
405
406 Scalar pc_; // The capillary pressure
407 PermeabilityType permeability_; // Effective permeability within the control volume
408
409 // Relative permeability within the control volume
410 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
411
412 // Effective diffusion coefficients for the phases
413 DiffusionCoefficients effectiveDiffCoeff_;
414
415protected:
416 const Impl &asImp_() const { return *static_cast<const Impl*>(this); }
417 Impl &asImp_() { return *static_cast<Impl*>(this); }
418};
425template <class Traits, bool useConstraintSolver>
426class TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>
427: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>>
428{
431
432 using Scalar = typename Traits::PrimaryVariables::value_type;
433 using ModelTraits = typename Traits::ModelTraits;
434
435 static constexpr int numFluidComps = ParentType::numFluidComponents();
436 // component indices
437 enum
438 {
439 comp0Idx = Traits::FluidSystem::comp0Idx,
440 comp1Idx = Traits::FluidSystem::comp1Idx,
441 phase0Idx = Traits::FluidSystem::phase0Idx,
442 phase1Idx = Traits::FluidSystem::phase1Idx
443 };
444
445 // phase presence indices
446 enum
447 {
448 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
449 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
450 bothPhases = ModelTraits::Indices::bothPhases
451 };
452
453 // primary variable indices
454 enum
455 {
456 switchIdx = ModelTraits::Indices::switchIdx,
457 pressureIdx = ModelTraits::Indices::pressureIdx
458 };
459
460 // formulations
461 static constexpr auto formulation = ModelTraits::priVarFormulation();
462
465public:
467 using FluidState = typename Traits::FluidState;
469 using FluidSystem = typename Traits::FluidSystem;
471 using SolidState = typename Traits::SolidState;
473 using SolidSystem = typename Traits::SolidSystem;
476
478 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
480 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
481
482 // check for permissive combinations
483 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
484
485 // The computations in the explicit composition update most probably assume a liquid-gas interface with
486 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
487 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
488 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
489
503 template<class ElemSol, class Problem, class Element, class Scv>
504 void completeFluidState(const ElemSol& elemSol,
505 const Problem& problem,
506 const Element& element,
507 const Scv& scv,
508 FluidState& fluidState,
509 SolidState& solidState)
510 {
511 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
512
513 const auto& priVars = elemSol[scv.localDofIndex()];
514 const auto phasePresence = priVars.state();
515
516 // calculate the phase compositions
517 typename FluidSystem::ParameterCache paramCache;
518
519 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
520 if(!useConstraintSolver)
521 {
522 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
523 {
524 assert(FluidSystem::isIdealMixture(phaseIdx));
525 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
526 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
527 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
528 }
529 }
530 }
531
532 // now comes the tricky part: calculate phase compositions
533 const Scalar p0 = fluidState.pressure(phase0Idx);
534 const Scalar p1 = fluidState.pressure(phase1Idx);
535 if (phasePresence == bothPhases)
536 {
537 // both phases are present, phase compositions are a result
538 // of the equilibrium between the phases. This is the job
539 // of the "MiscibleMultiPhaseComposition" constraint solver
540 if(useConstraintSolver)
542 paramCache);
543 // ... or calculated explicitly this way ...
544 else
545 {
546 // get the partial pressure of the main component of the first phase within the
547 // second phase == vapor pressure due to equilibrium. Note that in this case the
548 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
549 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0;
550
551 // get the partial pressure of the main component of the gas phase
552 const Scalar partPressGas = p1 - partPressLiquid;
553
554 // calculate the mole fractions of the components within the nonwetting phase
555 const Scalar xnn = partPressGas / p1;
556 const Scalar xnw = partPressLiquid / p1;
557
558 // calculate the mole fractions of the components within the wetting phase
559 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
560 // (see implementation in respective fluidsystem)
561 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0);
562 const Scalar xww = 1.0 - xwn;
563
564 // set all mole fractions
565 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
566 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
567 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
568 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
569 }
570 }
571 else if (phasePresence == secondPhaseOnly)
572 {
573 // only the second phase is present, composition is stored explicitly.
574 if( useMoles() )
575 {
576 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1 - priVars[switchIdx]);
577 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
578 }
579 // setMassFraction() has only to be called 1-numComponents times
580 else
581 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
582
583 // calculate the composition of the remaining phases (as
584 // well as the densities of all phases). This is the job
585 // of the "ComputeFromReferencePhase" constraint solver
586 if (useConstraintSolver)
588 paramCache,
589 phase1Idx);
590 // ... or calculated explicitly this way ...
591 else
592 {
593 // note that the water phase is actually not existing!
594 // thus, this is used as phase switch criterion
595 const Scalar xnw = priVars[switchIdx];
596 const Scalar xnn = 1.0 - xnw;
597
598 // first, xww:
599 // xnw * pn = "actual" (hypothetical) vapor pressure
600 // fugacityCoefficient * pw = vapor pressure given by thermodynamic conditions
601 // Here, xww is not actually the mole fraction of water in the wetting phase
602 // xww is only the ratio of "actual" vapor pressure / "thermodynamic" vapor pressure
603 // If xww > 1 : gas is over-saturated with water vapor,
604 // condensation takes place (see switch criterion in model)
605 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 );
606
607 // second, xwn:
608 // partialPressure / xwn = Henry
609 // partialPressure = xnn * pn
610 // xwn = xnn * pn / Henry
611 // Henry = fugacityCoefficient * pw
612 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 );
613
614 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
615 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
616 }
617 }
618 else if (phasePresence == firstPhaseOnly)
619 {
620 // only the wetting phase is present, i.e. wetting phase
621 // composition is stored explicitly.
622 if( useMoles() ) // mole-fraction formulation
623 {
624 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
625 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
626 }
627 // setMassFraction() has only to be called 1-numComponents times
628 else // mass-fraction formulation
629 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
630
631 // calculate the composition of the remaining phases (as
632 // well as the densities of all phases). This is the job
633 // of the "ComputeFromReferencePhase" constraint solver
634 if (useConstraintSolver)
636 paramCache,
637 phase0Idx);
638 // ... or calculated explicitly this way ...
639 else
640 {
641 // note that the gas phase is actually not existing!
642 // thus, this is used as phase switch criterion
643 const Scalar xwn = priVars[switchIdx];
644
645 // first, xnw:
646 // psteam = xnw * pn = partial pressure of water in gas phase
647 // psteam = fugacityCoefficient * pw
648 const Scalar xnw = ( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 )/p1;
649
650 // second, xnn:
651 // xwn = partialPressure / Henry
652 // partialPressure = pn * xnn
653 // xwn = pn * xnn / Henry
654 // xnn = xwn * Henry / pn
655 // Henry = fugacityCoefficient * pw
656 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 )/p1;
657
658 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
659 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
660 }
661 }
662
663 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
664 {
665 // set the viscosity and desity here if constraintsolver is not used
666 if(!useConstraintSolver)
667 {
668 paramCache.updateComposition(fluidState, phaseIdx);
669 const Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
670 fluidState.setDensity(phaseIdx, rho);
671 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
672 fluidState.setMolarDensity(phaseIdx, rhoMolar);
673 }
674
675 // compute and set the enthalpy
676 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
677 fluidState.setViscosity(phaseIdx,mu);
678 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
679 fluidState.setEnthalpy(phaseIdx, h);
680 }
681 }
682
683};
691template <class Traits, bool useConstraintSolver>
692class TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>
693: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>>
694{
697
698 using Scalar = typename Traits::PrimaryVariables::value_type;
699 using ModelTraits = typename Traits::ModelTraits;
700
701 static constexpr int numFluidComps = ParentType::numFluidComponents();
702 // component indices
703 enum
704 {
705 comp0Idx = Traits::FluidSystem::comp0Idx,
706 comp1Idx = Traits::FluidSystem::comp1Idx,
707 phase0Idx = Traits::FluidSystem::phase0Idx,
708 phase1Idx = Traits::FluidSystem::phase1Idx
709 };
710
711 // phase presence indices
712 enum
713 {
714 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
715 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
716 bothPhases = ModelTraits::Indices::bothPhases
717 };
718
719 // primary variable indices
720 enum
721 {
722 switchIdx = ModelTraits::Indices::switchIdx,
723 pressureIdx = ModelTraits::Indices::pressureIdx
724 };
725
726 // formulations
727 static constexpr auto formulation = ModelTraits::priVarFormulation();
728
729 using PermeabilityType = typename Traits::PermeabilityType;
732
733public:
735 using FluidState = typename Traits::FluidState;
737 using FluidSystem = typename Traits::FluidSystem;
739 using SolidState = typename Traits::SolidState;
741 using SolidSystem = typename Traits::SolidSystem;
744
746 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
748 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
749
750 // check for permissive combinations
751 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
752 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
753
754 // The computations in the explicit composition update most probably assume a liquid-gas interface with
755 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
756 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
757 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
758
772 template<class ElemSol, class Problem, class Element, class Scv>
773 void completeFluidState(const ElemSol& elemSol,
774 const Problem& problem,
775 const Element& element,
776 const Scv& scv,
777 FluidState& fluidState,
778 SolidState& solidState)
779 {
780 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
781
782 const auto& priVars = elemSol[scv.localDofIndex()];
783
785 // set the fluid compositions
787 typename FluidSystem::ParameterCache paramCache;
788 paramCache.updateAll(fluidState);
789
790 updateMoleFraction(fluidState,
791 paramCache,
792 priVars);
793
794
795 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
796 {
797 // compute and set the enthalpy
798 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
799 fluidState.setViscosity(phaseIdx,mu);
800 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
801 fluidState.setEnthalpy(phaseIdx, h);
802 }
803 }
804
812 void updateMoleFraction(FluidState &actualFluidState,
813 typename Traits::FluidSystem::ParameterCache & paramCache,
814 const typename Traits::PrimaryVariables& priVars)
815 {
816 const auto phasePresence = priVars.state();
817
818 Scalar xwnNonEquil = 0.0;
819 Scalar xwwNonEquil = 0.0;
820 Scalar xnwNonEquil = 0.0;
821 Scalar xnnNonEquil = 0.0;
822
823 if (phasePresence == bothPhases)
824 {
825 xwnNonEquil = priVars[ModelTraits::numFluidPhases()];
826 xwwNonEquil = 1-xwnNonEquil;
827 xnwNonEquil = priVars[ModelTraits::numFluidPhases()+comp1Idx];
828
829 //if the capillary pressure is high, we need to use Kelvin's equation to reduce the vapor pressure
830 if (actualFluidState.saturation(phase0Idx) < 0.01)
831 {
832 const Scalar p1 = actualFluidState.pressure(phase1Idx);
833 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
834 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
835 }
836 xnnNonEquil = 1- xnwNonEquil;
837
838 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
839 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
840 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
841 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
842
843 //check if we need this
844 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
845 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
846 {
847 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
848 paramCache,
849 phaseIdx,
850 compIdx);
851 actualFluidState.setFugacityCoefficient(phaseIdx,
852 compIdx,
853 phi);
854 }
855
856 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
857 equilFluidState.assign(actualFluidState) ;
858 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
859 if(!useConstraintSolver)
860 {
861 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
862 {
863 assert(FluidSystem::isIdealMixture(phaseIdx));
864 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
865 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
866 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
867 }
868 }
869 }
870
871 // now comes the tricky part: calculate phase compositions
872 const Scalar p0 = equilFluidState.pressure(phase0Idx);
873 const Scalar p1 = equilFluidState.pressure(phase1Idx);
874 // both phases are present, phase compositions are a result
875 // of the equilibrium between the phases. This is the job
876 // of the "MiscibleMultiPhaseComposition" constraint solver
877 if(useConstraintSolver)
879 paramCache);
880 // ... or calculated explicitly this way ...
881 else
882 {
883 // get the partial pressure of the main component of the first phase within the
884 // second phase == vapor pressure due to equilibrium. Note that in this case the
885 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
886 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
887
888 // get the partial pressure of the main component of the gas phase
889 const Scalar partPressGas = p1 - partPressLiquid;
890
891 // calculate the mole fractions of the components within the nonwetting phase
892 const Scalar xnn = partPressGas / p1;
893 const Scalar xnw = partPressLiquid / p1;
894
895 // calculate the mole fractions of the components within the wetting phase
896 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
897 // (see implementation in respective fluidsystem)
898 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
899 const Scalar xww = 1.0 - xwn;
900
901 // set all mole fractions
902 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
903 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
904 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
905 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
906 }
907
908 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
909 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
910 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
911 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
912 }
913 }
914 }
915 else
916 {
917 DUNE_THROW(Dune::InvalidStateException, "nonequilibrium is only possible for 2 phases present ");
918 }
919 // compute densities of all phases
920 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
921 {
922 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
923 actualFluidState.setDensity(phaseIdx, rho);
924 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
925 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
926 }
927 }
928
936 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
937 {
938 return xEquil_[phaseIdx][compIdx] ;
939 }
940
941private:
942 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
943};
944
945} // end namespace Dumux
946
947#endif
The available discretization methods in Dumux.
Computes all quantities of a generic fluid state if a reference phase has been specified.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.
Defines an enumeration for the formulations accepted by the two-phase model.
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition: formulation.hh:35
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
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 phasePresence() noexcept
I/O name of phase presence.
Definition: name.hh:147
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
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:77
static void solve(FluidState &fluidState, ParameterCache &paramCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:112
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:60
static void solve(FluidState &fluidState, ParameterCache &paramCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:82
Definition: porousmediumflow/2p2c/volumevariables.hh:44
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition: porousmediumflow/2p2c/volumevariables.hh:65
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/2p2c/volumevariables.hh:258
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:181
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:370
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:364
const Impl & asImp_() const
Definition: porousmediumflow/2p2c/volumevariables.hh:416
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:109
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:351
Impl & asImp_()
Definition: porousmediumflow/2p2c/volumevariables.hh:417
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:342
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/2p2c/volumevariables.hh:399
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:377
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:137
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:393
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:333
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:107
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:358
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:383
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:314
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:244
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:120
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:296
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:111
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:267
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:118
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:113
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/2p2c/volumevariables.hh:287
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:250
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:323
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:305
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/2p2c/volumevariables.hh:277
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:467
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:471
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:504
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:469
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:473
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:478
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:480
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:773
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:737
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:739
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:748
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/2p2c/volumevariables.hh:936
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:735
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:746
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:741
void updateMoleFraction(FluidState &actualFluidState, typename Traits::FluidSystem::ParameterCache &paramCache, const typename Traits::PrimaryVariables &priVars)
Updates composition of all phases from the primary variables.
Definition: porousmediumflow/2p2c/volumevariables.hh:812
The primary variable switch controlling the phase presence state variable.
Definition: 2pnc/primaryvariableswitch.hh:41
Definition: porousmediumflow/nonisothermal/volumevariables.hh:75
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:40
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition: porousmediumflow/volumevariables.hh:52
const PrimaryVariables & priVars() const
Returns the vector of primary variables.
Definition: porousmediumflow/volumevariables.hh:76
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/volumevariables.hh:64
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.
The primary variable switch for the 2pnc model.