version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_2P2C_VOLUME_VARIABLES_HH
15#define DUMUX_2P2C_VOLUME_VARIABLES_HH
16
20
27
28namespace Dumux {
29
30// forward declaration
31template <class Traits, bool enableChemicalNonEquilibrium, bool useConstraintSolver>
33
34
40template <class Traits, bool useConstraintSolver = true>
41using TwoPTwoCVolumeVariables = TwoPTwoCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium(), useConstraintSolver>;
42
49template <class Traits, class Impl>
52, public EnergyVolumeVariables<Traits, TwoPTwoCVolumeVariables<Traits> >
53{
56
57 using Scalar = typename Traits::PrimaryVariables::value_type;
58 using ModelTraits = typename Traits::ModelTraits;
59
60 static constexpr int numFluidComps = ParentType::numFluidComponents();
61 // component indices
62 enum
63 {
64 comp0Idx = Traits::FluidSystem::comp0Idx,
65 comp1Idx = Traits::FluidSystem::comp1Idx,
66 phase0Idx = Traits::FluidSystem::phase0Idx,
67 phase1Idx = Traits::FluidSystem::phase1Idx
68 };
69
70 // phase presence indices
71 enum
72 {
73 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
74 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
75 bothPhases = ModelTraits::Indices::bothPhases
76 };
77
78 // primary variable indices
79 enum
80 {
81 switchIdx = ModelTraits::Indices::switchIdx,
82 pressureIdx = ModelTraits::Indices::pressureIdx
83 };
84
85 // formulations
86 static constexpr auto formulation = ModelTraits::priVarFormulation();
87
88 using PermeabilityType = typename Traits::PermeabilityType;
91 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
92 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
93public:
95 using FluidState = typename Traits::FluidState;
97 using FluidSystem = typename Traits::FluidSystem;
99 using Indices = typename ModelTraits::Indices;
101 using SolidState = typename Traits::SolidState;
103 using SolidSystem = typename Traits::SolidSystem;
106
108 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
110 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
111
112 // check for permissive combinations
113 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
114 static_assert(ModelTraits::numFluidComponents() == 2, "NumComponents set in the model is not two!");
115 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
116
126 template<class ElemSol, class Problem, class Element, class Scv>
127 void update(const ElemSol& elemSol, const Problem& problem, const Element& element, const Scv& scv)
128 {
129 ParentType::update(elemSol, problem, element, scv);
130 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
131
132 const auto& spatialParams = problem.spatialParams();
133 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
134
135 const int wPhaseIdx = fluidState_.wettingPhase();
136 const int nPhaseIdx = 1 - wPhaseIdx;
137
138 // relative permeabilities -> require wetting phase saturation as parameter!
139 relativePermeability_[wPhaseIdx] = fluidMatrixInteraction.krw(saturation(wPhaseIdx));
140 relativePermeability_[nPhaseIdx] = fluidMatrixInteraction.krn(saturation(wPhaseIdx));
141
142 // porosity & permeabilty
143 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
144 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
145 permeability_ = spatialParams.permeability(element, scv, elemSol);
146
147 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
148 {
149 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
150 };
151
152 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
153
154 EnergyVolVars::updateEffectiveThermalConductivity();
155 }
156
170 template<class ElemSol, class Problem, class Element, class Scv>
171 void completeFluidState(const ElemSol& elemSol,
172 const Problem& problem,
173 const Element& element,
174 const Scv& scv,
177 {
178 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
179
180 const auto& priVars = elemSol[scv.localDofIndex()];
181 const auto phasePresence = priVars.state();
182
183 const auto& spatialParams = problem.spatialParams();
184 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
185 const auto wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
186 fluidState.setWettingPhase(wPhaseIdx);
187
188 // set the saturations
189 if (phasePresence == firstPhaseOnly)
190 {
191 fluidState.setSaturation(phase0Idx, 1.0);
192 fluidState.setSaturation(phase1Idx, 0.0);
193 }
194 else if (phasePresence == secondPhaseOnly)
195 {
196 fluidState.setSaturation(phase0Idx, 0.0);
197 fluidState.setSaturation(phase1Idx, 1.0);
198 }
199 else if (phasePresence == bothPhases)
200 {
201 if (formulation == TwoPFormulation::p0s1)
202 {
203 fluidState.setSaturation(phase1Idx, priVars[switchIdx]);
204 fluidState.setSaturation(phase0Idx, 1 - priVars[switchIdx]);
205 }
206 else
207 {
208 fluidState.setSaturation(phase0Idx, priVars[switchIdx]);
209 fluidState.setSaturation(phase1Idx, 1 - priVars[switchIdx]);
210 }
211 }
212 else
213 DUNE_THROW(Dune::InvalidStateException, "Invalid phase presence.");
214
215 // set pressures of the fluid phases
216 pc_ = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
217 if (formulation == TwoPFormulation::p0s1)
218 {
219 fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
220 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] + pc_
221 : priVars[pressureIdx] - pc_);
222 }
223 else
224 {
225 fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
226 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] - pc_
227 : priVars[pressureIdx] + pc_);
228 }
229 }
230
234 const FluidState &fluidState() const
235 { return fluidState_; }
236
240 const SolidState &solidState() const
241 { return solidState_; }
242
248 Scalar averageMolarMass(int phaseIdx) const
249 { return fluidState_.averageMolarMass(phaseIdx); }
250
257 Scalar saturation(const int phaseIdx) const
258 { return fluidState_.saturation(phaseIdx); }
259
267 Scalar massFraction(const int phaseIdx, const int compIdx) const
268 { return fluidState_.massFraction(phaseIdx, compIdx); }
269
277 Scalar moleFraction(const int phaseIdx, const int compIdx) const
278 { return fluidState_.moleFraction(phaseIdx, compIdx); }
279
286 Scalar density(const int phaseIdx) const
287 { return fluidState_.density(phaseIdx); }
288
295 Scalar viscosity(const int phaseIdx) const
296 { return fluidState_.viscosity(phaseIdx); }
297
304 Scalar molarDensity(const int phaseIdx) const
305 { return fluidState_.molarDensity(phaseIdx); }
306
313 Scalar pressure(const int phaseIdx) const
314 { return fluidState_.pressure(phaseIdx); }
315
323 Scalar temperature() const
324 { return fluidState_.temperature(/*phaseIdx=*/0); }
325
332 Scalar relativePermeability(const int phaseIdx) const
333 { return relativePermeability_[phaseIdx]; }
334
341 Scalar mobility(const int phaseIdx) const
342 { return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
343
348 Scalar capillaryPressure() const
349 { return pc_; }
350
354 Scalar porosity() const
355 { return solidState_.porosity(); }
356
360 const PermeabilityType& permeability() const
361 { return permeability_; }
362
366 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
367 {
368 typename FluidSystem::ParameterCache paramCache;
369 paramCache.updatePhase(fluidState_, phaseIdx);
370 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
371 }
372
376 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
377 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
378
382 int wettingPhase() const
383 { return fluidState_.wettingPhase(); }
384
385private:
386 FluidState fluidState_;
387 SolidState solidState_;
388
389 Scalar pc_; // The capillary pressure
390 PermeabilityType permeability_; // Effective permeability within the control volume
391
392 // Relative permeability within the control volume
393 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
394
395 // Effective diffusion coefficients for the phases
396 DiffusionCoefficients effectiveDiffCoeff_;
397
398protected:
399 const Impl &asImp_() const { return *static_cast<const Impl*>(this); }
400 Impl &asImp_() { return *static_cast<Impl*>(this); }
401};
408template <class Traits, bool useConstraintSolver>
409class TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>
410: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>>
411{
414
415 using Scalar = typename Traits::PrimaryVariables::value_type;
416 using ModelTraits = typename Traits::ModelTraits;
417
418 static constexpr int numFluidComps = ParentType::numFluidComponents();
419 // component indices
420 enum
421 {
422 comp0Idx = Traits::FluidSystem::comp0Idx,
423 comp1Idx = Traits::FluidSystem::comp1Idx,
424 phase0Idx = Traits::FluidSystem::phase0Idx,
425 phase1Idx = Traits::FluidSystem::phase1Idx
426 };
427
428 // phase presence indices
429 enum
430 {
431 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
432 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
433 bothPhases = ModelTraits::Indices::bothPhases
434 };
435
436 // primary variable indices
437 enum
438 {
439 switchIdx = ModelTraits::Indices::switchIdx,
440 pressureIdx = ModelTraits::Indices::pressureIdx
441 };
442
443 // formulations
444 static constexpr auto formulation = ModelTraits::priVarFormulation();
445
448public:
450 using FluidState = typename Traits::FluidState;
452 using FluidSystem = typename Traits::FluidSystem;
454 using SolidState = typename Traits::SolidState;
456 using SolidSystem = typename Traits::SolidSystem;
459
461 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
463 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
464
465 // check for permissive combinations
466 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
467
468 // The computations in the explicit composition update most probably assume a liquid-gas interface with
469 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
470 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
471 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
472
486 template<class ElemSol, class Problem, class Element, class Scv>
487 void completeFluidState(const ElemSol& elemSol,
488 const Problem& problem,
489 const Element& element,
490 const Scv& scv,
491 FluidState& fluidState,
492 SolidState& solidState)
493 {
494 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
495
496 const auto& priVars = elemSol[scv.localDofIndex()];
497 const auto phasePresence = priVars.state();
498
499 // calculate the phase compositions
500 typename FluidSystem::ParameterCache paramCache;
501
502 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
503 if(!useConstraintSolver)
504 {
505 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
506 {
507 assert(FluidSystem::isIdealMixture(phaseIdx));
508 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
509 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
510 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
511 }
512 }
513 }
514
515 // now comes the tricky part: calculate phase compositions
516 const Scalar p0 = fluidState.pressure(phase0Idx);
517 const Scalar p1 = fluidState.pressure(phase1Idx);
518 if (phasePresence == bothPhases)
519 {
520 // both phases are present, phase compositions are a result
521 // of the equilibrium between the phases. This is the job
522 // of the "MiscibleMultiPhaseComposition" constraint solver
523 if(useConstraintSolver)
525 paramCache);
526 // ... or calculated explicitly this way ...
527 else
528 {
529 // get the partial pressure of the main component of the first phase within the
530 // second phase == vapor pressure due to equilibrium. Note that in this case the
531 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
532 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0;
533
534 // get the partial pressure of the main component of the gas phase
535 const Scalar partPressGas = p1 - partPressLiquid;
536
537 // calculate the mole fractions of the components within the nonwetting phase
538 const Scalar xnn = partPressGas / p1;
539 const Scalar xnw = partPressLiquid / p1;
540
541 // calculate the mole fractions of the components within the wetting phase
542 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
543 // (see implementation in respective fluidsystem)
544 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0);
545 const Scalar xww = 1.0 - xwn;
546
547 // set all mole fractions
548 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
549 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
550 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
551 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
552 }
553 }
554 else if (phasePresence == secondPhaseOnly)
555 {
556 // only the second phase is present, composition is stored explicitly.
557 if( useMoles() )
558 {
559 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1 - priVars[switchIdx]);
560 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
561 }
562 // setMassFraction() has only to be called 1-numComponents times
563 else
564 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
565
566 // calculate the composition of the remaining phases (as
567 // well as the densities of all phases). This is the job
568 // of the "ComputeFromReferencePhase" constraint solver
569 if (useConstraintSolver)
571 paramCache,
572 phase1Idx);
573 // ... or calculated explicitly this way ...
574 else
575 {
576 // note that the water phase is actually not existing!
577 // thus, this is used as phase switch criterion
578 const Scalar xnw = priVars[switchIdx];
579 const Scalar xnn = 1.0 - xnw;
580
581 // first, xww:
582 // xnw * pn = "actual" (hypothetical) vapor pressure
583 // fugacityCoefficient * pw = vapor pressure given by thermodynamic conditions
584 // Here, xww is not actually the mole fraction of water in the wetting phase
585 // xww is only the ratio of "actual" vapor pressure / "thermodynamic" vapor pressure
586 // If xww > 1 : gas is over-saturated with water vapor,
587 // condensation takes place (see switch criterion in model)
588 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 );
589
590 // second, xwn:
591 // partialPressure / xwn = Henry
592 // partialPressure = xnn * pn
593 // xwn = xnn * pn / Henry
594 // Henry = fugacityCoefficient * pw
595 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 );
596
597 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
598 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
599 }
600 }
601 else if (phasePresence == firstPhaseOnly)
602 {
603 // only the wetting phase is present, i.e. wetting phase
604 // composition is stored explicitly.
605 if( useMoles() ) // mole-fraction formulation
606 {
607 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
608 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
609 }
610 // setMassFraction() has only to be called 1-numComponents times
611 else // mass-fraction formulation
612 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
613
614 // calculate the composition of the remaining phases (as
615 // well as the densities of all phases). This is the job
616 // of the "ComputeFromReferencePhase" constraint solver
617 if (useConstraintSolver)
619 paramCache,
620 phase0Idx);
621 // ... or calculated explicitly this way ...
622 else
623 {
624 // note that the gas phase is actually not existing!
625 // thus, this is used as phase switch criterion
626 const Scalar xwn = priVars[switchIdx];
627
628 // first, xnw:
629 // psteam = xnw * pn = partial pressure of water in gas phase
630 // psteam = fugacityCoefficient * pw
631 const Scalar xnw = ( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 )/p1;
632
633 // second, xnn:
634 // xwn = partialPressure / Henry
635 // partialPressure = pn * xnn
636 // xwn = pn * xnn / Henry
637 // xnn = xwn * Henry / pn
638 // Henry = fugacityCoefficient * pw
639 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 )/p1;
640
641 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
642 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
643 }
644 }
645
646 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
647 {
648 // set the viscosity and desity here if constraintsolver is not used
649 if(!useConstraintSolver)
650 {
651 paramCache.updateComposition(fluidState, phaseIdx);
652 const Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
653 fluidState.setDensity(phaseIdx, rho);
654 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
655 fluidState.setMolarDensity(phaseIdx, rhoMolar);
656 }
657
658 // compute and set the enthalpy
659 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
660 fluidState.setViscosity(phaseIdx,mu);
661 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
662 fluidState.setEnthalpy(phaseIdx, h);
663 }
664 }
665
666};
674template <class Traits, bool useConstraintSolver>
675class TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>
676: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>>
677{
680
681 using Scalar = typename Traits::PrimaryVariables::value_type;
682 using ModelTraits = typename Traits::ModelTraits;
683
684 static constexpr int numFluidComps = ParentType::numFluidComponents();
685 // component indices
686 enum
687 {
688 comp0Idx = Traits::FluidSystem::comp0Idx,
689 comp1Idx = Traits::FluidSystem::comp1Idx,
690 phase0Idx = Traits::FluidSystem::phase0Idx,
691 phase1Idx = Traits::FluidSystem::phase1Idx
692 };
693
694 // phase presence indices
695 enum
696 {
697 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
698 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
699 bothPhases = ModelTraits::Indices::bothPhases
700 };
701
702 // primary variable indices
703 enum
704 {
705 switchIdx = ModelTraits::Indices::switchIdx,
706 pressureIdx = ModelTraits::Indices::pressureIdx
707 };
708
709 // formulations
710 static constexpr auto formulation = ModelTraits::priVarFormulation();
711
712 using PermeabilityType = typename Traits::PermeabilityType;
715
716public:
718 using FluidState = typename Traits::FluidState;
720 using FluidSystem = typename Traits::FluidSystem;
722 using SolidState = typename Traits::SolidState;
724 using SolidSystem = typename Traits::SolidSystem;
727
729 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
731 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
732
733 // check for permissive combinations
734 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
735 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
736
737 // The computations in the explicit composition update most probably assume a liquid-gas interface with
738 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
739 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
740 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
741
755 template<class ElemSol, class Problem, class Element, class Scv>
756 void completeFluidState(const ElemSol& elemSol,
757 const Problem& problem,
758 const Element& element,
759 const Scv& scv,
760 FluidState& fluidState,
761 SolidState& solidState)
762 {
763 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
764
765 const auto& priVars = elemSol[scv.localDofIndex()];
766
768 // set the fluid compositions
770 typename FluidSystem::ParameterCache paramCache;
771 paramCache.updateAll(fluidState);
772
773 updateMoleFraction(fluidState,
774 paramCache,
775 priVars);
776
777
778 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
779 {
780 // compute and set the enthalpy
781 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
782 fluidState.setViscosity(phaseIdx,mu);
783 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
784 fluidState.setEnthalpy(phaseIdx, h);
785 }
786 }
787
795 void updateMoleFraction(FluidState &actualFluidState,
796 typename Traits::FluidSystem::ParameterCache & paramCache,
797 const typename Traits::PrimaryVariables& priVars)
798 {
799 const auto phasePresence = priVars.state();
800
801 Scalar xwnNonEquil = 0.0;
802 Scalar xwwNonEquil = 0.0;
803 Scalar xnwNonEquil = 0.0;
804 Scalar xnnNonEquil = 0.0;
805
806 if (phasePresence == bothPhases)
807 {
808 xwnNonEquil = priVars[ModelTraits::numFluidPhases()];
809 xwwNonEquil = 1-xwnNonEquil;
810 xnwNonEquil = priVars[ModelTraits::numFluidPhases()+comp1Idx];
811
812 //if the capillary pressure is high, we need to use Kelvin's equation to reduce the vapor pressure
813 if (actualFluidState.saturation(phase0Idx) < 0.01)
814 {
815 const Scalar p1 = actualFluidState.pressure(phase1Idx);
816 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
817 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
818 }
819 xnnNonEquil = 1- xnwNonEquil;
820
821 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
822 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
823 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
824 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
825
826 //check if we need this
827 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
828 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
829 {
830 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
831 paramCache,
832 phaseIdx,
833 compIdx);
834 actualFluidState.setFugacityCoefficient(phaseIdx,
835 compIdx,
836 phi);
837 }
838
839 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
840 equilFluidState.assign(actualFluidState) ;
841 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
842 if(!useConstraintSolver)
843 {
844 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
845 {
846 assert(FluidSystem::isIdealMixture(phaseIdx));
847 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
848 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
849 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
850 }
851 }
852 }
853
854 // now comes the tricky part: calculate phase compositions
855 const Scalar p0 = equilFluidState.pressure(phase0Idx);
856 const Scalar p1 = equilFluidState.pressure(phase1Idx);
857 // both phases are present, phase compositions are a result
858 // of the equilibrium between the phases. This is the job
859 // of the "MiscibleMultiPhaseComposition" constraint solver
860 if(useConstraintSolver)
862 paramCache);
863 // ... or calculated explicitly this way ...
864 else
865 {
866 // get the partial pressure of the main component of the first phase within the
867 // second phase == vapor pressure due to equilibrium. Note that in this case the
868 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
869 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
870
871 // get the partial pressure of the main component of the gas phase
872 const Scalar partPressGas = p1 - partPressLiquid;
873
874 // calculate the mole fractions of the components within the nonwetting phase
875 const Scalar xnn = partPressGas / p1;
876 const Scalar xnw = partPressLiquid / p1;
877
878 // calculate the mole fractions of the components within the wetting phase
879 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
880 // (see implementation in respective fluidsystem)
881 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
882 const Scalar xww = 1.0 - xwn;
883
884 // set all mole fractions
885 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
886 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
887 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
888 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
889 }
890
891 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
892 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
893 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
894 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
895 }
896 }
897 }
898 else
899 {
900 DUNE_THROW(Dune::InvalidStateException, "nonequilibrium is only possible for 2 phases present ");
901 }
902 // compute densities of all phases
903 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
904 {
905 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
906 actualFluidState.setDensity(phaseIdx, rho);
907 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
908 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
909 }
910 }
911
919 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
920 {
921 return xEquil_[phaseIdx][compIdx] ;
922 }
923
924private:
925 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
926};
927
928} // end namespace Dumux
929
930#endif
The primary variable switch for the 2pnc model.
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:64
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:99
Definition: porousmediumflow/nonisothermal/volumevariables.hh:63
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:47
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:69
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:28
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition: porousmediumflow/volumevariables.hh:40
const PrimaryVariables & priVars() const
Returns the vector of primary variables.
Definition: porousmediumflow/volumevariables.hh:64
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:52
The primary variable switch controlling the phase presence state variable.
Definition: 2pnc/primaryvariableswitch.hh:29
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition: porousmediumflow/2p2c/volumevariables.hh:53
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/2p2c/volumevariables.hh:248
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:171
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:360
typename ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/2p2c/volumevariables.hh:99
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:354
const Impl & asImp_() const
Definition: porousmediumflow/2p2c/volumevariables.hh:399
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:97
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:341
Impl & asImp_()
Definition: porousmediumflow/2p2c/volumevariables.hh:400
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:332
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/2p2c/volumevariables.hh:382
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:127
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:376
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:323
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:95
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:348
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:366
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:304
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:234
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:110
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:286
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:101
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:257
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:108
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:103
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:277
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:240
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:313
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:295
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:267
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:450
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:454
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:487
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:452
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:456
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:461
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:463
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:756
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:720
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:722
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:731
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:919
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:718
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:729
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:724
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:795
Definition: porousmediumflow/2p2c/volumevariables.hh:32
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Computes all quantities of a generic fluid state if a reference phase has been specified.
Defines an enumeration for the formulations accepted by the two-phase 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:24
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition: formulation.hh:23
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
The available discretization methods in Dumux.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
std::string phasePresence() noexcept
I/O name of phase presence.
Definition: name.hh:135
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:62
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
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.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.