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