version 3.10-dev
porousmediumflow/3p3c/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//
13#ifndef DUMUX_3P3C_VOLUME_VARIABLES_HH
14#define DUMUX_3P3C_VOLUME_VARIABLES_HH
15
20
25
27
28namespace Dumux {
29
30namespace Detail {
31// helper struct and function detecting if the fluid matrix interaction features a adsorptionModel() function
32template <class FluidMatrixInteraction>
33using AdsorptionModelDetector = decltype(std::declval<FluidMatrixInteraction>().adsorptionModel());
34
35template<class FluidMatrixInteraction>
36static constexpr bool hasAdsorptionModel()
37{ return Dune::Std::is_detected<AdsorptionModelDetector, FluidMatrixInteraction>::value; }
38
39}
40
46template <class Traits>
49, public EnergyVolumeVariables<Traits, ThreePThreeCVolumeVariables<Traits> >
50{
53
54 using Scalar = typename Traits::PrimaryVariables::value_type;
55 using PermeabilityType = typename Traits::PermeabilityType;
56
57 using FS = typename Traits::FluidSystem;
60
61 using ModelTraits = typename Traits::ModelTraits;
62 using Idx = typename ModelTraits::Indices;
63 static constexpr int numFluidComps = ParentType::numFluidComponents();
64 enum {
65 wCompIdx = FS::wCompIdx,
66 gCompIdx = FS::gCompIdx,
67 nCompIdx = FS::nCompIdx,
68
69 wPhaseIdx = FS::wPhaseIdx,
70 gPhaseIdx = FS::gPhaseIdx,
71 nPhaseIdx = FS::nPhaseIdx,
72
73 switch1Idx = Idx::switch1Idx,
74 switch2Idx = Idx::switch2Idx,
75 pressureIdx = Idx::pressureIdx
76 };
77
78 // present phases
79 enum {
80 threePhases = Idx::threePhases,
81 wPhaseOnly = Idx::wPhaseOnly,
82 gnPhaseOnly = Idx::gnPhaseOnly,
83 wnPhaseOnly = Idx::wnPhaseOnly,
84 gPhaseOnly = Idx::gPhaseOnly,
85 wgPhaseOnly = Idx::wgPhaseOnly
86 };
87
88 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
89 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
90
91public:
93 using FluidState = typename Traits::FluidState;
95 using FluidSystem = typename Traits::FluidSystem;
97 using Indices = typename ModelTraits::Indices;
99 using SolidState = typename Traits::SolidState;
101 using SolidSystem = typename Traits::SolidSystem;
104
114 template<class ElemSol, class Problem, class Element, class Scv>
115 void update(const ElemSol &elemSol,
116 const Problem &problem,
117 const Element &element,
118 const Scv& scv)
119 {
120 ParentType::update(elemSol, problem, element, scv);
121 const auto& priVars = elemSol[scv.localDofIndex()];
122 const auto phasePresence = priVars.state();
123
124 constexpr bool useConstraintSolver = ModelTraits::useConstraintSolver();
125
126 EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState_, solidState_);
127
128 /* first the saturations */
129 if (phasePresence == threePhases)
130 {
131 sw_ = priVars[switch1Idx];
132 sn_ = priVars[switch2Idx];
133 sg_ = 1. - sw_ - sn_;
134 }
135 else if (phasePresence == wPhaseOnly)
136 {
137 sw_ = 1.;
138 sn_ = 0.;
139 sg_ = 0.;
140 }
141 else if (phasePresence == gnPhaseOnly)
142 {
143 sw_ = 0.;
144 sn_ = priVars[switch2Idx];
145 sg_ = 1. - sn_;
146 }
147 else if (phasePresence == wnPhaseOnly)
148 {
149 sn_ = priVars[switch2Idx];
150 sw_ = 1. - sn_;
151 sg_ = 0.;
152 }
153 else if (phasePresence == gPhaseOnly)
154 {
155 sw_ = 0.;
156 sn_ = 0.;
157 sg_ = 1.;
158 }
159 else if (phasePresence == wgPhaseOnly)
160 {
161 sw_ = priVars[switch1Idx];
162 sn_ = 0.;
163 sg_ = 1. - sw_;
164 }
165 else
166 DUNE_THROW(Dune::InvalidStateException, "phasePresence: " << phasePresence << " is invalid.");
167
168 fluidState_.setSaturation(wPhaseIdx, sw_);
169 fluidState_.setSaturation(gPhaseIdx, sg_);
170 fluidState_.setSaturation(nPhaseIdx, sn_);
171
172 /* now the pressures */
173 pg_ = priVars[pressureIdx];
174
175 // calculate capillary pressures
176
177 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
178 Scalar pcgw = fluidMatrixInteraction.pcgw(sw_, sn_);
179 Scalar pcnw = fluidMatrixInteraction.pcnw(sw_, sn_);
180 Scalar pcgn = fluidMatrixInteraction.pcgn(sw_, sn_);
181
182 const Scalar pcAlpha = fluidMatrixInteraction.pcAlpha(sw_, sn_);
183 const Scalar pcNW1 = 0.0; // TODO: this should be possible to assign in the problem file
184
185 pn_ = pg_- pcAlpha * pcgn - (1.-pcAlpha)*(pcgw - pcNW1);
186 pw_ = pn_ - pcAlpha * pcnw - (1.-pcAlpha)*pcNW1;
187
188 fluidState_.setPressure(wPhaseIdx, pw_);
189 fluidState_.setPressure(gPhaseIdx, pg_);
190 fluidState_.setPressure(nPhaseIdx, pn_);
191
192 // calculate and set all fugacity coefficients. this is
193 // possible because we require all phases to be an ideal
194 // mixture, i.e. fugacity coefficients are not supposed to
195 // depend on composition!
196 typename FluidSystem::ParameterCache paramCache;
197 // assert(FluidSystem::isIdealGas(gPhaseIdx));
198 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++ phaseIdx) {
199 assert(FluidSystem::isIdealMixture(phaseIdx));
200
201 for (int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++ compIdx) {
202 Scalar phi = FluidSystem::fugacityCoefficient(fluidState_, paramCache, phaseIdx, compIdx);
203 fluidState_.setFugacityCoefficient(phaseIdx, compIdx, phi);
204 }
205 }
206
207 // now comes the tricky part: calculate phase composition
208 if (phasePresence == threePhases) {
209 // all phases are present, phase compositions are a
210 // result of the the gas <-> liquid equilibrium. This is
211 // the job of the "MiscibleMultiPhaseComposition"
212 // constraint solver ...
213 if (useConstraintSolver) {
215 paramCache);
216 }
217 // ... or calculated explicitly this way ...
218 // please note that we experienced some problems with un-regularized
219 // partial pressures due to their calculation from fugacity coefficients -
220 // that's why they are regularized below "within physically meaningful bounds"
221 else {
222 Scalar partPressH2O = FluidSystem::fugacityCoefficient(fluidState_,
223 wPhaseIdx,
224 wCompIdx) * pw_;
225 if (partPressH2O > pg_) partPressH2O = pg_;
226 Scalar partPressNAPL = FluidSystem::fugacityCoefficient(fluidState_,
227 nPhaseIdx,
228 nCompIdx) * pn_;
229 if (partPressNAPL > pg_) partPressNAPL = pg_;
230 Scalar partPressAir = pg_ - partPressH2O - partPressNAPL;
231
232 Scalar xgn = partPressNAPL/pg_;
233 Scalar xgw = partPressH2O/pg_;
234 Scalar xgg = partPressAir/pg_;
235
236 // actually, it's nothing else than Henry coefficient
237 Scalar xwn = partPressNAPL
238 / (FluidSystem::fugacityCoefficient(fluidState_,
239 wPhaseIdx,nCompIdx)
240 * pw_);
241 Scalar xwg = partPressAir
242 / (FluidSystem::fugacityCoefficient(fluidState_,
243 wPhaseIdx,gCompIdx)
244 * pw_);
245 Scalar xww = 1.-xwg-xwn;
246
247 Scalar xnn = 1.-2.e-10;
248 Scalar xna = 1.e-10;
249 Scalar xnw = 1.e-10;
250
251 fluidState_.setMoleFraction(wPhaseIdx, wCompIdx, xww);
252 fluidState_.setMoleFraction(wPhaseIdx, gCompIdx, xwg);
253 fluidState_.setMoleFraction(wPhaseIdx, nCompIdx, xwn);
254 fluidState_.setMoleFraction(gPhaseIdx, wCompIdx, xgw);
255 fluidState_.setMoleFraction(gPhaseIdx, gCompIdx, xgg);
256 fluidState_.setMoleFraction(gPhaseIdx, nCompIdx, xgn);
257 fluidState_.setMoleFraction(nPhaseIdx, wCompIdx, xnw);
258 fluidState_.setMoleFraction(nPhaseIdx, gCompIdx, xna);
259 fluidState_.setMoleFraction(nPhaseIdx, nCompIdx, xnn);
260
261 Scalar rhoW = FluidSystem::density(fluidState_, wPhaseIdx);
262 Scalar rhoG = FluidSystem::density(fluidState_, gPhaseIdx);
263 Scalar rhoN = FluidSystem::density(fluidState_, nPhaseIdx);
264 Scalar rhoWMolar = FluidSystem::molarDensity(fluidState_, wPhaseIdx);
265 Scalar rhoGMolar = FluidSystem::molarDensity(fluidState_, gPhaseIdx);
266 Scalar rhoNMolar = FluidSystem::molarDensity(fluidState_, nPhaseIdx);
267
268 fluidState_.setDensity(wPhaseIdx, rhoW);
269 fluidState_.setDensity(gPhaseIdx, rhoG);
270 fluidState_.setDensity(nPhaseIdx, rhoN);
271 fluidState_.setMolarDensity(wPhaseIdx, rhoWMolar);
272 fluidState_.setMolarDensity(gPhaseIdx, rhoGMolar);
273 fluidState_.setMolarDensity(nPhaseIdx, rhoNMolar);
274 }
275 }
276 else if (phasePresence == wPhaseOnly) {
277 // only the water phase is present, water phase composition is
278 // stored explicitly.
279
280 // extract mole fractions in the water phase
281 Scalar xwg = priVars[switch1Idx];
282 Scalar xwn = priVars[switch2Idx];
283 Scalar xww = 1 - xwg - xwn;
284
285 // write water mole fractions in the fluid state
286 fluidState_.setMoleFraction(wPhaseIdx, wCompIdx, xww);
287 fluidState_.setMoleFraction(wPhaseIdx, gCompIdx, xwg);
288 fluidState_.setMoleFraction(wPhaseIdx, nCompIdx, xwn);
289
290 // calculate the composition of the remaining phases (as
291 // well as the densities of all phases). this is the job
292 // of the "ComputeFromReferencePhase" constraint solver ...
293 if (useConstraintSolver)
294 {
296 paramCache,
297 wPhaseIdx);
298 }
299 // ... or calculated explicitly this way ...
300 else {
301 // note that the gas phase is actually not existing!
302 // thus, this is used as phase switch criterion
303 Scalar xgg = xwg * FluidSystem::fugacityCoefficient(fluidState_,
304 wPhaseIdx,gCompIdx)
305 * pw_ / pg_;
306 Scalar xgn = xwn * FluidSystem::fugacityCoefficient(fluidState_,
307 wPhaseIdx,nCompIdx)
308 * pw_ / pg_;
309 Scalar xgw = FluidSystem::fugacityCoefficient(fluidState_,
310 wPhaseIdx,wCompIdx)
311 * pw_ / pg_;
312
313
314 // note that the gas phase is actually not existing!
315 // thus, this is used as phase switch criterion
316 Scalar xnn = xwn * FluidSystem::fugacityCoefficient(fluidState_,
317 wPhaseIdx,nCompIdx)
318 * pw_;
319 Scalar xna = 1.e-10;
320 Scalar xnw = 1.e-10;
321
322 fluidState_.setMoleFraction(gPhaseIdx, wCompIdx, xgw);
323 fluidState_.setMoleFraction(gPhaseIdx, gCompIdx, xgg);
324 fluidState_.setMoleFraction(gPhaseIdx, nCompIdx, xgn);
325 fluidState_.setMoleFraction(nPhaseIdx, wCompIdx, xnw);
326 fluidState_.setMoleFraction(nPhaseIdx, gCompIdx, xna);
327 fluidState_.setMoleFraction(nPhaseIdx, nCompIdx, xnn);
328
329 Scalar rhoW = FluidSystem::density(fluidState_, wPhaseIdx);
330 Scalar rhoG = FluidSystem::density(fluidState_, gPhaseIdx);
331 Scalar rhoN = FluidSystem::density(fluidState_, nPhaseIdx);
332 Scalar rhoWMolar = FluidSystem::molarDensity(fluidState_, wPhaseIdx);
333 Scalar rhoGMolar = FluidSystem::molarDensity(fluidState_, gPhaseIdx);
334 Scalar rhoNMolar = FluidSystem::molarDensity(fluidState_, nPhaseIdx);
335
336 fluidState_.setDensity(wPhaseIdx, rhoW);
337 fluidState_.setDensity(gPhaseIdx, rhoG);
338 fluidState_.setDensity(nPhaseIdx, rhoN);
339 fluidState_.setMolarDensity(wPhaseIdx, rhoWMolar);
340 fluidState_.setMolarDensity(gPhaseIdx, rhoGMolar);
341 fluidState_.setMolarDensity(nPhaseIdx, rhoNMolar);
342 }
343 }
344 else if (phasePresence == gnPhaseOnly) {
345 // only gas and NAPL phases are present
346 // we have all (partly hypothetical) phase pressures
347 // and temperature and the mole fraction of water in
348 // the gas phase
349
350 // we have all (partly hypothetical) phase pressures
351 // and temperature and the mole fraction of water in
352 // the gas phase
353 Scalar partPressNAPL = fluidState_.fugacityCoefficient(nPhaseIdx, nCompIdx)*pn_;
354 if (partPressNAPL > pg_) partPressNAPL = pg_;
355
356 Scalar xgw = priVars[switch1Idx];
357 Scalar xgn = partPressNAPL/pg_;
358 Scalar xgg = 1.-xgw-xgn;
359
360 // write mole fractions in the fluid state
361 fluidState_.setMoleFraction(gPhaseIdx, wCompIdx, xgw);
362 fluidState_.setMoleFraction(gPhaseIdx, gCompIdx, xgg);
363 fluidState_.setMoleFraction(gPhaseIdx, nCompIdx, xgn);
364
365 // calculate the composition of the remaining phases (as
366 // well as the densities of all phases). this is the job
367 // of the "ComputeFromReferencePhase" constraint solver
369 paramCache,
370 gPhaseIdx);
371 }
372 else if (phasePresence == wnPhaseOnly) {
373 // only water and NAPL phases are present
374 Scalar partPressNAPL = fluidState_.fugacityCoefficient(nPhaseIdx,nCompIdx)*pn_;
375 if (partPressNAPL > pg_) partPressNAPL = pg_;
376 Scalar henryC = fluidState_.fugacityCoefficient(wPhaseIdx,nCompIdx)*pw_;
377
378 Scalar xwg = priVars[switch1Idx];
379 Scalar xwn = partPressNAPL/henryC;
380 Scalar xww = 1.-xwg-xwn;
381
382 // write mole fractions in the fluid state
383 fluidState_.setMoleFraction(wPhaseIdx, wCompIdx, xww);
384 fluidState_.setMoleFraction(wPhaseIdx, gCompIdx, xwg);
385 fluidState_.setMoleFraction(wPhaseIdx, nCompIdx, xwn);
386
387 // calculate the composition of the remaining phases (as
388 // well as the densities of all phases). this is the job
389 // of the "ComputeFromReferencePhase" constraint solver
391 paramCache,
392 wPhaseIdx);
393 }
394 else if (phasePresence == gPhaseOnly) {
395 // only the gas phase is present, gas phase composition is
396 // stored explicitly here below.
397
398 const Scalar xgw = priVars[switch1Idx];
399 const Scalar xgn = priVars[switch2Idx];
400 Scalar xgg = 1 - xgw - xgn;
401
402 // write mole fractions in the fluid state
403 fluidState_.setMoleFraction(gPhaseIdx, wCompIdx, xgw);
404 fluidState_.setMoleFraction(gPhaseIdx, gCompIdx, xgg);
405 fluidState_.setMoleFraction(gPhaseIdx, nCompIdx, xgn);
406
407 // calculate the composition of the remaining phases (as
408 // well as the densities of all phases). this is the job
409 // of the "ComputeFromReferencePhase" constraint solver ...
410 if (useConstraintSolver)
411 {
413 paramCache,
414 gPhaseIdx);
415 }
416 // ... or calculated explicitly this way ...
417 else {
418
419 // note that the water phase is actually not existing!
420 // thus, this is used as phase switch criterion
421 Scalar xww = xgw * pg_
422 / (FluidSystem::fugacityCoefficient(fluidState_,
423 wPhaseIdx,wCompIdx)
424 * pw_);
425 Scalar xwn = 1.e-10;
426 Scalar xwg = 1.e-10;
427
428 // note that the NAPL phase is actually not existing!
429 // thus, this is used as phase switch criterion
430 Scalar xnn = xgn * pg_
431 / (FluidSystem::fugacityCoefficient(fluidState_,
432 nPhaseIdx,nCompIdx)
433 * pn_);
434 Scalar xna = 1.e-10;
435 Scalar xnw = 1.e-10;
436
437 fluidState_.setMoleFraction(wPhaseIdx, wCompIdx, xww);
438 fluidState_.setMoleFraction(wPhaseIdx, gCompIdx, xwg);
439 fluidState_.setMoleFraction(wPhaseIdx, nCompIdx, xwn);
440 fluidState_.setMoleFraction(nPhaseIdx, wCompIdx, xnw);
441 fluidState_.setMoleFraction(nPhaseIdx, gCompIdx, xna);
442 fluidState_.setMoleFraction(nPhaseIdx, nCompIdx, xnn);
443
444 Scalar rhoW = FluidSystem::density(fluidState_, wPhaseIdx);
445 Scalar rhoG = FluidSystem::density(fluidState_, gPhaseIdx);
446 Scalar rhoN = FluidSystem::density(fluidState_, nPhaseIdx);
447 Scalar rhoWMolar = FluidSystem::molarDensity(fluidState_, wPhaseIdx);
448 Scalar rhoGMolar = FluidSystem::molarDensity(fluidState_, gPhaseIdx);
449 Scalar rhoNMolar = FluidSystem::molarDensity(fluidState_, nPhaseIdx);
450
451 fluidState_.setDensity(wPhaseIdx, rhoW);
452 fluidState_.setDensity(gPhaseIdx, rhoG);
453 fluidState_.setDensity(nPhaseIdx, rhoN);
454 fluidState_.setMolarDensity(wPhaseIdx, rhoWMolar);
455 fluidState_.setMolarDensity(gPhaseIdx, rhoGMolar);
456 fluidState_.setMolarDensity(nPhaseIdx, rhoNMolar);
457 }
458 }
459 else if (phasePresence == wgPhaseOnly) {
460 // only water and gas phases are present
461 Scalar xgn = priVars[switch2Idx];
462 Scalar partPressH2O = fluidState_.fugacityCoefficient(wPhaseIdx, wCompIdx)*pw_;
463 if (partPressH2O > pg_) partPressH2O = pg_;
464
465 Scalar xgw = partPressH2O/pg_;
466 Scalar xgg = 1.-xgn-xgw;
467
468 // write mole fractions in the fluid state
469 fluidState_.setMoleFraction(gPhaseIdx, wCompIdx, xgw);
470 fluidState_.setMoleFraction(gPhaseIdx, gCompIdx, xgg);
471 fluidState_.setMoleFraction(gPhaseIdx, nCompIdx, xgn);
472
473 // calculate the composition of the remaining phases (as
474 // well as the densities of all phases). this is the job
475 // of the "ComputeFromReferencePhase" constraint solver ...
476 if (useConstraintSolver)
477 {
479 paramCache,
480 gPhaseIdx);
481 }
482 // ... or calculated explicitly this way ...
483 else {
484 // actually, it's nothing else than Henry coefficient
485 Scalar xwn = xgn * pg_
486 / (FluidSystem::fugacityCoefficient(fluidState_,
487 wPhaseIdx,nCompIdx)
488 * pw_);
489 Scalar xwg = xgg * pg_
490 / (FluidSystem::fugacityCoefficient(fluidState_,
491 wPhaseIdx,gCompIdx)
492 * pw_);
493 Scalar xww = 1.-xwg-xwn;
494
495 // note that the NAPL phase is actually not existing!
496 // thus, this is used as phase switch criterion
497 Scalar xnn = xgn * pg_
498 / (FluidSystem::fugacityCoefficient(fluidState_,
499 nPhaseIdx,nCompIdx)
500 * pn_);
501 Scalar xna = 1.e-10;
502 Scalar xnw = 1.e-10;
503
504 fluidState_.setMoleFraction(wPhaseIdx, wCompIdx, xww);
505 fluidState_.setMoleFraction(wPhaseIdx, gCompIdx, xwg);
506 fluidState_.setMoleFraction(wPhaseIdx, nCompIdx, xwn);
507 fluidState_.setMoleFraction(nPhaseIdx, wCompIdx, xnw);
508 fluidState_.setMoleFraction(nPhaseIdx, gCompIdx, xna);
509 fluidState_.setMoleFraction(nPhaseIdx, nCompIdx, xnn);
510
511 Scalar rhoW = FluidSystem::density(fluidState_, wPhaseIdx);
512 Scalar rhoG = FluidSystem::density(fluidState_, gPhaseIdx);
513 Scalar rhoN = FluidSystem::density(fluidState_, nPhaseIdx);
514 Scalar rhoWMolar = FluidSystem::molarDensity(fluidState_, wPhaseIdx);
515 Scalar rhoGMolar = FluidSystem::molarDensity(fluidState_, gPhaseIdx);
516 Scalar rhoNMolar = FluidSystem::molarDensity(fluidState_, nPhaseIdx);
517
518 fluidState_.setDensity(wPhaseIdx, rhoW);
519 fluidState_.setDensity(gPhaseIdx, rhoG);
520 fluidState_.setDensity(nPhaseIdx, rhoN);
521 fluidState_.setMolarDensity(wPhaseIdx, rhoWMolar);
522 fluidState_.setMolarDensity(gPhaseIdx, rhoGMolar);
523 fluidState_.setMolarDensity(nPhaseIdx, rhoNMolar);
524 }
525 }
526 else
527 DUNE_THROW(Dune::InvalidStateException, "phasePresence: " << phasePresence << " is invalid.");
528
529 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
530 {
531 // mobilities
532 const Scalar mu =
534 paramCache,
535 phaseIdx);
536 fluidState_.setViscosity(phaseIdx,mu);
537
538 const Scalar kr = fluidMatrixInteraction.kr(phaseIdx,
539 fluidState_.saturation(wPhaseIdx),
540 fluidState_.saturation(nPhaseIdx));
541 mobility_[phaseIdx] = kr / mu;
542 }
543
544 // material dependent parameters for NAPL adsorption (only if law is provided)
545 if constexpr (Detail::hasAdsorptionModel<std::decay_t<decltype(fluidMatrixInteraction)>>())
546 bulkDensTimesAdsorpCoeff_ = fluidMatrixInteraction.adsorptionModel().bulkDensTimesAdsorpCoeff();
547
548 /* compute the diffusion coefficient
549 * \note This is the part of the diffusion coefficient determined by the fluid state, e.g.
550 * important if they are tabularized. In the diffusive flux computation (e.g. Fick's law)
551 * this gets converted into an efficient coefficient depending on saturation and porosity.
552 * We can then add a normalized tensorial component
553 * e.g. obtained from DTI from the spatial params (currently not implemented)
554 */
555
556 auto getEffectiveDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
557 {
558 return EffDiffModel::effectiveDiffusionCoefficient(*this, phaseIdx, compIIdx, compJIdx);
559 };
560
561 // porosity & permeabilty
562 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
563
564 effectiveDiffCoeff_.update(getEffectiveDiffusionCoefficient);
565
566 EnergyVolVars::updateSolidEnergyParams(elemSol, problem, element, scv, solidState_);
567 permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
568 EnergyVolVars::updateEffectiveThermalConductivity();
569
570 // compute and set the enthalpy
571 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
572 {
573 Scalar h = EnergyVolVars::enthalpy(fluidState_, paramCache, phaseIdx);
574 fluidState_.setEnthalpy(phaseIdx, h);
575 }
576 }
577
581 const FluidState &fluidState() const
582 { return fluidState_; }
583
587 const SolidState &solidState() const
588 { return solidState_; }
589
595 Scalar averageMolarMass(int phaseIdx) const
596 { return fluidState_.averageMolarMass(phaseIdx); }
597
604 Scalar saturation(const int phaseIdx) const
605 { return fluidState_.saturation(phaseIdx); }
606
614 Scalar massFraction(const int phaseIdx, const int compIdx) const
615 { return fluidState_.massFraction(phaseIdx, compIdx); }
616
624 Scalar moleFraction(const int phaseIdx, const int compIdx) const
625 { return fluidState_.moleFraction(phaseIdx, compIdx); }
626
633 Scalar density(const int phaseIdx) const
634 { return fluidState_.density(phaseIdx); }
635
642 Scalar molarDensity(const int phaseIdx) const
643 { return fluidState_.molarDensity(phaseIdx); }
644
651 Scalar pressure(const int phaseIdx) const
652 { return fluidState_.pressure(phaseIdx); }
653
661 Scalar temperature() const
662 { return fluidState_.temperature(/*phaseIdx=*/0); }
663
670 Scalar mobility(const int phaseIdx) const
671 { return mobility_[phaseIdx]; }
672
676 Scalar capillaryPressure() const
677 { return fluidState_.capillaryPressure(); }
678
682 Scalar porosity() const
683 { return solidState_.porosity(); }
684
689 {
690 if (bulkDensTimesAdsorpCoeff_)
691 return bulkDensTimesAdsorpCoeff_.value();
692 else
693 DUNE_THROW(Dune::NotImplemented, "Your spatialParams do not provide an adsorption model");
694 }
695
699 const PermeabilityType& permeability() const
700 { return permeability_; }
701
702 /*
703 * \brief Returns the binary diffusion coefficients for a phase in \f$[m^2/s]\f$.
704 */
705 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
706 {
707 typename FluidSystem::ParameterCache paramCache;
708 paramCache.updatePhase(fluidState_, phaseIdx);
709 return FluidSystem::diffusionCoefficient(fluidState_, paramCache, phaseIdx, compJIdx);
710 }
711
715 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
716 { return effectiveDiffCoeff_(phaseIdx, compIIdx, compJIdx); }
717
718protected:
721
722
723private:
724 Scalar sw_, sg_, sn_, pg_, pw_, pn_;
725
726 Scalar moleFrac_[ModelTraits::numFluidPhases()][ModelTraits::numFluidComponents()];
727 Scalar massFrac_[ModelTraits::numFluidPhases()][ModelTraits::numFluidComponents()];
728
729 PermeabilityType permeability_;
730 Scalar mobility_[ModelTraits::numFluidPhases()];
731 OptionalScalar<Scalar> bulkDensTimesAdsorpCoeff_;
732
733 DiffusionCoefficients effectiveDiffCoeff_;
734};
735
736} // end namespace Dumux
737
738#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: 3p3c/primaryvariableswitch.hh:28
Contains the quantities which are are constant within a finite volume in the three-phase three-compon...
Definition: porousmediumflow/3p3c/volumevariables.hh:50
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/3p3c/volumevariables.hh:614
typename ModelTraits::Indices Indices
export the indices
Definition: porousmediumflow/3p3c/volumevariables.hh:97
Scalar averageMolarMass(int phaseIdx) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/3p3c/volumevariables.hh:595
FluidState fluidState_
Definition: porousmediumflow/3p3c/volumevariables.hh:719
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Update all quantities for a given control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:115
SolidState solidState_
Definition: porousmediumflow/3p3c/volumevariables.hh:720
Scalar saturation(const int phaseIdx) const
Returns the effective saturation of a given phase within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:604
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/3p3c/volumevariables.hh:624
Scalar molarDensity(const int phaseIdx) const
Returns the molar density of a given phase within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:642
typename Traits::SolidState SolidState
export type of solid state
Definition: porousmediumflow/3p3c/volumevariables.hh:99
Scalar pressure(const int phaseIdx) const
Returns the effective pressure of a given phase within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:651
const SolidState & solidState() const
Returns the phase state for the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:587
const PermeabilityType & permeability() const
Returns the average permeability within the control volume in .
Definition: porousmediumflow/3p3c/volumevariables.hh:699
typename Traits::FluidState FluidState
export fluid state type
Definition: porousmediumflow/3p3c/volumevariables.hh:93
Scalar temperature() const
Returns temperature inside the sub-control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:661
Scalar capillaryPressure() const
Returns the effective capillary pressure within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:676
Scalar mobility(const int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:670
const FluidState & fluidState() const
Returns the phase state for the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:581
Scalar density(const int phaseIdx) const
Returns the mass density of a given phase within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:633
typename Traits::FluidSystem FluidSystem
export fluid system type
Definition: porousmediumflow/3p3c/volumevariables.hh:95
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/3p3c/volumevariables.hh:715
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Definition: porousmediumflow/3p3c/volumevariables.hh:705
Scalar bulkDensTimesAdsorpCoeff() const
Returns the adsorption information.
Definition: porousmediumflow/3p3c/volumevariables.hh:688
Scalar porosity() const
Returns the average porosity within the control volume.
Definition: porousmediumflow/3p3c/volumevariables.hh:682
typename Traits::SolidSystem SolidSystem
export type of solid system
Definition: porousmediumflow/3p3c/volumevariables.hh:101
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.
A central place for various physical constants occurring in some equations.
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
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
decltype(std::declval< FluidMatrixInteraction >().adsorptionModel()) AdsorptionModelDetector
Definition: porousmediumflow/3p3c/volumevariables.hh:33
static constexpr bool hasAdsorptionModel()
Definition: porousmediumflow/3p3c/volumevariables.hh:36
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
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
A wrapper that can either contain a valid Scalar or NaN.
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.
T value() const
Definition: optionalscalar.hh:36
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.