3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/2p2c/volumevariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_2P2C_VOLUME_VARIABLES_HH
27#define DUMUX_2P2C_VOLUME_VARIABLES_HH
28
32
39
40namespace Dumux {
41
42// forward declaration
43template <class Traits, bool enableChemicalNonEquilibrium, bool useConstraintSolver>
45
46
52template <class Traits, bool useConstraintSolver = true>
53using TwoPTwoCVolumeVariables = TwoPTwoCVolumeVariablesImplementation<Traits, Traits::ModelTraits::enableChemicalNonEquilibrium(), useConstraintSolver>;
54
61template <class Traits, class Impl>
64, public EnergyVolumeVariables<Traits, TwoPTwoCVolumeVariables<Traits> >
65{
68
69 using Scalar = typename Traits::PrimaryVariables::value_type;
70 using ModelTraits = typename Traits::ModelTraits;
71
72 static constexpr int numFluidComps = ParentType::numFluidComponents();
73 // component indices
74 enum
75 {
76 comp0Idx = Traits::FluidSystem::comp0Idx,
77 comp1Idx = Traits::FluidSystem::comp1Idx,
78 phase0Idx = Traits::FluidSystem::phase0Idx,
79 phase1Idx = Traits::FluidSystem::phase1Idx
80 };
81
82 // phase presence indices
83 enum
84 {
85 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
86 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
87 bothPhases = ModelTraits::Indices::bothPhases
88 };
89
90 // primary variable indices
91 enum
92 {
93 switchIdx = ModelTraits::Indices::switchIdx,
94 pressureIdx = ModelTraits::Indices::pressureIdx
95 };
96
97 // formulations
98 static constexpr auto formulation = ModelTraits::priVarFormulation();
99
100 using PermeabilityType = typename Traits::PermeabilityType;
103public:
105 using FluidState = typename Traits::FluidState;
107 using FluidSystem = typename Traits::FluidSystem;
109 using SolidState = typename Traits::SolidState;
111 using SolidSystem = typename Traits::SolidSystem;
114
116 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
118 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
119
120 // check for permissive combinations
121 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
122 static_assert(ModelTraits::numFluidComponents() == 2, "NumComponents set in the model is not two!");
123 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
124
134 template<class ElemSol, class Problem, class Element, class Scv>
135 void update(const ElemSol& elemSol, const Problem& problem, const Element& element, const Scv& scv)
136 {
137 ParentType::update(elemSol, problem, element, scv);
138 asImp_().completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
139
140 // Second instance of a parameter cache. Could be avoided if
141 // diffusion coefficients also became part of the fluid state.
142 typename FluidSystem::ParameterCache paramCache;
143 paramCache.updateAll(fluidState_);
144
145 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
146 const auto& matParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
147
148 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
149 const int nPhaseIdx = 1 - wPhaseIdx;
150
151 // relative permeabilities -> require wetting phase saturation as parameter!
152 relativePermeability_[wPhaseIdx] = MaterialLaw::krw(matParams, saturation(wPhaseIdx));
153 relativePermeability_[nPhaseIdx] = MaterialLaw::krn(matParams, saturation(wPhaseIdx));
154
155 // binary diffusion coefficients
156 diffCoeff_[phase0Idx] = FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phase0Idx, comp0Idx, comp1Idx);
157 diffCoeff_[phase1Idx] = FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phase1Idx, comp0Idx, comp1Idx);
158
159 // porosity & permeabilty
160 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
161 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
162 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
163 }
164
178 template<class ElemSol, class Problem, class Element, class Scv>
179 void completeFluidState(const ElemSol& elemSol,
180 const Problem& problem,
181 const Element& element,
182 const Scv& scv,
185 {
186 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
187
188 const auto& priVars = elemSol[scv.localDofIndex()];
189 const auto phasePresence = priVars.state();
190
191 using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
192 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
193 const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
194 fluidState.setWettingPhase(wPhaseIdx);
195
196 // set the saturations
197 if (phasePresence == firstPhaseOnly)
198 {
199 fluidState.setSaturation(phase0Idx, 1.0);
200 fluidState.setSaturation(phase1Idx, 0.0);
201 }
202 else if (phasePresence == secondPhaseOnly)
203 {
204 fluidState.setSaturation(phase0Idx, 0.0);
205 fluidState.setSaturation(phase1Idx, 1.0);
206 }
207 else if (phasePresence == bothPhases)
208 {
209 if (formulation == TwoPFormulation::p0s1)
210 {
211 fluidState.setSaturation(phase1Idx, priVars[switchIdx]);
212 fluidState.setSaturation(phase0Idx, 1 - priVars[switchIdx]);
213 }
214 else
215 {
216 fluidState.setSaturation(phase0Idx, priVars[switchIdx]);
217 fluidState.setSaturation(phase1Idx, 1 - priVars[switchIdx]);
218 }
219 }
220 else
221 DUNE_THROW(Dune::InvalidStateException, "Invalid phase presence.");
222
223 // set pressures of the fluid phases
224 pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
225 if (formulation == TwoPFormulation::p0s1)
226 {
227 fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
228 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] + pc_
229 : priVars[pressureIdx] - pc_);
230 }
231 else
232 {
233 fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
234 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] - pc_
235 : priVars[pressureIdx] + pc_);
236 }
237 }
238
242 const FluidState &fluidState() const
243 { return fluidState_; }
244
248 const SolidState &solidState() const
249 { return solidState_; }
250
256 Scalar averageMolarMass(int phaseIdx) const
257 { return fluidState_.averageMolarMass(phaseIdx); }
258
265 Scalar saturation(const int phaseIdx) const
266 { return fluidState_.saturation(phaseIdx); }
267
275 Scalar massFraction(const int phaseIdx, const int compIdx) const
276 { return fluidState_.massFraction(phaseIdx, compIdx); }
277
285 Scalar moleFraction(const int phaseIdx, const int compIdx) const
286 { return fluidState_.moleFraction(phaseIdx, compIdx); }
287
294 Scalar density(const int phaseIdx) const
295 { return fluidState_.density(phaseIdx); }
296
303 Scalar viscosity(const int phaseIdx) const
304 { return fluidState_.viscosity(phaseIdx); }
305
312 Scalar molarDensity(const int phaseIdx) const
313 { return fluidState_.molarDensity(phaseIdx); }
314
321 Scalar pressure(const int phaseIdx) const
322 { return fluidState_.pressure(phaseIdx); }
323
331 Scalar temperature() const
332 { return fluidState_.temperature(/*phaseIdx=*/0); }
333
340 Scalar relativePermeability(const int phaseIdx) const
341 { return relativePermeability_[phaseIdx]; }
342
349 Scalar mobility(const int phaseIdx) const
350 { return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
351
356 Scalar capillaryPressure() const
357 { return pc_; }
358
362 Scalar porosity() const
363 { return solidState_.porosity(); }
364
368 const PermeabilityType& permeability() const
369 { return permeability_; }
370
374 Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
375 {
376 if(phaseIdx == compIdx)
377 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficient called for phaseIdx = compIdx");
378 else
379 return diffCoeff_[phaseIdx];
380 }
381
382private:
383 FluidState fluidState_;
384 SolidState solidState_;
385
386 Scalar pc_; // The capillary pressure
387 PermeabilityType permeability_; // Effective permeability within the control volume
388
389 // Relative permeability within the control volume
390 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
391
392 // Binary diffusion coefficients for the phases
393 std::array<Scalar, ModelTraits::numFluidPhases()> diffCoeff_;
394
395protected:
396 const Impl &asImp_() const { return *static_cast<const Impl*>(this); }
397 Impl &asImp_() { return *static_cast<Impl*>(this); }
398};
405template <class Traits, bool useConstraintSolver>
406class TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>
407: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, false, useConstraintSolver>>
408{
411
412 using Scalar = typename Traits::PrimaryVariables::value_type;
413 using ModelTraits = typename Traits::ModelTraits;
414
415 static constexpr int numFluidComps = ParentType::numFluidComponents();
416 // component indices
417 enum
418 {
419 comp0Idx = Traits::FluidSystem::comp0Idx,
420 comp1Idx = Traits::FluidSystem::comp1Idx,
421 phase0Idx = Traits::FluidSystem::phase0Idx,
422 phase1Idx = Traits::FluidSystem::phase1Idx
423 };
424
425 // phase presence indices
426 enum
427 {
428 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
429 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
430 bothPhases = ModelTraits::Indices::bothPhases
431 };
432
433 // primary variable indices
434 enum
435 {
436 switchIdx = ModelTraits::Indices::switchIdx,
437 pressureIdx = ModelTraits::Indices::pressureIdx
438 };
439
440 // formulations
441 static constexpr auto formulation = ModelTraits::priVarFormulation();
442
445public:
447 using FluidState = typename Traits::FluidState;
449 using FluidSystem = typename Traits::FluidSystem;
451 using SolidState = typename Traits::SolidState;
453 using SolidSystem = typename Traits::SolidSystem;
456
458 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
460 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
461
462 // check for permissive combinations
463 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
464
465 // The computations in the explicit composition update most probably assume a liquid-gas interface with
466 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
467 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
468 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
469
483 template<class ElemSol, class Problem, class Element, class Scv>
484 void completeFluidState(const ElemSol& elemSol,
485 const Problem& problem,
486 const Element& element,
487 const Scv& scv,
488 FluidState& fluidState,
489 SolidState& solidState)
490 {
491 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
492
493 const auto& priVars = elemSol[scv.localDofIndex()];
494 const auto phasePresence = priVars.state();
495
496 // calculate the phase compositions
497 typename FluidSystem::ParameterCache paramCache;
498
499 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
500 if(!useConstraintSolver)
501 {
502 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
503 {
504 assert(FluidSystem::isIdealMixture(phaseIdx));
505 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
506 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
507 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
508 }
509 }
510 }
511
512 // now comes the tricky part: calculate phase compositions
513 const Scalar p0 = fluidState.pressure(phase0Idx);
514 const Scalar p1 = fluidState.pressure(phase1Idx);
515 if (phasePresence == bothPhases)
516 {
517 // both phases are present, phase compositions are a result
518 // of the equilibrium between the phases. This is the job
519 // of the "MiscibleMultiPhaseComposition" constraint solver
520 if(useConstraintSolver)
522 paramCache);
523 // ... or calculated explicitly this way ...
524 else
525 {
526 // get the partial pressure of the main component of the first phase within the
527 // second phase == vapor pressure due to equilibrium. Note that in this case the
528 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
529 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0;
530
531 // get the partial pressure of the main component of the gas phase
532 const Scalar partPressGas = p1 - partPressLiquid;
533
534 // calculate the mole fractions of the components within the nonwetting phase
535 const Scalar xnn = partPressGas / p1;
536 const Scalar xnw = partPressLiquid / p1;
537
538 // calculate the mole fractions of the components within the wetting phase
539 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
540 // (see implementation in respective fluidsystem)
541 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0);
542 const Scalar xww = 1.0 - xwn;
543
544 // set all mole fractions
545 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
546 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
547 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
548 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
549 }
550 }
551 else if (phasePresence == secondPhaseOnly)
552 {
553 // only the second phase is present, composition is stored explicitly.
554 if( useMoles() )
555 {
556 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1 - priVars[switchIdx]);
557 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
558 }
559 // setMassFraction() has only to be called 1-numComponents times
560 else
561 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
562
563 // calculate the composition of the remaining phases (as
564 // well as the densities of all phases). This is the job
565 // of the "ComputeFromReferencePhase" constraint solver
566 if (useConstraintSolver)
568 paramCache,
569 phase1Idx);
570 // ... or calculated explicitly this way ...
571 else
572 {
573 // note that the water phase is actually not existing!
574 // thus, this is used as phase switch criterion
575 const Scalar xnw = priVars[switchIdx];
576 const Scalar xnn = 1.0 - xnw;
577
578 // first, xww:
579 // xnw * pn = "actual" (hypothetical) vapor pressure
580 // fugacityCoefficient * pw = vapor pressure given by thermodynamic conditions
581 // Here, xww is not actually the mole fraction of water in the wetting phase
582 // xww is only the ratio of "actual" vapor pressure / "thermodynamic" vapor pressure
583 // If xww > 1 : gas is over-saturated with water vapor,
584 // condensation takes place (see switch criterion in model)
585 const Scalar xww = xnw*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 );
586
587 // second, xwn:
588 // partialPressure / xwn = Henry
589 // partialPressure = xnn * pn
590 // xwn = xnn * pn / Henry
591 // Henry = fugacityCoefficient * pw
592 const Scalar xwn = xnn*p1/( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 );
593
594 fluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
595 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
596 }
597 }
598 else if (phasePresence == firstPhaseOnly)
599 {
600 // only the wetting phase is present, i.e. wetting phase
601 // composition is stored explicitly.
602 if( useMoles() ) // mole-fraction formulation
603 {
604 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
605 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
606 }
607 // setMassFraction() has only to be called 1-numComponents times
608 else // mass-fraction formulation
609 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
610
611 // calculate the composition of the remaining phases (as
612 // well as the densities of all phases). This is the job
613 // of the "ComputeFromReferencePhase" constraint solver
614 if (useConstraintSolver)
616 paramCache,
617 phase0Idx);
618 // ... or calculated explicitly this way ...
619 else
620 {
621 // note that the gas phase is actually not existing!
622 // thus, this is used as phase switch criterion
623 const Scalar xwn = priVars[switchIdx];
624
625 // first, xnw:
626 // psteam = xnw * pn = partial pressure of water in gas phase
627 // psteam = fugacityCoefficient * pw
628 const Scalar xnw = ( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx)*p0 )/p1;
629
630 // second, xnn:
631 // xwn = partialPressure / Henry
632 // partialPressure = pn * xnn
633 // xwn = pn * xnn / Henry
634 // xnn = xwn * Henry / pn
635 // Henry = fugacityCoefficient * pw
636 const Scalar xnn = xwn*( FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx)*p0 )/p1;
637
638 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
639 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
640 }
641 }
642
643 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
644 {
645 // set the viscosity and desity here if constraintsolver is not used
646 if(!useConstraintSolver)
647 {
648 paramCache.updateComposition(fluidState, phaseIdx);
649 const Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
650 fluidState.setDensity(phaseIdx, rho);
651 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
652 fluidState.setMolarDensity(phaseIdx, rhoMolar);
653 }
654
655 // compute and set the enthalpy
656 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
657 fluidState.setViscosity(phaseIdx,mu);
658 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
659 fluidState.setEnthalpy(phaseIdx, h);
660 }
661 }
662
663};
671template <class Traits, bool useConstraintSolver>
672class TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>
673: public TwoPTwoCVolumeVariablesBase<Traits, TwoPTwoCVolumeVariablesImplementation<Traits, true, useConstraintSolver>>
674{
677
678 using Scalar = typename Traits::PrimaryVariables::value_type;
679 using ModelTraits = typename Traits::ModelTraits;
680
681 static constexpr int numFluidComps = ParentType::numFluidComponents();
682 // component indices
683 enum
684 {
685 comp0Idx = Traits::FluidSystem::comp0Idx,
686 comp1Idx = Traits::FluidSystem::comp1Idx,
687 phase0Idx = Traits::FluidSystem::phase0Idx,
688 phase1Idx = Traits::FluidSystem::phase1Idx
689 };
690
691 // phase presence indices
692 enum
693 {
694 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
695 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
696 bothPhases = ModelTraits::Indices::bothPhases
697 };
698
699 // primary variable indices
700 enum
701 {
702 switchIdx = ModelTraits::Indices::switchIdx,
703 pressureIdx = ModelTraits::Indices::pressureIdx
704 };
705
706 // formulations
707 static constexpr auto formulation = ModelTraits::priVarFormulation();
708
709 using PermeabilityType = typename Traits::PermeabilityType;
712
713public:
715 using FluidState = typename Traits::FluidState;
717 using FluidSystem = typename Traits::FluidSystem;
719 using SolidState = typename Traits::SolidState;
721 using SolidSystem = typename Traits::SolidSystem;
724
726 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
728 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
729
730 // check for permissive combinations
731 static_assert(useMoles() || (!useMoles() && useConstraintSolver), "if !UseMoles, UseConstraintSolver has to be set to true");
732 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
733
734 // The computations in the explicit composition update most probably assume a liquid-gas interface with
735 // liquid as first phase. TODO: is this really needed? The constraint solver does the job anyway, doesn't it?
736 static_assert(useConstraintSolver || (!FluidSystem::isGas(phase0Idx) && FluidSystem::isGas(phase1Idx)),
737 "Explicit composition calculation has to be re-checked for NON-liquid-gas equilibria");
738
752 template<class ElemSol, class Problem, class Element, class Scv>
753 void completeFluidState(const ElemSol& elemSol,
754 const Problem& problem,
755 const Element& element,
756 const Scv& scv,
757 FluidState& fluidState,
758 SolidState& solidState)
759 {
760 ParentType::completeFluidState(elemSol, problem, element, scv, fluidState, solidState);
761
762 const auto& priVars = elemSol[scv.localDofIndex()];
763
765 // set the fluid compositions
767 typename FluidSystem::ParameterCache paramCache;
768 paramCache.updateAll(fluidState);
769
770 updateMoleFraction(fluidState,
771 paramCache,
772 priVars);
773
774
775 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
776 {
777 // compute and set the enthalpy
778 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
779 fluidState.setViscosity(phaseIdx,mu);
780 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
781 fluidState.setEnthalpy(phaseIdx, h);
782 }
783 }
784
792 void updateMoleFraction(FluidState &actualFluidState,
793 typename Traits::FluidSystem::ParameterCache & paramCache,
794 const typename Traits::PrimaryVariables& priVars)
795 {
796 const auto phasePresence = priVars.state();
797
798 Scalar xwnNonEquil = 0.0;
799 Scalar xwwNonEquil = 0.0;
800 Scalar xnwNonEquil = 0.0;
801 Scalar xnnNonEquil = 0.0;
802
803 if (phasePresence == bothPhases)
804 {
805 xwnNonEquil = priVars[ModelTraits::numFluidPhases()];
806 xwwNonEquil = 1-xwnNonEquil;
807 xnwNonEquil = priVars[ModelTraits::numFluidPhases()+comp1Idx];
808
809 //if the capillary pressure is high, we need to use Kelvin's equation to reduce the vapor pressure
810 if (actualFluidState.saturation(phase0Idx) < 0.01)
811 {
812 const Scalar p1 = actualFluidState.pressure(phase1Idx);
813 const Scalar partPressLiquid = FluidSystem::vaporPressure(actualFluidState, comp0Idx);
814 xnwNonEquil =std::min(partPressLiquid/p1, xnwNonEquil);
815 }
816 xnnNonEquil = 1- xnwNonEquil;
817
818 actualFluidState.setMoleFraction(phase0Idx, comp0Idx, xwwNonEquil);
819 actualFluidState.setMoleFraction(phase0Idx, comp1Idx, xwnNonEquil);
820 actualFluidState.setMoleFraction(phase1Idx, comp0Idx, xnwNonEquil);
821 actualFluidState.setMoleFraction(phase1Idx, comp1Idx, xnnNonEquil);
822
823 //check if we need this
824 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
825 for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
826 {
827 const Scalar phi = FluidSystem::fugacityCoefficient(actualFluidState,
828 paramCache,
829 phaseIdx,
830 compIdx);
831 actualFluidState.setFugacityCoefficient(phaseIdx,
832 compIdx,
833 phi);
834 }
835
836 FluidState equilFluidState; // the fluidState *on the interface* i.e. chemical equilibrium
837 equilFluidState.assign(actualFluidState) ;
838 // If constraint solver is not used, get the phase pressures and set the fugacity coefficients here
839 if(!useConstraintSolver)
840 {
841 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx)
842 {
843 assert(FluidSystem::isIdealMixture(phaseIdx));
844 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
845 Scalar phi = FluidSystem::fugacityCoefficient(equilFluidState, paramCache, phaseIdx, compIdx);
846 equilFluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
847 }
848 }
849 }
850
851 // now comes the tricky part: calculate phase compositions
852 const Scalar p0 = equilFluidState.pressure(phase0Idx);
853 const Scalar p1 = equilFluidState.pressure(phase1Idx);
854 // both phases are present, phase compositions are a result
855 // of the equilibrium between the phases. This is the job
856 // of the "MiscibleMultiPhaseComposition" constraint solver
857 if(useConstraintSolver)
859 paramCache);
860 // ... or calculated explicitly this way ...
861 else
862 {
863 // get the partial pressure of the main component of the first phase within the
864 // second phase == vapor pressure due to equilibrium. Note that in this case the
865 // fugacityCoefficient * p is the vapor pressure (see implementation in respective fluidsystem)
866 const Scalar partPressLiquid = FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp0Idx)*p0;
867
868 // get the partial pressure of the main component of the gas phase
869 const Scalar partPressGas = p1 - partPressLiquid;
870
871 // calculate the mole fractions of the components within the nonwetting phase
872 const Scalar xnn = partPressGas / p1;
873 const Scalar xnw = partPressLiquid / p1;
874
875 // calculate the mole fractions of the components within the wetting phase
876 // note that in this case the fugacityCoefficient * p is the Henry Coefficient
877 // (see implementation in respective fluidsystem)
878 const Scalar xwn = partPressGas / (FluidSystem::fugacityCoefficient(equilFluidState, phase0Idx, comp1Idx)*p0);
879 const Scalar xww = 1.0 - xwn;
880
881 // set all mole fractions
882 equilFluidState.setMoleFraction(phase0Idx, comp0Idx, xww);
883 equilFluidState.setMoleFraction(phase0Idx, comp1Idx, xwn);
884 equilFluidState.setMoleFraction(phase1Idx, comp0Idx, xnw);
885 equilFluidState.setMoleFraction(phase1Idx, comp1Idx, xnn);
886 }
887
888 // Setting the equilibrium composition (in a kinetic model not necessarily the same as the actual mole fraction)
889 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx){
890 for (int compIdx=0; compIdx< numFluidComps; ++ compIdx){
891 xEquil_[phaseIdx][compIdx] = equilFluidState.moleFraction(phaseIdx, compIdx);
892 }
893 }
894 }
895 else
896 {
897 DUNE_THROW(Dune::InvalidStateException, "nonequilibrium is only possible for 2 phases present ");
898 }
899 // compute densities of all phases
900 for(int phaseIdx=0; phaseIdx<ModelTraits::numFluidPhases(); ++phaseIdx)
901 {
902 const Scalar rho = FluidSystem::density(actualFluidState, paramCache, phaseIdx);
903 actualFluidState.setDensity(phaseIdx, rho);
904 const Scalar rhoMolar = FluidSystem::molarDensity(actualFluidState, paramCache, phaseIdx);
905 actualFluidState.setMolarDensity(phaseIdx, rhoMolar);
906 }
907 }
908
916 const Scalar xEquil(const unsigned int phaseIdx, const unsigned int compIdx) const
917 {
918 return xEquil_[phaseIdx][compIdx] ;
919 }
920
921private:
922 std::array<std::array<Scalar, numFluidComps>, ModelTraits::numFluidPhases()> xEquil_;
923};
924
925} // end namespace Dumux
926
927#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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string phasePresence() noexcept
I/O name of phase presence.
Definition: name.hh:147
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:77
static void solve(FluidState &fluidState, ParameterCache &paramCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:112
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:60
static void solve(FluidState &fluidState, ParameterCache &paramCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:82
Definition: porousmediumflow/2p2c/volumevariables.hh:44
Contains the quantities which are constant within a finite volume in the two-phase two-component mode...
Definition: porousmediumflow/2p2c/volumevariables.hh:65
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/2p2c/volumevariables.hh:256
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:179
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:368
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:362
const Impl & asImp_() const
Definition: porousmediumflow/2p2c/volumevariables.hh:396
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:107
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:349
Impl & asImp_()
Definition: porousmediumflow/2p2c/volumevariables.hh:397
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:340
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2p2c/volumevariables.hh:374
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:135
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:331
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:105
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:356
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:312
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:242
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:118
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:294
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:109
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:265
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:116
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:111
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:285
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2p2c/volumevariables.hh:248
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:321
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/2p2c/volumevariables.hh:303
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:275
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:447
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:451
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:484
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:449
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:453
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:458
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:460
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:753
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:717
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2p2c/volumevariables.hh:719
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2p2c/volumevariables.hh:728
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:916
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/2p2c/volumevariables.hh:715
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2p2c/volumevariables.hh:726
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2p2c/volumevariables.hh:721
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:792
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.