3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/nonequilibrium/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 *****************************************************************************/
30#ifndef DUMUX_NONEQUILIBRIUM_VOLUME_VARIABLES_HH
31#define DUMUX_NONEQUILIBRIUM_VOLUME_VARIABLES_HH
32
33#include <cassert>
34#include <array>
35
38
39namespace Dumux {
40
47template<class Traits, class EquilibriumVolumeVariables, bool enableChemicalNonEquilibrium,
48 bool enableThermalNonEquilibrium, int numEnergyEqFluid>
50
51template<class Traits, class EquilibriumVolumeVariables>
54 EquilibriumVolumeVariables,
55 Traits::ModelTraits::enableChemicalNonEquilibrium(),
56 Traits::ModelTraits::enableThermalNonEquilibrium(),
57 Traits::ModelTraits::numEnergyEqFluid()>;
58
60// specialization for the case of NO kinetic mass but kinetic energy transfer of two fluids and a solid
62template<class Traits, class EquilibriumVolumeVariables>
64 EquilibriumVolumeVariables,
65 false/*chemicalNonEquilibrium?*/,
66 true/*thermalNonEquilibrium?*/,
67 2>
68: public EquilibriumVolumeVariables
69{
70 using ParentType = EquilibriumVolumeVariables;
71 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
72 using Scalar = typename Traits::PrimaryVariables::value_type;
73
74 using ModelTraits = typename Traits::ModelTraits;
75
76 using FS = typename Traits::FluidSystem;
77 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
78 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
79
80 static constexpr auto phase0Idx = FS::phase0Idx;
81 static constexpr auto phase1Idx = FS::phase1Idx;
82 static constexpr auto sPhaseIdx = FS::numPhases;
83
85
86 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
87 using InterfacialAreasArray = std::array<std::array<Scalar, ModelTraits::numFluidPhases()+numEnergyEqSolid>,
88 ModelTraits::numFluidPhases()+numEnergyEqSolid>;
89
90public:
91 using FluidState = typename Traits::FluidState;
92 using Indices = typename ModelTraits::Indices;
101 template<class ElemSol, class Problem, class Element, class Scv>
102 void update(const ElemSol &elemSol,
103 const Problem &problem,
104 const Element &element,
105 const Scv& scv)
106 {
107 // Update parent type (also completes the fluid state)
108 ParentType::update(elemSol, problem, element, scv);
109
110 ParameterCache paramCache;
111 paramCache.updateAll(this->fluidState());
112 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
113 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
114 }
115
126 template<class ElemSol, class Problem, class Element, class Scv>
127 void updateDimLessNumbers(const ElemSol& elemSol,
128 const FluidState& fluidState,
129 const ParameterCache& paramCache,
130 const Problem& problem,
131 const Element& element,
132 const Scv& scv)
133 {
134 factorEnergyTransfer_ = problem.spatialParams().factorEnergyTransfer(element, scv, elemSol);
135 characteristicLength_ = problem.spatialParams().characteristicLength(element, scv, elemSol);
136
137 // set the dimensionless numbers and obtain respective quantities
138 const unsigned int vIdxGlobal = scv.dofIndex();
139 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
140 {
141 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
142 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
143 const auto density = fluidState.density(phaseIdx);
144 const auto kinematicViscosity = dynamicViscosity/density;
145
146 using FluidSystem = typename Traits::FluidSystem;
147 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
148 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
149 const auto porosity = this->porosity();
150
151 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_,kinematicViscosity);
152 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
153 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
154 prandtlNumber_[phaseIdx],
155 porosity,
156 ModelTraits::nusseltFormulation());
157
158 }
159 }
160
171 template<class ElemSol, class Problem, class Element, class Scv>
172 void updateInterfacialArea(const ElemSol& elemSol,
173 const FluidState& fluidState,
174 const ParameterCache& paramCache,
175 const Problem& problem,
176 const Element& element,
177 const Scv& scv)
178 {
179 // obtain (standard) material parameters (needed for the residual saturations)
180 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
181
182 //obtain parameters for interfacial area constitutive relations
183 const auto& aWettingNonWettingSurfaceParams = problem.spatialParams().aWettingNonWettingSurfaceParams(element, scv, elemSol);
184
185 const Scalar pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
186 const Scalar Sw = fluidState.saturation(phase0Idx);
187
188 using AwnSurface = typename Problem::SpatialParams::AwnSurface;
189 const auto awn = AwnSurface::interfacialArea(aWettingNonWettingSurfaceParams, materialParams, Sw, pc);
190 interfacialArea_[phase0Idx][phase1Idx] = awn;
191 interfacialArea_[phase1Idx][phase0Idx] = interfacialArea_[phase0Idx][phase1Idx];
192 interfacialArea_[phase0Idx][phase0Idx] = 0.;
193
194 using AnsSurface = typename Problem::SpatialParams::AnsSurface;
195 const auto& aNonWettingSolidSurfaceParams = problem.spatialParams().aNonWettingSolidSurfaceParams(element, scv, elemSol);
196 const auto ans = AnsSurface::interfacialArea(aNonWettingSolidSurfaceParams, materialParams, Sw, pc);
197
198 // Switch for using a a_{wn} relations that has some "maximum capillary pressure" as parameter
199 // That value is obtained by regularization of the pc(Sw) function.
200 static const bool computeAwsFromAnsAndPcMax = getParam<bool>("SpatialParams.ComputeAwsFromAnsAndPcMax", true);
201 if (computeAwsFromAnsAndPcMax)
202 {
203 // I know the solid surface from the pore network. But it is more consistent to use the fit value.
204 const Scalar pcMax = aWettingNonWettingSurfaceParams.pcMax();
205 const auto solidSurface = AnsSurface::interfacialArea(aNonWettingSolidSurfaceParams, materialParams, /*Sw=*/0., pcMax);
206 interfacialArea_[phase0Idx][sPhaseIdx] = solidSurface - ans;
207 }
208 else
209 {
210 using AwsSurface = typename Problem::SpatialParams::AwsSurface;
211 const auto& aWettingSolidSurfaceParams = problem.spatialParams().aWettingSolidSurfaceParams(element, scv, elemSol);
212 interfacialArea_[phase0Idx][sPhaseIdx] = AwsSurface::interfacialArea(aWettingSolidSurfaceParams, materialParams, Sw, pc);
213 }
214
215 interfacialArea_[sPhaseIdx][phase0Idx] = interfacialArea_[phase0Idx][sPhaseIdx];
216 interfacialArea_[sPhaseIdx][sPhaseIdx] = 0.;
217
218 interfacialArea_[phase1Idx][sPhaseIdx] = ans;
219 interfacialArea_[sPhaseIdx][phase1Idx] = interfacialArea_[phase1Idx][sPhaseIdx];
220 interfacialArea_[phase1Idx][phase1Idx] = 0.;
221 }
222
227 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
228 {
229 // there is no interfacial area between a phase and itself
230 assert(phaseIIdx not_eq phaseJIdx);
231 return interfacialArea_[phaseIIdx][phaseJIdx];
232 }
233
235 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
236 { return reynoldsNumber_[phaseIdx]; }
237
239 const Scalar prandtlNumber(const unsigned int phaseIdx) const
240 { return prandtlNumber_[phaseIdx]; }
241
243 const Scalar nusseltNumber(const unsigned int phaseIdx) const
244 { return nusseltNumber_[phaseIdx]; }
245
247 const Scalar characteristicLength() const
248 { return characteristicLength_; }
249
251 const Scalar factorEnergyTransfer() const
252 { return factorEnergyTransfer_; }
253
254private:
256 NumFluidPhasesArray reynoldsNumber_;
257 NumFluidPhasesArray prandtlNumber_;
258 NumFluidPhasesArray nusseltNumber_;
259
260 Scalar characteristicLength_;
261 Scalar factorEnergyTransfer_;
262 InterfacialAreasArray interfacialArea_;
263};
264
266// specialization for the case of NO kinetic mass but kinetic energy transfer of a fluid mixture and a solid
268template<class Traits, class EquilibriumVolumeVariables>
270 EquilibriumVolumeVariables,
271 false/*chemicalNonEquilibrium?*/,
272 true/*thermalNonEquilibrium?*/,
273 1>
274: public EquilibriumVolumeVariables
275{
276 using ParentType = EquilibriumVolumeVariables;
277 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
278 using Scalar = typename Traits::PrimaryVariables::value_type;
279
280 using ModelTraits = typename Traits::ModelTraits;
281 using Indices = typename ModelTraits::Indices;
282 using FS = typename Traits::FluidSystem;
283 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
284 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
285
287
288 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
289
290public:
291 using FluidState = typename Traits::FluidState;
301 template<class ElemSol, class Problem, class Element, class Scv>
302 void update(const ElemSol &elemSol,
303 const Problem &problem,
304 const Element &element,
305 const Scv& scv)
306 {
307 // Update parent type (also completes the fluid state)
308 ParentType::update(elemSol, problem, element, scv);
309
310 ParameterCache paramCache;
311 paramCache.updateAll(this->fluidState());
312 //only update of DimLessNumbers is necessary here, as interfacial area
313 // is easy due to only one fluid with a solid and is directly computed in localresidual
314 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
315 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
316 }
317
328 template<class ElemSol, class Problem, class Element, class Scv>
329 void updateDimLessNumbers(const ElemSol& elemSol,
330 const FluidState& fluidState,
331 const ParameterCache& paramCache,
332 const Problem& problem,
333 const Element& element,
334 const Scv& scv)
335 {
336 factorEnergyTransfer_ = problem.spatialParams().factorEnergyTransfer(element, scv, elemSol);
337 characteristicLength_ = problem.spatialParams().characteristicLength(element, scv, elemSol);
338
339 // set the dimensionless numbers and obtain respective quantities
340 const unsigned int vIdxGlobal = scv.dofIndex();
341 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
342 {
343 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
344 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
345 const auto density = fluidState.density(phaseIdx);
346 const auto kinematicViscosity = dynamicViscosity/density;
347
348 using FluidSystem = typename Traits::FluidSystem;
349 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
350 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
351 const auto porosity = this->porosity();
352
353 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_,kinematicViscosity);
354 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
355 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
356 prandtlNumber_[phaseIdx],
357 porosity,
358 ModelTraits::nusseltFormulation());
359 }
360 }
361
362
373 template<class ElemSol, class Problem, class Element, class Scv>
374 void updateInterfacialArea(const ElemSol& elemSol,
375 const FluidState& fluidState,
376 const ParameterCache& paramCache,
377 const Problem& problem,
378 const Element& element,
379 const Scv& scv)
380 {
381 using FluidSolidInterfacialAreaFormulation = typename Problem::SpatialParams::FluidSolidInterfacialAreaFormulation;
382 interfacialArea_ = FluidSolidInterfacialAreaFormulation::fluidSolidInterfacialArea(this->porosity(), characteristicLength());
383 }
384
386 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
387 { return reynoldsNumber_[phaseIdx]; }
388
390 const Scalar prandtlNumber(const unsigned int phaseIdx) const
391 { return prandtlNumber_[phaseIdx]; }
392
394 const Scalar nusseltNumber(const unsigned int phaseIdx) const
395 { return nusseltNumber_[phaseIdx]; }
396
398 const Scalar characteristicLength() const
399 { return characteristicLength_; }
400
402 const Scalar factorEnergyTransfer() const
403 { return factorEnergyTransfer_; }
404
405 const Scalar fluidSolidInterfacialArea() const
406 {return interfacialArea_;}
407
408private:
410 NumFluidPhasesArray reynoldsNumber_;
411 NumFluidPhasesArray prandtlNumber_;
412 NumFluidPhasesArray nusseltNumber_;
413
414 Scalar characteristicLength_;
415 Scalar factorEnergyTransfer_;
416 Scalar interfacialArea_ ;
417};
418
420// specialization for the case of (only) kinetic mass transfer. Be careful, we do not test this!
422template<class Traits, class EquilibriumVolumeVariables>
424 EquilibriumVolumeVariables,
425 true/*chemicalNonEquilibrium?*/,
426 false/*thermalNonEquilibrium?*/,
427 0>
428: public EquilibriumVolumeVariables
429{
430 using ParentType = EquilibriumVolumeVariables;
431 using FluidState = typename Traits::FluidState;
432 using FS = typename Traits::FluidSystem;
433 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
434 using Scalar = typename Traits::PrimaryVariables::value_type;
435
436 using ModelTraits = typename Traits::ModelTraits;
437
438 static constexpr auto phase0Idx = FS::phase0Idx;
439 static constexpr auto phase1Idx = FS::phase1Idx;
440 static constexpr auto wCompIdx = FS::comp0Idx;
441 static constexpr auto nCompIdx = FS::comp1Idx;
442
444
445 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
446
447public:
448 using Indices = typename ModelTraits::Indices;
458 template<class ElemSol, class Problem, class Element, class Scv>
459 void update(const ElemSol &elemSol,
460 const Problem &problem,
461 const Element &element,
462 const Scv& scv)
463 {
464 // Update parent type (also completes the fluid state)
465 ParentType::update(elemSol, problem, element, scv);
466
467 ParameterCache paramCache;
468 paramCache.updateAll(this->fluidState());
469 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
470 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
471 }
472
483 template<class ElemSol, class Problem, class Element, class Scv>
484 void updateDimLessNumbers(const ElemSol& elemSol,
485 const FluidState& fluidState,
486 const ParameterCache& paramCache,
487 const Problem& problem,
488 const Element& element,
489 const Scv& scv)
490 {
491 factorMassTransfer_ = problem.spatialParams().factorMassTransfer(element, scv, elemSol);
492 characteristicLength_ = problem.spatialParams().characteristicLength(element, scv, elemSol);
493
494 // set the dimensionless numbers and obtain respective quantities.
495 const unsigned int vIdxGlobal = scv.dofIndex();
496 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
497 {
498 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
499 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
500 const auto density = fluidState.density(phaseIdx);
501 const auto kinematicViscosity = dynamicViscosity/density;
502
503 // diffusion coefficient of non-wetting component in wetting phase
504 using FluidSystem = typename Traits::FluidSystem;
505 const auto diffCoeff = FluidSystem::binaryDiffusionCoefficient(fluidState,
506 paramCache,
507 phaseIdx,
508 wCompIdx,
509 nCompIdx);
510
511 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_, kinematicViscosity);
512 schmidtNumber_[phaseIdx] = DimLessNum::schmidtNumber(dynamicViscosity, density, diffCoeff);
513 sherwoodNumber_[phaseIdx] = DimLessNum::sherwoodNumber(reynoldsNumber_[phaseIdx],
514 schmidtNumber_[phaseIdx],
515 ModelTraits::sherwoodFormulation());
516 }
517 }
518
529 template<class ElemSol, class Problem, class Element, class Scv>
530 void updateInterfacialArea(const ElemSol& elemSol,
531 const FluidState& fluidState,
532 const ParameterCache& paramCache,
533 const Problem& problem,
534 const Element& element,
535 const Scv& scv)
536 {
537 // obtain parameters for awnsurface and material law
538 const auto& awnSurfaceParams = problem.spatialParams().aWettingNonWettingSurfaceParams(element, scv, elemSol) ;
539 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol) ;
540
541 const auto Sw = fluidState.saturation(phase0Idx) ;
542 const auto pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
543
544 // when we only consider chemical non-equilibrium there is only mass transfer between
545 // the fluid phases, so in 2p only interfacial area between wetting and non-wetting
546 using AwnSurface = typename Problem::SpatialParams::AwnSurface;
547 interfacialArea_ = AwnSurface::interfacialArea(awnSurfaceParams, materialParams, Sw, pc);
548 }
549
553 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
554 {
555 // so far there is only a model for kinetic mass transfer between fluid phases
556 assert( (phaseIIdx == phase1Idx && phaseJIdx == phase0Idx)
557 || (phaseIIdx == phase0Idx && phaseJIdx == phase1Idx) );
558 return interfacialArea_;
559 }
560
562 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
563 { return reynoldsNumber_[phaseIdx]; }
564
566 const Scalar schmidtNumber(const unsigned int phaseIdx) const
567 { return schmidtNumber_[phaseIdx]; }
568
570 const Scalar sherwoodNumber(const unsigned int phaseIdx) const
571 { return sherwoodNumber_[phaseIdx]; }
572
574 const Scalar characteristicLength() const
575 { return characteristicLength_; }
576
578 const Scalar factorMassTransfer() const
579 { return factorMassTransfer_; }
580
581private:
582 Scalar characteristicLength_;
583 Scalar factorMassTransfer_;
584 Scalar interfacialArea_ ;
585 NumFluidPhasesArray sherwoodNumber_;
586 NumFluidPhasesArray schmidtNumber_;
587 NumFluidPhasesArray reynoldsNumber_;
588};
589
590// Specialization where everything is in non-equilibrium.
591template<class Traits, class EquilibriumVolumeVariables>
593 EquilibriumVolumeVariables,
594 true/*chemicalNonEquilibrium?*/,
595 true/*thermalNonEquilibrium?*/,
596 2>
597: public EquilibriumVolumeVariables
598{
599 using ParentType = EquilibriumVolumeVariables;
600 using FluidState = typename Traits::FluidState;
601 using FS = typename Traits::FluidSystem;
602 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
603 using Scalar = typename Traits::PrimaryVariables::value_type;
604
605 using ModelTraits = typename Traits::ModelTraits;
606 using Indices = typename ModelTraits::Indices;
607 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
608 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
609
610 static constexpr auto phase0Idx = FS::phase0Idx;
611 static constexpr auto phase1Idx = FS::phase1Idx;
612 static constexpr auto sPhaseIdx = FS::numPhases;
613 static constexpr auto wCompIdx = FS::comp0Idx;
614 static constexpr auto nCompIdx = FS::comp1Idx;
615
617 static_assert((numEnergyEqFluid > 1), "This model only deals with energy transfer between two fluids and one solid phase");
618
619 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
620 using InterfacialAreasArray = std::array<std::array<Scalar, ModelTraits::numFluidPhases()+numEnergyEqSolid>,
621 ModelTraits::numFluidPhases()+numEnergyEqSolid>;
622
623public:
633 template<class ElemSol, class Problem, class Element, class Scv>
634 void update(const ElemSol &elemSol,
635 const Problem &problem,
636 const Element &element,
637 const Scv& scv)
638 {
639 // Update parent type (also completes the fluid state)
640 ParentType::update(elemSol, problem, element, scv);
641
642 ParameterCache paramCache;
643 paramCache.updateAll(this->fluidState());
644 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
645 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
646 }
647
658 template<class ElemSol, class Problem, class Element, class Scv>
659 void updateDimLessNumbers(const ElemSol& elemSol,
660 const FluidState& fluidState,
661 const ParameterCache& paramCache,
662 const Problem& problem,
663 const Element& element,
664 const Scv& scv)
665 {
666 factorMassTransfer_ = problem.spatialParams().factorMassTransfer(element, scv, elemSol);
667 factorEnergyTransfer_ = problem.spatialParams().factorEnergyTransfer(element, scv, elemSol);
668 characteristicLength_ = problem.spatialParams().characteristicLength(element, scv, elemSol);
669
670 const auto vIdxGlobal = scv.dofIndex();
671 using FluidSystem = typename Traits::FluidSystem;
672 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
673 {
674 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
675 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
676 const auto density = fluidState.density(phaseIdx);
677 const auto kinematicViscosity = dynamicViscosity/density;
678 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
679 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
680
681 // diffusion coefficient of non-wetting component in wetting phase
682 const auto porosity = this->porosity();
683 const auto diffCoeff = FluidSystem::binaryDiffusionCoefficient(fluidState,
684 paramCache,
685 phaseIdx,
686 wCompIdx,
687 nCompIdx);
688
689 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_, kinematicViscosity);
690 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
691 schmidtNumber_[phaseIdx] = DimLessNum::schmidtNumber(dynamicViscosity, density, diffCoeff);
692 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
693 prandtlNumber_[phaseIdx],
694 porosity,
695 ModelTraits::nusseltFormulation());
696 // If Diffusion is not enabled, Sherwood is divided by zero
697 sherwoodNumber_[phaseIdx] = DimLessNum::sherwoodNumber(reynoldsNumber_[phaseIdx],
698 schmidtNumber_[phaseIdx],
699 ModelTraits::sherwoodFormulation());
700 }
701 }
702
713 template<class ElemSol, class Problem, class Element, class Scv>
714 void updateInterfacialArea(const ElemSol& elemSol,
715 const FluidState& fluidState,
716 const ParameterCache& paramCache,
717 const Problem& problem,
718 const Element& element,
719 const Scv& scv)
720 {
721 // obtain (standard) material parameters (needed for the residual saturations)
722 const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
723
724 //obtain parameters for interfacial area constitutive relations
725 const auto& aWettingNonWettingSurfaceParams = problem.spatialParams().aWettingNonWettingSurfaceParams(element, scv, elemSol);
726
727 const Scalar pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
728 const Scalar Sw = fluidState.saturation(phase0Idx);
729
730 using AwnSurface = typename Problem::SpatialParams::AwnSurface;
731 const auto awn = AwnSurface::interfacialArea(aWettingNonWettingSurfaceParams, materialParams, Sw, pc);
732 interfacialArea_[phase0Idx][phase1Idx] = awn;
733 interfacialArea_[phase1Idx][phase0Idx] = interfacialArea_[phase0Idx][phase1Idx];
734 interfacialArea_[phase0Idx][phase0Idx] = 0.;
735
736 using AnsSurface = typename Problem::SpatialParams::AnsSurface;
737 const auto& aNonWettingSolidSurfaceParams = problem.spatialParams().aNonWettingSolidSurfaceParams(element, scv, elemSol);
738 const auto ans = AnsSurface::interfacialArea(aNonWettingSolidSurfaceParams, materialParams, Sw, pc);
739
740 // Switch for using a a_{wn} relations that has some "maximum capillary pressure" as parameter.
741 // That value is obtained by regularization of the pc(Sw) function.
742 static const bool computeAwsFromAnsAndPcMax = getParam<bool>("SpatialParams.ComputeAwsFromAnsAndPcMax", true);
743 if (computeAwsFromAnsAndPcMax)
744 {
745 // I know the solid surface from the pore network. But it is more consistent to use the fit value.
746 const Scalar pcMax = aWettingNonWettingSurfaceParams.pcMax();
747 const auto solidSurface = AnsSurface::interfacialArea(aNonWettingSolidSurfaceParams, materialParams, /*Sw=*/0., pcMax);
748 interfacialArea_[phase0Idx][sPhaseIdx] = solidSurface - ans;
749 }
750 else
751 {
752 using AwsSurface = typename Problem::SpatialParams::AwsSurface;
753 const auto& aWettingSolidSurfaceParams = problem.spatialParams().aWettingSolidSurfaceParams(element, scv, elemSol);
754 interfacialArea_[phase0Idx][sPhaseIdx] = AwsSurface::interfacialArea(aWettingSolidSurfaceParams, materialParams, Sw, pc);
755 }
756
757 interfacialArea_[sPhaseIdx][phase0Idx] = interfacialArea_[phase0Idx][sPhaseIdx];
758 interfacialArea_[sPhaseIdx][sPhaseIdx] = 0.;
759
760 interfacialArea_[phase1Idx][sPhaseIdx] = ans;
761 interfacialArea_[sPhaseIdx][phase1Idx] = interfacialArea_[phase1Idx][sPhaseIdx];
762 interfacialArea_[phase1Idx][phase1Idx] = 0.;
763 }
764
769 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
770 {
771 // there is no interfacial area between a phase and itself
772 assert(phaseIIdx not_eq phaseJIdx);
773 return interfacialArea_[phaseIIdx][phaseJIdx];
774 }
775
777 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
778 { return reynoldsNumber_[phaseIdx]; }
779
781 const Scalar prandtlNumber(const unsigned int phaseIdx) const
782 { return prandtlNumber_[phaseIdx]; }
783
785 const Scalar nusseltNumber(const unsigned int phaseIdx) const
786 { return nusseltNumber_[phaseIdx]; }
787
789 const Scalar schmidtNumber(const unsigned int phaseIdx) const
790 { return schmidtNumber_[phaseIdx]; }
791
793 const Scalar sherwoodNumber(const unsigned int phaseIdx) const
794 { return sherwoodNumber_[phaseIdx]; }
795
797 const Scalar characteristicLength() const
798 { return characteristicLength_; }
799
801 const Scalar factorEnergyTransfer() const
802 { return factorEnergyTransfer_; }
803
805 const Scalar factorMassTransfer() const
806 { return factorMassTransfer_; }
807
808private:
810 NumFluidPhasesArray reynoldsNumber_;
811 NumFluidPhasesArray prandtlNumber_;
812 NumFluidPhasesArray nusseltNumber_;
813 NumFluidPhasesArray schmidtNumber_;
814 NumFluidPhasesArray sherwoodNumber_;
815 Scalar characteristicLength_;
816 Scalar factorEnergyTransfer_;
817 Scalar factorMassTransfer_;
818 InterfacialAreasArray interfacialArea_;
819};
820
821} // namespace Dumux
822
823#endif
Collection of functions, calculating dimensionless numbers.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: adapt.hh:29
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
Collection of functions which calculate dimensionless numbers. Each number has it's own function....
Definition: dimensionlessnumbers.hh:64
This class contains the volume variables required for the modules which require the specific interfac...
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:49
void updateInterfacialArea(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the phases.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:172
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:251
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:239
const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
The specific interfacial area between two fluid phases [m^2 / m^3].
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:227
typename ModelTraits::Indices Indices
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:92
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Update all quantities for a given control volume.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:102
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:235
void updateDimLessNumbers(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates dimensionless numbers.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:127
const Scalar nusseltNumber(const unsigned int phaseIdx) const
access function Nusselt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:243
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:247
typename Traits::FluidState FluidState
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:91
const Scalar fluidSolidInterfacialArea() const
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:405
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:398
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:386
void updateDimLessNumbers(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates dimensionless numbers.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:329
void updateInterfacialArea(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the solid and the fluid phase.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:374
const Scalar nusseltNumber(const unsigned int phaseIdx) const
access function Nusselt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:394
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:390
typename Traits::FluidState FluidState
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:291
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:402
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Update all quantities for a given control volume.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:302
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Update all quantities for a given control volume.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:459
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:562
void updateDimLessNumbers(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the phases.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:484
typename ModelTraits::Indices Indices
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:448
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:574
const Scalar sherwoodNumber(const unsigned int phaseIdx) const
access function Sherwood Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:570
const Scalar factorMassTransfer() const
access function pre factor mass transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:578
const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
The specific interfacial area between two fluid phases [m^2 / m^3].
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:553
void updateInterfacialArea(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the phases.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:530
const Scalar schmidtNumber(const unsigned int phaseIdx) const
access function Schmidt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:566
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:797
const Scalar sherwoodNumber(const unsigned int phaseIdx) const
access function Sherwood Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:793
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Update all quantities for a given control volume.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:634
const Scalar factorMassTransfer() const
access function pre factor mass transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:805
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:801
const Scalar schmidtNumber(const unsigned int phaseIdx) const
access function Schmidt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:789
const Scalar nusseltNumber(const unsigned int phaseIdx) const
access function Nusselt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:785
const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
The specific interfacial area between two fluid phases [m^2 / m^3].
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:769
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:781
void updateDimLessNumbers(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the phases.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:659
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:777
void updateInterfacialArea(const ElemSol &elemSol, const FluidState &fluidState, const ParameterCache &paramCache, const Problem &problem, const Element &element, const Scv &scv)
Updates the volume specific interfacial area [m^2 / m^3] between the phases.
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:714