version 3.10-dev
porousmediumflow/2pnc/volumevariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_2PNC_VOLUME_VARIABLES_HH
15#define DUMUX_2PNC_VOLUME_VARIABLES_HH
16
17#include <iostream>
18#include <vector>
19
20#include <dumux/common/math.hh>
22
25
30
32
34
35namespace Dumux {
36
42template <class Traits>
45, public EnergyVolumeVariables<Traits, TwoPNCVolumeVariables<Traits> >
46{
49 using Scalar = typename Traits::PrimaryVariables::value_type;
50 using PermeabilityType = typename Traits::PermeabilityType;
51 using FS = typename Traits::FluidSystem;
52 using ModelTraits = typename Traits::ModelTraits;
53 static constexpr int numFluidComps = ParentType::numFluidComponents();
54 enum
55 {
56 numMajorComponents = ModelTraits::numFluidPhases(),
57
58 // phase indices
59 phase0Idx = FS::phase0Idx,
60 phase1Idx = FS::phase1Idx,
61
62 // component indices
63 comp0Idx = FS::comp0Idx,
64 comp1Idx = FS::comp1Idx,
65
66 // phase presence enums
67 secondPhaseOnly = ModelTraits::Indices::secondPhaseOnly,
68 firstPhaseOnly = ModelTraits::Indices::firstPhaseOnly,
69 bothPhases = ModelTraits::Indices::bothPhases,
70
71 // primary variable indices
72 pressureIdx = ModelTraits::Indices::pressureIdx,
73 switchIdx = ModelTraits::Indices::switchIdx
74 };
75
76 static constexpr auto formulation = ModelTraits::priVarFormulation();
77 static constexpr bool setFirstPhaseMoleFractions = ModelTraits::setMoleFractionsForFirstPhase();
78
81 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
82 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
83
84public:
86 using FluidState = typename Traits::FluidState;
88 using FluidSystem = typename Traits::FluidSystem;
90 using Indices = typename ModelTraits::Indices;
92 using SolidState = typename Traits::SolidState;
94 using SolidSystem = typename Traits::SolidSystem;
97
99 static constexpr bool useMoles() { return Traits::ModelTraits::useMoles(); }
101 static constexpr TwoPFormulation priVarFormulation() { return formulation; }
102
103 // check for permissive specifications
104 static_assert(useMoles(), "use moles has to be set true in the 2pnc model");
105 static_assert(ModelTraits::numFluidPhases() == 2, "NumPhases set in the model is not two!");
106 static_assert((formulation == TwoPFormulation::p0s1 || formulation == TwoPFormulation::p1s0), "Chosen TwoPFormulation not supported!");
107
117 template<class ElemSol, class Problem, class Element, class Scv>
118 void update(const ElemSol &elemSol,
119 const Problem &problem,
120 const Element &element,
121 const Scv& scv)
122 {
123 ParentType::update(elemSol, problem, element, scv);
124
125 completeFluidState(elemSol, problem, element, scv, fluidState_, solidState_);
126
128 // calculate the remaining quantities
130
131 const auto& spatialParams = problem.spatialParams();
132 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
133
134 const int wPhaseIdx = fluidState_.wettingPhase();
135 const int nPhaseIdx = 1 - wPhaseIdx;
136
137 // mobilities -> require wetting phase saturation as parameter!
138 mobility_[wPhaseIdx] = fluidMatrixInteraction.krw(saturation(wPhaseIdx))/fluidState_.viscosity(wPhaseIdx);
139 mobility_[nPhaseIdx] = fluidMatrixInteraction.krn(saturation(wPhaseIdx))/fluidState_.viscosity(nPhaseIdx);
140
141 //update porosity before calculating the effective properties depending on it
142 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
143
144 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
145 {
146 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
147 };
148
149 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
150
151 // calculate the remaining quantities
152 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
153 permeability_ = spatialParams.permeability(element, scv, elemSol);
154 EnergyVolVars::updateEffectiveThermalConductivity();
155 }
156
169 template<class ElemSol, class Problem, class Element, class Scv>
170 void completeFluidState(const ElemSol& elemSol,
171 const Problem& problem,
172 const Element& element,
173 const Scv& scv,
176 {
177 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
178
179 const auto& priVars = elemSol[scv.localDofIndex()];
180 const auto phasePresence = priVars.state();
181
182 const auto& spatialParams = problem.spatialParams();
183 const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol);
184 const auto wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, scv, elemSol);
185 fluidState.setWettingPhase(wPhaseIdx);
186
187 // set the saturations
188 if (phasePresence == secondPhaseOnly)
189 {
190 fluidState.setSaturation(phase0Idx, 0.0);
191 fluidState.setSaturation(phase1Idx, 1.0);
192 }
193 else if (phasePresence == firstPhaseOnly)
194 {
195 fluidState.setSaturation(phase0Idx, 1.0);
196 fluidState.setSaturation(phase1Idx, 0.0);
197 }
198 else if (phasePresence == bothPhases)
199 {
200 if (formulation == TwoPFormulation::p0s1)
201 {
202 fluidState.setSaturation(phase1Idx, priVars[switchIdx]);
203 fluidState.setSaturation(phase0Idx, 1 - priVars[switchIdx]);
204 }
205 else
206 {
207 fluidState.setSaturation(phase0Idx, priVars[switchIdx]);
208 fluidState.setSaturation(phase1Idx, 1 - priVars[switchIdx]);
209 }
210 }
211 else
212 DUNE_THROW(Dune::InvalidStateException, "phasePresence: " << phasePresence << " is invalid.");
213
214 // set pressures of the fluid phases
215 pc_ = fluidMatrixInteraction.pc(fluidState.saturation(wPhaseIdx));
216 if (formulation == TwoPFormulation::p0s1)
217 {
218 fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
219 fluidState.setPressure(phase1Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] + pc_
220 : priVars[pressureIdx] - pc_);
221 }
222 else
223 {
224 fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
225 fluidState.setPressure(phase0Idx, (wPhaseIdx == phase0Idx) ? priVars[pressureIdx] - pc_
226 : priVars[pressureIdx] + pc_);
227 }
228
229 // calculate the phase compositions
230 typename FluidSystem::ParameterCache paramCache;
231
232 // now comes the tricky part: calculate phase composition
233 if (phasePresence == bothPhases)
234 {
235 // both phases are present, phase composition results from
236 // the first <-> second phase equilibrium. This is the job
237 // of the "MiscibleMultiPhaseComposition" constraint solver
238
239 // set the known mole fractions in the fluidState so that they
240 // can be used by the MiscibleMultiPhaseComposition constraint solver
241
242 const int knownPhaseIdx = setFirstPhaseMoleFractions ? phase0Idx : phase1Idx;
243 for (int compIdx = numMajorComponents; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
244 fluidState.setMoleFraction(knownPhaseIdx, compIdx, priVars[compIdx]);
245
247 paramCache,
248 knownPhaseIdx);
249 }
250 else if (phasePresence == secondPhaseOnly)
251 {
252 Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()> moleFrac;
253
254 moleFrac[comp0Idx] = priVars[switchIdx];
255 Scalar sumMoleFracOtherComponents = moleFrac[comp0Idx];
256
257 for (int compIdx = numMajorComponents; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
258 {
259 moleFrac[compIdx] = priVars[compIdx];
260 sumMoleFracOtherComponents += moleFrac[compIdx];
261 }
262
263 moleFrac[comp1Idx] = 1 - sumMoleFracOtherComponents;
264
265 // Set fluid state mole fractions
266 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
267 fluidState.setMoleFraction(phase1Idx, compIdx, moleFrac[compIdx]);
268
269 // calculate the composition of the remaining phases (as
270 // well as the densities of all phases). this is the job
271 // of the "ComputeFromReferencePhase" constraint solver
273 paramCache,
274 phase1Idx);
275 }
276 else if (phasePresence == firstPhaseOnly)
277 {
278 // only the first phase is present, i.e. first phase composition
279 // is stored explicitly. extract _mass_ fractions in the second phase
280 Dune::FieldVector<Scalar, ModelTraits::numFluidComponents()> moleFrac;
281
282 moleFrac[comp1Idx] = priVars[switchIdx];
283 Scalar sumMoleFracOtherComponents = moleFrac[comp1Idx];
284 for (int compIdx = numMajorComponents; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
285 {
286 moleFrac[compIdx] = priVars[compIdx];
287
288 sumMoleFracOtherComponents += moleFrac[compIdx];
289 }
290
291 moleFrac[comp0Idx] = 1 - sumMoleFracOtherComponents;
292
293 // convert mass to mole fractions and set the fluid state
294 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
295 fluidState.setMoleFraction(phase0Idx, compIdx, moleFrac[compIdx]);
296
297 // calculate the composition of the remaining phases (as
298 // well as the densities of all phases). this is the job
299 // of the "ComputeFromReferencePhase" constraint solver
301 paramCache,
302 phase0Idx);
303 }
304 paramCache.updateAll(fluidState);
305 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
306 {
307 Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
308 Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, phaseIdx);
309 fluidState.setViscosity(phaseIdx, mu);
310 fluidState.setEnthalpy(phaseIdx, h);
311 }
312 }
313
317 const FluidState &fluidState() const
318 { return fluidState_; }
319
323 const SolidState &solidState() const
324 { return solidState_; }
325
331 Scalar averageMolarMass(int phaseIdx) const
332 { return fluidState_.averageMolarMass(phaseIdx); }
333
340 Scalar saturation(int phaseIdx) const
341 { return fluidState_.saturation(phaseIdx); }
342
349 Scalar density(int phaseIdx) const
350 { return fluidState_.density(phaseIdx); }
351
358 Scalar viscosity(int phaseIdx) const
359 { return fluidState_.viscosity(phaseIdx); }
360
367 Scalar molarDensity(int phaseIdx) const
368 {
369 if (phaseIdx < ModelTraits::numFluidPhases())
370 return fluidState_.molarDensity(phaseIdx);
371
372 else
373 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
374 }
375
382 Scalar pressure(int phaseIdx) const
383 { return fluidState_.pressure(phaseIdx); }
384
392 Scalar temperature() const
393 { return fluidState_.temperature(/*phaseIdx=*/0); }
394
401 Scalar mobility(int phaseIdx) const
402 { return mobility_[phaseIdx]; }
403
408 Scalar capillaryPressure() const
409 { return pc_; }
410
414 Scalar porosity() const
415 { return solidState_.porosity(); }
416
420 const PermeabilityType& permeability() const
421 { return permeability_; }
422
426 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
427 {
428 typename FluidSystem::ParameterCache paramCache;
429 paramCache.updatePhase(fluidState_, phaseIdx);
430 return FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, compIIdx, compJIdx);
431 }
432
436 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
437 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
438
445 Scalar massFraction(int phaseIdx, int compIdx) const
446 { return fluidState_.massFraction(phaseIdx, compIdx); }
447
454 Scalar moleFraction(int phaseIdx, int compIdx) const
455 { return fluidState_.moleFraction(phaseIdx, compIdx); }
456
460 int wettingPhase() const
461 { return fluidState_.wettingPhase(); }
462
463protected:
466
467private:
468
469 Scalar pc_; // The capillary pressure
470 Scalar porosity_; // Effective porosity within the control volume
471 PermeabilityType permeability_; // Effective permeability within the control volume
472 Scalar mobility_[ModelTraits::numFluidPhases()]; // Effective mobility within the control volume
473 DiffusionCoefficients effectiveDiffCoeff_;
474
475};
476
477} // end namespace Dumux
478
479#endif
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:64
static void solve(FluidState &fluidState, ParameterCache &paramCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:99
Definition: porousmediumflow/nonisothermal/volumevariables.hh:63
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:47
static void solve(FluidState &fluidState, ParameterCache &paramCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:69
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:28
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition: porousmediumflow/volumevariables.hh:40
const PrimaryVariables & priVars() const
Returns the vector of primary variables.
Definition: porousmediumflow/volumevariables.hh:64
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/volumevariables.hh:52
The primary variable switch controlling the phase presence state variable.
Definition: 2pnc/primaryvariableswitch.hh:29
Contains the quantities which are are constant within a finite volume in the two-phase,...
Definition: porousmediumflow/2pnc/volumevariables.hh:46
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:414
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/2pnc/volumevariables.hh:426
Scalar molarDensity(int phaseIdx) const
Returns the molar density of a given phase within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:367
static constexpr TwoPFormulation priVarFormulation()
Return the two-phase formulation used here.
Definition: porousmediumflow/2pnc/volumevariables.hh:101
Scalar viscosity(int phaseIdx) const
Returns the dynamic viscosity of a given phase within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:358
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume in .
Definition: porousmediumflow/2pnc/volumevariables.hh:408
Scalar density(int phaseIdx) const
Returns the mass density of a given phase within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:349
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:118
const SolidState & solidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:323
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the mole fraction of a component in the phase.
Definition: porousmediumflow/2pnc/volumevariables.hh:454
typename ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/2pnc/volumevariables.hh:90
typename Traits::SolidState SolidState
Export type of solid state.
Definition: porousmediumflow/2pnc/volumevariables.hh:92
Scalar saturation(int phaseIdx) const
Returns the saturation of a given phase within the control volume in .
Definition: porousmediumflow/2pnc/volumevariables.hh:340
typename Traits::SolidSystem SolidSystem
Export type of solid system.
Definition: porousmediumflow/2pnc/volumevariables.hh:94
const PermeabilityType & permeability() const
Returns the permeability within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:420
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/2pnc/volumevariables.hh:436
Scalar temperature() const
Returns temperature inside the sub-control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:392
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/2pnc/volumevariables.hh:331
SolidState solidState_
Definition: porousmediumflow/2pnc/volumevariables.hh:465
FluidState fluidState_
Definition: porousmediumflow/2pnc/volumevariables.hh:464
typename Traits::FluidSystem FluidSystem
Export fluid system type.
Definition: porousmediumflow/2pnc/volumevariables.hh:88
typename Traits::FluidState FluidState
Export fluid state type.
Definition: porousmediumflow/2pnc/volumevariables.hh:86
void completeFluidState(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, FluidState &fluidState, SolidState &solidState)
Sets complete fluid state.
Definition: porousmediumflow/2pnc/volumevariables.hh:170
int wettingPhase() const
Returns the wetting phase index.
Definition: porousmediumflow/2pnc/volumevariables.hh:460
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of a component in the phase.
Definition: porousmediumflow/2pnc/volumevariables.hh:445
Scalar pressure(int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:382
static constexpr bool useMoles()
Return whether moles or masses are balanced.
Definition: porousmediumflow/2pnc/volumevariables.hh:99
Scalar mobility(int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:401
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: porousmediumflow/2pnc/volumevariables.hh:317
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Computes all quantities of a generic fluid state if a reference phase has been specified.
Defines an enumeration for the formulations accepted by the two-phase model.
void updateSolidVolumeFractions(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, SolidState &solidState, const int solidVolFracOffset)
update the solid volume fractions (inert and reacitve) and set them in the solidstate
Definition: updatesolidvolumefractions.hh:24
TwoPFormulation
Enumerates the formulations which the two-phase model accepts.
Definition: formulation.hh:23
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
Define some often used mathematical functions.
The available discretization methods in Dumux.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
std::string phasePresence() noexcept
I/O name of phase presence.
Definition: name.hh:135
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:62
Definition: adapt.hh:17
Base class for the model specific class which provides access to all volume averaged quantities.
Base class for the model specific class which provides access to all volume averaged quantities.
The primary variable switch for the extended Richards model.
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.