3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/co2/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_CO2_VOLUME_VARIABLES_HH
27#define DUMUX_CO2_VOLUME_VARIABLES_HH
28
29#include <array>
30
31#include <dune/common/exceptions.hh>
32
37
39
40namespace Dumux {
41
47template <class Traits>
50, public EnergyVolumeVariables<Traits, TwoPTwoCCO2VolumeVariables<Traits> >
51{
54
55 using Scalar = typename Traits::PrimaryVariables::value_type;
56 using ModelTraits = typename Traits::ModelTraits;
57 static constexpr int numFluidComps = ParentType::numFluidComponents();
58
59 // component indices
60 enum
61 {
62 comp0Idx = Traits::FluidSystem::comp0Idx,
63 comp1Idx = Traits::FluidSystem::comp1Idx,
64 phase0Idx = Traits::FluidSystem::phase0Idx,
65 phase1Idx = Traits::FluidSystem::phase1Idx
66 };
67
68 // phase presence indices
69 enum
70 {
71 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
72 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
73 bothPhases = ModelTraits::Indices::bothPhases
74 };
75
76 // primary variable indices
77 enum
78 {
79 switchIdx = ModelTraits::Indices::switchIdx,
80 pressureIdx = ModelTraits::Indices::pressureIdx
81 };
82
83 // formulation
84 static constexpr auto formulation = ModelTraits::priVarFormulation();
85
86 // type used for the permeability
87 using PermeabilityType = typename Traits::PermeabilityType;
88
89 // type used for the diffusion coefficients
90 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
91 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
92
93public:
95 using FluidState = typename Traits::FluidState;
97 using FluidSystem = typename Traits::FluidSystem;
99 using Indices = typename ModelTraits::Indices;
101 using SolidState = typename Traits::SolidState;
103 using SolidSystem = typename Traits::SolidSystem;
106
107
109 static constexpr bool useMoles() { return ModelTraits::useMoles(); }
111 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
112
113 // check for permissive combinations
114 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
115 static_assert(ModelTraits::numFluidComponents() == 2, "NumComponents set in the model is not two!");
116 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
117
127 template<class ElemSol, class Problem, class Element, class Scv>
128 void update(const ElemSol& elemSol, const Problem& problem, const Element& element, const Scv& scv)
129 {
130 ParentType::update(elemSol, problem, element, scv);
131 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
132
133 // Second instance of a parameter cache. Could be avoided if
134 // diffusion coefficients also became part of the fluid state.
135 typename FluidSystem::ParameterCache paramCache;
136 paramCache.updateAll(fluidState_);
137
138 const auto& spatialParams = problem.spatialParams();
139 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
140
141 const int wPhaseIdx = fluidState_.wettingPhase();
142 const int nPhaseIdx = 1 - wPhaseIdx;
143
144 // relative permeabilities -> require wetting phase saturation as parameter!
145 relativePermeability_[wPhaseIdx] = fluidMatrixInteraction.krw(saturation(wPhaseIdx));
146 relativePermeability_[nPhaseIdx] = fluidMatrixInteraction.krn(saturation(wPhaseIdx));
147
148 // porosity & permeabilty
149 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
150 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
151 permeability_ = spatialParams.permeability(element, scv, elemSol);
152
153 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
154 {
155 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
156 };
157
158 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
159
160 EnergyVolVars::updateEffectiveThermalConductivity();
161 }
162
177 template<class ElemSol, class Problem, class Element, class Scv>
178 void completeFluidState(const ElemSol& elemSol,
179 const Problem& problem,
180 const Element& element,
181 const Scv& scv,
184 {
185 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
186
187 const auto& priVars = elemSol[scv.localDofIndex()];
188 const auto phasePresence = priVars.state();
189
190 const auto& spatialParams = problem.spatialParams();
191 const auto wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
192 fluidState.setWettingPhase(wPhaseIdx);
193
194 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
195
196 // set the saturations
197 if (phasePresence == secondPhaseOnly)
198 {
199 fluidState.setSaturation(phase0Idx, 0.0);
200 fluidState.setSaturation(phase1Idx, 1.0);
201 }
202 else if (phasePresence == firstPhaseOnly)
203 {
204 fluidState.setSaturation(phase0Idx, 1.0);
205 fluidState.setSaturation(phase1Idx, 0.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_ = fluidMatrixInteraction.pc(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 // calculate the phase compositions
239 typename FluidSystem::ParameterCache paramCache;
240 // both phases are present
241 if (phasePresence == bothPhases)
242 {
243 //Get the equilibrium mole fractions from the FluidSystem and set them in the fluidState
244 //xCO2 = equilibrium mole fraction of CO2 in the liquid phase
245 //yH2O = equilibrium mole fraction of H2O in the gas phase
246 const auto xwCO2 = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase0Idx);
247 const auto xgH2O = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase1Idx);
248 const auto xwH2O = 1 - xwCO2;
249 const auto xgCO2 = 1 - xgH2O;
250 fluidState.setMoleFraction(phase0Idx, comp0Idx, xwH2O);
251 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwCO2);
252 fluidState.setMoleFraction(phase1Idx, comp0Idx, xgH2O);
253 fluidState.setMoleFraction(phase1Idx, comp1Idx, xgCO2);
254 }
255
256 // only the nonwetting phase is present, i.e. nonwetting phase
257 // composition is stored explicitly.
258 else if (phasePresence == secondPhaseOnly)
259 {
260 if( useMoles() ) // mole-fraction formulation
261 {
262 // set the fluid state
263 fluidState.setMoleFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
264 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1-priVars[switchIdx]);
265 // TODO give values for non-existing wetting phase
266 const auto xwCO2 = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase0Idx);
267 const auto xwH2O = 1 - xwCO2;
268 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwCO2);
269 fluidState.setMoleFraction(phase0Idx, comp0Idx, xwH2O);
270 }
271 else // mass-fraction formulation
272 {
273 // setMassFraction() has only to be called 1-numComponents times
274 fluidState.setMassFraction(phase1Idx, comp0Idx, priVars[switchIdx]);
275 // TODO give values for non-existing wetting phase
276 const auto xwCO2 = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase0Idx);
277 const auto xwH2O = 1 - xwCO2;
278 fluidState.setMoleFraction(phase0Idx, comp1Idx, xwCO2);
279 fluidState.setMoleFraction(phase0Idx, comp0Idx, xwH2O);
280 }
281 }
282
283 // only the wetting phase is present, i.e. wetting phase
284 // composition is stored explicitly.
285 else if (phasePresence == firstPhaseOnly)
286 {
287 if( useMoles() ) // mole-fraction formulation
288 {
289 // convert mass to mole fractions and set the fluid state
290 fluidState.setMoleFraction(phase0Idx, comp0Idx, 1-priVars[switchIdx]);
291 fluidState.setMoleFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
292 // TODO give values for non-existing nonwetting phase
293 Scalar xnH2O = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase1Idx);
294 Scalar xnCO2 = 1 - xnH2O;
295 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnCO2);
296 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnH2O);
297 }
298 else // mass-fraction formulation
299 {
300 // setMassFraction() has only to be called 1-numComponents times
301 fluidState.setMassFraction(phase0Idx, comp1Idx, priVars[switchIdx]);
302 // TODO give values for non-existing nonwetting phase
303 Scalar xnH2O = FluidSystem::equilibriumMoleFraction(fluidState, paramCache, phase1Idx);
304 Scalar xnCO2 = 1 - xnH2O;
305 fluidState.setMoleFraction(phase1Idx, comp1Idx, xnCO2);
306 fluidState.setMoleFraction(phase1Idx, comp0Idx, xnH2O);
307 }
308 }
309
310 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
311 {
312 // set the viscosity and desity here if constraintsolver is not used
313 paramCache.updateComposition(fluidState, phaseIdx);
314 const Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
315 fluidState.setDensity(phaseIdx, rho);
316 const Scalar rhoMolar = FluidSystem::molarDensity(fluidState, phaseIdx);
317 fluidState.setMolarDensity(phaseIdx, rhoMolar);
318 const Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
319 fluidState.setViscosity(phaseIdx,mu);
320
321 // compute and set the enthalpy
322 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
323 fluidState.setEnthalpy(phaseIdx, h);
324 }
325 }
326
330 const FluidState &fluidState() const
331 { return fluidState_; }
332
336 const SolidState &solidState() const
337 { return solidState_; }
338
344 Scalar averageMolarMass(int phaseIdx) const
345 { return fluidState_.averageMolarMass(phaseIdx); }
346
353 Scalar saturation(const int phaseIdx) const
354 { return fluidState_.saturation(phaseIdx); }
355
363 Scalar massFraction(const int phaseIdx, const int compIdx) const
364 { return fluidState_.massFraction(phaseIdx, compIdx); }
365
373 Scalar moleFraction(const int phaseIdx, const int compIdx) const
374 { return fluidState_.moleFraction(phaseIdx, compIdx); }
375
382 Scalar density(const int phaseIdx) const
383 { return fluidState_.density(phaseIdx); }
384
391 Scalar viscosity(const int phaseIdx) const
392 { return fluidState_.viscosity(phaseIdx); }
393
400 Scalar molarDensity(const int phaseIdx) const
401 { return fluidState_.molarDensity(phaseIdx) ; }
402
409 Scalar pressure(const int phaseIdx) const
410 { return fluidState_.pressure(phaseIdx); }
411
419 Scalar temperature() const
420 { return fluidState_.temperature(/*phaseIdx=*/0); }
421
428 Scalar relativePermeability(const int phaseIdx) const
429 { return relativePermeability_[phaseIdx]; }
430
437 Scalar mobility(const int phaseIdx) const
438 { return relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx); }
439
444 Scalar capillaryPressure() const
445 { return fluidState_.pressure(phase1Idx) - fluidState_.pressure(phase0Idx); }
446
450 Scalar porosity() const
451 { return solidState_.porosity(); }
452
456 const PermeabilityType& permeability() const
457 { return permeability_; }
458
462 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
463 {
464 typename FluidSystem::ParameterCache paramCache;
465 paramCache.updatePhase(fluidState_, phaseIdx);
466 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
467 }
468
472 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
473 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
474
475
479 int wettingPhase() const
480 { return fluidState_.wettingPhase(); }
481
482private:
483 FluidState fluidState_;
484 SolidState solidState_;
485 Scalar pc_; // The capillary pressure
486 PermeabilityType permeability_; // Effective permeability within the control volume
487
488 // Relative permeability within the control volume
489 std::array<Scalar, ModelTraits::numFluidPhases()> relativePermeability_;
490
491 // Effective diffusion coefficients for the phases
492 DiffusionCoefficients effectiveDiffCoeff_;
493};
494
495} // end namespace Dumux
496
497#endif
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
The primary variable switch for the 2p2c-CO2 model controlling the phase presence state variable.
Definition: co2/primaryvariableswitch.hh:45
Contains the quantities which are are constant within a finite volume in the CO2 model.
Definition: porousmediumflow/co2/volumevariables.hh:51
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/co2/volumevariables.hh:363
const SolidState & solidState() const
Returns the phase state for the control volume.
Definition: porousmediumflow/co2/volumevariables.hh:336
typename Traits::FluidState FluidState
The type of the object returned by the fluidState() method.
Definition: porousmediumflow/co2/volumevariables.hh:95
Scalar saturation(const int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:353
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Completes the fluid state.
Definition: porousmediumflow/co2/volumevariables.hh:178
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:437
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/co2/volumevariables.hh:472
Scalar porosity() const
Returns the average porosity within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:450
typename ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/co2/volumevariables.hh:99
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:382
Scalar viscosity(const int phaseIdx) const
Returns the dynamic viscosity of the fluid within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:391
typename Traits::FluidSystem FluidSystem
The fluid system used here.
Definition: porousmediumflow/co2/volumevariables.hh:97
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/co2/volumevariables.hh:479
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/co2/volumevariables.hh:101
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/co2/volumevariables.hh:103
Scalar temperature() const
Returns temperature within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:419
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:409
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/co2/volumevariables.hh:109
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:456
const FluidState & fluidState() const
Returns the phase state within the control volume.
Definition: porousmediumflow/co2/volumevariables.hh:330
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/co2/volumevariables.hh:128
Scalar molarDensity(const int phaseIdx) const
Returns the mass density of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:400
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/co2/volumevariables.hh:344
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/co2/volumevariables.hh:373
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/co2/volumevariables.hh:111
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/co2/volumevariables.hh:462
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:444
Scalar relativePermeability(const int phaseIdx) const
Returns the relative permeability of a given phase within the control volume in .
Definition: porousmediumflow/co2/volumevariables.hh:428
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 extended Richards model.