3.4
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 Indices = typename ModelTraits::Indices;
113 using SolidState = typename Traits::SolidState;
115 using SolidSystem = typename Traits::SolidSystem;
118
120 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
122 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
123
124 // check for permissive combinations
125 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
126 static_assert(ModelTraits::numFluidComponents() == 2, "NumComponents set in the model is not two!");
127 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
128
138 template<class ElemSol, class Problem, class Element, class Scv>
139 void update(const ElemSol& elemSol, const Problem& problem, const Element& element, const Scv& scv)
140 {
141 ParentType::update(elemSol, problem, element, scv);
142 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
143
144 const auto& spatialParams = problem.spatialParams();
145 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
146
147 const int wPhaseIdx = fluidState_.wettingPhase();
148 const int nPhaseIdx = 1 - wPhaseIdx;
149
150 // relative permeabilities -> require wetting phase saturation as parameter!
151 relativePermeability_[wPhaseIdx] = fluidMatrixInteraction.krw(saturation(wPhaseIdx));
152 relativePermeability_[nPhaseIdx] = fluidMatrixInteraction.krn(saturation(wPhaseIdx));
153
154 // porosity & permeabilty
155 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
156 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
157 permeability_ = spatialParams.permeability(element, scv, elemSol);
158
159 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
160 {
161 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
162 };
163
164 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
165
166 EnergyVolVars::updateEffectiveThermalConductivity();
167 }
168
182 template<class ElemSol, class Problem, class Element, class Scv>
183 void completeFluidState(const ElemSol& elemSol,
184 const Problem& problem,
185 const Element& element,
186 const Scv& scv,
189 {
190 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
191
192 const auto& priVars = elemSol[scv.localDofIndex()];
193 const auto phasePresence = priVars.state();
194
195 const auto& spatialParams = problem.spatialParams();
196 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
197 const auto wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
198 fluidState.setWettingPhase(wPhaseIdx);
199
200 // set the saturations
201 if (phasePresence == firstPhaseOnly)
202 {
203 fluidState.setSaturation(phase0Idx, 1.0);
204 fluidState.setSaturation(phase1Idx, 0.0);
205 }
206 else if (phasePresence == secondPhaseOnly)
207 {
208 fluidState.setSaturation(phase0Idx, 0.0);
209 fluidState.setSaturation(phase1Idx, 1.0);
210 }
211 else if (phasePresence == bothPhases)
212 {
213 if (formulation == TwoPFormulation::p0s1)
214 {
215 fluidState.setSaturation(phase1Idx, priVars[switchIdx]);
216 fluidState.setSaturation(phase0Idx, 1 - priVars[switchIdx]);
217 }
218 else
219 {
220 fluidState.setSaturation(phase0Idx, priVars[switchIdx]);
221 fluidState.setSaturation(phase1Idx, 1 - priVars[switchIdx]);
222 }
223 }
224 else
225 DUNE_THROW(Dune::InvalidStateException, "Invalid phase presence.");
226
227 // set pressures of the fluid phases
228 pc_ = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
229 if (formulation == TwoPFormulation::p0s1)
230 {
231 fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
232 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] + pc_
233 : priVars[pressureIdx] - pc_);
234 }
235 else
236 {
237 fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
238 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] - pc_
239 : priVars[pressureIdx] + pc_);
240 }
241 }
242
246 const FluidState &fluidState() const
247 { return fluidState_; }
248
252 const SolidState &solidState() const
253 { return solidState_; }
254
260 Scalar averageMolarMass(int phaseIdx) const
261 { return fluidState_.averageMolarMass(phaseIdx); }
262
269 Scalar saturation(const int phaseIdx) const
270 { return fluidState_.saturation(phaseIdx); }
271
279 Scalar massFraction(const int phaseIdx, const int compIdx) const
280 { return fluidState_.massFraction(phaseIdx, compIdx); }
281
289 Scalar moleFraction(const int phaseIdx, const int compIdx) const
290 { return fluidState_.moleFraction(phaseIdx, compIdx); }
291
298 Scalar density(const int phaseIdx) const
299 { return fluidState_.density(phaseIdx); }
300
307 Scalar viscosity(const int phaseIdx) const
308 { return fluidState_.viscosity(phaseIdx); }
309
316 Scalar molarDensity(const int phaseIdx) const
317 { return fluidState_.molarDensity(phaseIdx); }
318
325 Scalar pressure(const int phaseIdx) const
326 { return fluidState_.pressure(phaseIdx); }
327
335 Scalar temperature() const
336 { return fluidState_.temperature(/*phaseIdx=*/0); }
337
344 Scalar relativePermeability(const int phaseIdx) const
345 { return relativePermeability_[phaseIdx]; }
346
353 Scalar mobility(const int phaseIdx) const
354 { return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
355
360 Scalar capillaryPressure() const
361 { return pc_; }
362
366 Scalar porosity() const
367 { return solidState_.porosity(); }
368
372 const PermeabilityType& permeability() const
373 { return permeability_; }
374
378 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
379 {
380 typename FluidSystem::ParameterCache paramCache;
381 paramCache.updatePhase(fluidState_, phaseIdx);
382 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
383 }
384
388 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
389 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
390
394 int wettingPhase() const
395 { return fluidState_.wettingPhase(); }
396
397private:
398 FluidState fluidState_;
399 SolidState solidState_;
400
401 Scalar pc_; // The capillary pressure
402 PermeabilityType permeability_; // Effective permeability within the control volume
403
404 // Relative permeability within the control volume
405 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
406
407 // Effective diffusion coefficients for the phases
408 DiffusionCoefficients effectiveDiffCoeff_;
409
410protected:
411 const Impl &asImp_() const { return *static_cast<const Impl*>(this); }
412 Impl &asImp_() { return *static_cast<Impl*>(this); }
413};
420template <class Traits, bool useConstraintSolver>
421class TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>
422: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>>
423{
426
427 using Scalar = typename Traits::PrimaryVariables::value_type;
428 using ModelTraits = typename Traits::ModelTraits;
429
430 static constexpr int numFluidComps = ParentType::numFluidComponents();
431 // component indices
432 enum
433 {
434 comp0Idx = Traits::FluidSystem::comp0Idx,
435 comp1Idx = Traits::FluidSystem::comp1Idx,
436 phase0Idx = Traits::FluidSystem::phase0Idx,
437 phase1Idx = Traits::FluidSystem::phase1Idx
438 };
439
440 // phase presence indices
441 enum
442 {
443 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
444 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
445 bothPhases = ModelTraits::Indices::bothPhases
446 };
447
448 // primary variable indices
449 enum
450 {
451 switchIdx = ModelTraits::Indices::switchIdx,
452 pressureIdx = ModelTraits::Indices::pressureIdx
453 };
454
455 // formulations
456 static constexpr auto formulation = ModelTraits::priVarFormulation();
457
460public:
462 using FluidState = typename Traits::FluidState;
464 using FluidSystem = typename Traits::FluidSystem;
466 using SolidState = typename Traits::SolidState;
468 using SolidSystem = typename Traits::SolidSystem;
471
473 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
475 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
476
477 // check for permissive combinations
478 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
479
480 // The computations in the explicit composition update most probably assume a liquid-gas interface with
481 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
482 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
483 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
484
498 template<class ElemSol, class Problem, class Element, class Scv>
499 void completeFluidState(const ElemSol& elemSol,
500 const Problem& problem,
501 const Element& element,
502 const Scv& scv,
503 FluidState& fluidState,
504 SolidState& solidState)
505 {
506 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
507
508 const auto& priVars = elemSol[scv.localDofIndex()];
509 const auto phasePresence = priVars.state();
510
511 // calculate the phase compositions
512 typename FluidSystem::ParameterCache paramCache;
513
514 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
515 if(!useConstraintSolver)
516 {
517 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
518 {
519 assert(FluidSystem::isIdealMixture(phaseIdx));
520 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
521 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
522 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
523 }
524 }
525 }
526
527 // now comes the tricky part: calculate phase compositions
528 const Scalar p0 = fluidState.pressure(phase0Idx);
529 const Scalar p1 = fluidState.pressure(phase1Idx);
530 if (phasePresence == bothPhases)
531 {
532 // both phases are present, phase compositions are a result
533 // of the equilibrium between the phases. This is the job
534 // of the "MiscibleMultiPhaseComposition" constraint solver
535 if(useConstraintSolver)
537 paramCache);
538 // ... or calculated explicitly this way ...
539 else
540 {
541 // get the partial pressure of the main component of the first phase within the
542 // second phase == vapor pressure due to equilibrium. Note that in this case the
543 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
544 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0;
545
546 // get the partial pressure of the main component of the gas phase
547 const Scalar partPressGas = p1 - partPressLiquid;
548
549 // calculate the mole fractions of the components within the nonwetting phase
550 const Scalar xnn = partPressGas / p1;
551 const Scalar xnw = partPressLiquid / p1;
552
553 // calculate the mole fractions of the components within the wetting phase
554 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
555 // (see implementation in respective fluidsystem)
556 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0);
557 const Scalar xww = 1.0 - xwn;
558
559 // set all mole fractions
560 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
561 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
562 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
563 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
564 }
565 }
566 else if (phasePresence == secondPhaseOnly)
567 {
568 // only the second phase is present, composition is stored explicitly.
569 if( useMoles() )
570 {
571 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1 - priVars[switchIdx]);
572 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
573 }
574 // setMassFraction() has only to be called 1-numComponents times
575 else
576 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
577
578 // calculate the composition of the remaining phases (as
579 // well as the densities of all phases). This is the job
580 // of the "ComputeFromReferencePhase" constraint solver
581 if (useConstraintSolver)
583 paramCache,
584 phase1Idx);
585 // ... or calculated explicitly this way ...
586 else
587 {
588 // note that the water phase is actually not existing!
589 // thus, this is used as phase switch criterion
590 const Scalar xnw = priVars[switchIdx];
591 const Scalar xnn = 1.0 - xnw;
592
593 // first, xww:
594 // xnw * pn = "actual" (hypothetical) vapor pressure
595 // fugacityCoefficient * pw = vapor pressure given by thermodynamic conditions
596 // Here, xww is not actually the mole fraction of water in the wetting phase
597 // xww is only the ratio of "actual" vapor pressure / "thermodynamic" vapor pressure
598 // If xww > 1 : gas is over-saturated with water vapor,
599 // condensation takes place (see switch criterion in model)
600 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 );
601
602 // second, xwn:
603 // partialPressure / xwn = Henry
604 // partialPressure = xnn * pn
605 // xwn = xnn * pn / Henry
606 // Henry = fugacityCoefficient * pw
607 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 );
608
609 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
610 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
611 }
612 }
613 else if (phasePresence == firstPhaseOnly)
614 {
615 // only the wetting phase is present, i.e. wetting phase
616 // composition is stored explicitly.
617 if( useMoles() ) // mole-fraction formulation
618 {
619 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
620 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
621 }
622 // setMassFraction() has only to be called 1-numComponents times
623 else // mass-fraction formulation
624 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
625
626 // calculate the composition of the remaining phases (as
627 // well as the densities of all phases). This is the job
628 // of the "ComputeFromReferencePhase" constraint solver
629 if (useConstraintSolver)
631 paramCache,
632 phase0Idx);
633 // ... or calculated explicitly this way ...
634 else
635 {
636 // note that the gas phase is actually not existing!
637 // thus, this is used as phase switch criterion
638 const Scalar xwn = priVars[switchIdx];
639
640 // first, xnw:
641 // psteam = xnw * pn = partial pressure of water in gas phase
642 // psteam = fugacityCoefficient * pw
643 const Scalar xnw = ( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 )/p1;
644
645 // second, xnn:
646 // xwn = partialPressure / Henry
647 // partialPressure = pn * xnn
648 // xwn = pn * xnn / Henry
649 // xnn = xwn * Henry / pn
650 // Henry = fugacityCoefficient * pw
651 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 )/p1;
652
653 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
654 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
655 }
656 }
657
658 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
659 {
660 // set the viscosity and desity here if constraintsolver is not used
661 if(!useConstraintSolver)
662 {
663 paramCache.updateComposition(fluidState, phaseIdx);
664 const Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
665 fluidState.setDensity(phaseIdx, rho);
666 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
667 fluidState.setMolarDensity(phaseIdx, rhoMolar);
668 }
669
670 // compute and set the enthalpy
671 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
672 fluidState.setViscosity(phaseIdx,mu);
673 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
674 fluidState.setEnthalpy(phaseIdx, h);
675 }
676 }
677
678};
686template <class Traits, bool useConstraintSolver>
687class TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>
688: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>>
689{
692
693 using Scalar = typename Traits::PrimaryVariables::value_type;
694 using ModelTraits = typename Traits::ModelTraits;
695
696 static constexpr int numFluidComps = ParentType::numFluidComponents();
697 // component indices
698 enum
699 {
700 comp0Idx = Traits::FluidSystem::comp0Idx,
701 comp1Idx = Traits::FluidSystem::comp1Idx,
702 phase0Idx = Traits::FluidSystem::phase0Idx,
703 phase1Idx = Traits::FluidSystem::phase1Idx
704 };
705
706 // phase presence indices
707 enum
708 {
709 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
710 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
711 bothPhases = ModelTraits::Indices::bothPhases
712 };
713
714 // primary variable indices
715 enum
716 {
717 switchIdx = ModelTraits::Indices::switchIdx,
718 pressureIdx = ModelTraits::Indices::pressureIdx
719 };
720
721 // formulations
722 static constexpr auto formulation = ModelTraits::priVarFormulation();
723
724 using PermeabilityType = typename Traits::PermeabilityType;
727
728public:
730 using FluidState = typename Traits::FluidState;
732 using FluidSystem = typename Traits::FluidSystem;
734 using SolidState = typename Traits::SolidState;
736 using SolidSystem = typename Traits::SolidSystem;
739
741 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
743 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
744
745 // check for permissive combinations
746 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
747 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
748
749 // The computations in the explicit composition update most probably assume a liquid-gas interface with
750 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
751 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
752 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
753
767 template<class ElemSol, class Problem, class Element, class Scv>
768 void completeFluidState(const ElemSol& elemSol,
769 const Problem& problem,
770 const Element& element,
771 const Scv& scv,
772 FluidState& fluidState,
773 SolidState& solidState)
774 {
775 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
776
777 const auto& priVars = elemSol[scv.localDofIndex()];
778
780 // set the fluid compositions
782 typename FluidSystem::ParameterCache paramCache;
783 paramCache.updateAll(fluidState);
784
785 updateMoleFraction(fluidState,
786 paramCache,
787 priVars);
788
789
790 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
791 {
792 // compute and set the enthalpy
793 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
794 fluidState.setViscosity(phaseIdx,mu);
795 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
796 fluidState.setEnthalpy(phaseIdx, h);
797 }
798 }
799
807 void updateMoleFraction(FluidState &actualFluidState,
808 typename Traits::FluidSystem::ParameterCache & paramCache,
809 const typename Traits::PrimaryVariables& priVars)
810 {
811 const auto phasePresence = priVars.state();
812
813 Scalar xwnNonEquil = 0.0;
814 Scalar xwwNonEquil = 0.0;
815 Scalar xnwNonEquil = 0.0;
816 Scalar xnnNonEquil = 0.0;
817
818 if (phasePresence == bothPhases)
819 {
820 xwnNonEquil = priVars[ModelTraits::numFluidPhases()];
821 xwwNonEquil = 1-xwnNonEquil;
822 xnwNonEquil = priVars[ModelTraits::numFluidPhases()+comp1Idx];
823
824 //if the capillary pressure is high, we need to use Kelvin's equation to reduce the vapor pressure
825 if (actualFluidState.saturation(phase0Idx) < 0.01)
826 {
827 const Scalar p1 = actualFluidState.pressure(phase1Idx);
828 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
829 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
830 }
831 xnnNonEquil = 1- xnwNonEquil;
832
833 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
834 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
835 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
836 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
837
838 //check if we need this
839 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
840 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
841 {
842 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
843 paramCache,
844 phaseIdx,
845 compIdx);
846 actualFluidState.setFugacityCoefficient(phaseIdx,
847 compIdx,
848 phi);
849 }
850
851 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
852 equilFluidState.assign(actualFluidState) ;
853 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
854 if(!useConstraintSolver)
855 {
856 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
857 {
858 assert(FluidSystem::isIdealMixture(phaseIdx));
859 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
860 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
861 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
862 }
863 }
864 }
865
866 // now comes the tricky part: calculate phase compositions
867 const Scalar p0 = equilFluidState.pressure(phase0Idx);
868 const Scalar p1 = equilFluidState.pressure(phase1Idx);
869 // both phases are present, phase compositions are a result
870 // of the equilibrium between the phases. This is the job
871 // of the "MiscibleMultiPhaseComposition" constraint solver
872 if(useConstraintSolver)
874 paramCache);
875 // ... or calculated explicitly this way ...
876 else
877 {
878 // get the partial pressure of the main component of the first phase within the
879 // second phase == vapor pressure due to equilibrium. Note that in this case the
880 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
881 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
882
883 // get the partial pressure of the main component of the gas phase
884 const Scalar partPressGas = p1 - partPressLiquid;
885
886 // calculate the mole fractions of the components within the nonwetting phase
887 const Scalar xnn = partPressGas / p1;
888 const Scalar xnw = partPressLiquid / p1;
889
890 // calculate the mole fractions of the components within the wetting phase
891 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
892 // (see implementation in respective fluidsystem)
893 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
894 const Scalar xww = 1.0 - xwn;
895
896 // set all mole fractions
897 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
898 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
899 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
900 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
901 }
902
903 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
904 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
905 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
906 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
907 }
908 }
909 }
910 else
911 {
912 DUNE_THROW(Dune::InvalidStateException, "nonequilibrium is only possible for 2 phases present ");
913 }
914 // compute densities of all phases
915 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
916 {
917 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
918 actualFluidState.setDensity(phaseIdx, rho);
919 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
920 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
921 }
922 }
923
931 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
932 {
933 return xEquil_[phaseIdx][compIdx] ;
934 }
935
936private:
937 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
938};
939
940} // end namespace Dumux
941
942#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:76
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:111
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:59
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:81
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:260
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:183
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:372
typename ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/2p2c/volumevariables.hh:111
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:366
const Impl & asImp_() const
Definition: porousmediumflow/2p2c/volumevariables.hh:411
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:353
Impl & asImp_()
Definition: porousmediumflow/2p2c/volumevariables.hh:412
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:344
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/2p2c/volumevariables.hh:394
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:139
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:388
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:335
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:360
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:378
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:316
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:246
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:122
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:298
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:113
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:269
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:120
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:115
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:289
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:252
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:325
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:307
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:279
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:462
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:466
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:499
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:464
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:468
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:473
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:475
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:768
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:732
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:734
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:743
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:931
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:730
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:741
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:736
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:807
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.