3.6-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 const auto& spatialParams = problem.spatialParams();
135 factorEnergyTransfer_ = spatialParams.factorEnergyTransfer(element, scv, elemSol);
136 characteristicLength_ = spatialParams.characteristicLength(element, scv, elemSol);
137
138 // set the dimensionless numbers and obtain respective quantities
139 const unsigned int vIdxGlobal = scv.dofIndex();
140 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
141 {
142 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
143 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
144 const auto density = fluidState.density(phaseIdx);
145 const auto kinematicViscosity = dynamicViscosity/density;
146
147 using FluidSystem = typename Traits::FluidSystem;
148 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
149 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
150 const auto porosity = this->porosity();
151
152 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_,kinematicViscosity);
153 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
154 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
155 prandtlNumber_[phaseIdx],
156 porosity,
157 ModelTraits::nusseltFormulation());
158
159 }
160 }
161
172 template<class ElemSol, class Problem, class Element, class Scv>
173 void updateInterfacialArea(const ElemSol& elemSol,
174 const FluidState& fluidState,
175 const ParameterCache& paramCache,
176 const Problem& problem,
177 const Element& element,
178 const Scv& scv)
179 {
180 const Scalar pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
181 const Scalar Sw = fluidState.saturation(phase0Idx);
182
183 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
184 const auto areaWN = fluidMatrixInteraction.wettingNonwettingInterface().area(Sw, pc);
185 const auto areaNS = fluidMatrixInteraction.nonwettingSolidInterface().area(Sw, pc);
186 interfacialArea_[phase0Idx][phase1Idx] = areaWN;
187 interfacialArea_[phase1Idx][phase0Idx] = interfacialArea_[phase0Idx][phase1Idx];
188 interfacialArea_[phase0Idx][phase0Idx] = 0.0;
189
190 // Switch for using a a_{wn} relations that has some "maximum capillary pressure" as parameter
191 // That value is obtained by regularization of the pc(Sw) function.
192 static const bool computeAwsFromAnsAndPcMax = getParam<bool>("SpatialParams.ComputeAwsFromAnsAndPcMax", true);
193 if (computeAwsFromAnsAndPcMax)
194 {
195 // I know the solid surface from the pore network. But it is more consistent to use the fit value.
196 const Scalar pcMax = fluidMatrixInteraction.wettingNonwettingInterface().basicParams().pcMax();
197 const auto solidSurface = fluidMatrixInteraction.nonwettingSolidInterface().area(/*Sw=*/0., pcMax);
198 interfacialArea_[phase0Idx][sPhaseIdx] = solidSurface - areaNS;
199 }
200 else
201 interfacialArea_[phase0Idx][sPhaseIdx] = fluidMatrixInteraction.wettingSolidInterface().area(Sw, pc);
202
203 interfacialArea_[sPhaseIdx][phase0Idx] = interfacialArea_[phase0Idx][sPhaseIdx];
204 interfacialArea_[sPhaseIdx][sPhaseIdx] = 0.0;
205
206 interfacialArea_[phase1Idx][sPhaseIdx] = areaNS;
207 interfacialArea_[sPhaseIdx][phase1Idx] = interfacialArea_[phase1Idx][sPhaseIdx];
208 interfacialArea_[phase1Idx][phase1Idx] = 0.0;
209 }
210
215 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
216 {
217 // there is no interfacial area between a phase and itself
218 assert(phaseIIdx not_eq phaseJIdx);
219 return interfacialArea_[phaseIIdx][phaseJIdx];
220 }
221
223 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
224 { return reynoldsNumber_[phaseIdx]; }
225
227 const Scalar prandtlNumber(const unsigned int phaseIdx) const
228 { return prandtlNumber_[phaseIdx]; }
229
231 const Scalar nusseltNumber(const unsigned int phaseIdx) const
232 { return nusseltNumber_[phaseIdx]; }
233
235 const Scalar characteristicLength() const
236 { return characteristicLength_; }
237
239 const Scalar factorEnergyTransfer() const
240 { return factorEnergyTransfer_; }
241
242private:
244 NumFluidPhasesArray reynoldsNumber_;
245 NumFluidPhasesArray prandtlNumber_;
246 NumFluidPhasesArray nusseltNumber_;
247
248 Scalar characteristicLength_;
249 Scalar factorEnergyTransfer_;
250 InterfacialAreasArray interfacialArea_;
251};
252
254// specialization for the case of NO kinetic mass but kinetic energy transfer of a fluid mixture and a solid
256template<class Traits, class EquilibriumVolumeVariables>
258 EquilibriumVolumeVariables,
259 false/*chemicalNonEquilibrium?*/,
260 true/*thermalNonEquilibrium?*/,
261 1>
262: public EquilibriumVolumeVariables
263{
264 using ParentType = EquilibriumVolumeVariables;
265 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
266 using Scalar = typename Traits::PrimaryVariables::value_type;
267
268 using ModelTraits = typename Traits::ModelTraits;
269 using FS = typename Traits::FluidSystem;
270 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
271 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
272
274
275 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
276
277public:
278 using Indices = typename ModelTraits::Indices;
279 using FluidState = typename Traits::FluidState;
280
290 template<class ElemSol, class Problem, class Element, class Scv>
291 void update(const ElemSol &elemSol,
292 const Problem &problem,
293 const Element &element,
294 const Scv& scv)
295 {
296 // Update parent type (also completes the fluid state)
297 ParentType::update(elemSol, problem, element, scv);
298
299 ParameterCache paramCache;
300 paramCache.updateAll(this->fluidState());
301 //only update of DimLessNumbers is necessary here, as interfacial area
302 // is easy due to only one fluid with a solid and is directly computed in localresidual
303 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
304 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
305 }
306
317 template<class ElemSol, class Problem, class Element, class Scv>
318 void updateDimLessNumbers(const ElemSol& elemSol,
319 const FluidState& fluidState,
320 const ParameterCache& paramCache,
321 const Problem& problem,
322 const Element& element,
323 const Scv& scv)
324 {
325 const auto& spatialParams = problem.spatialParams();
326 factorEnergyTransfer_ = spatialParams.factorEnergyTransfer(element, scv, elemSol);
327 characteristicLength_ = spatialParams.characteristicLength(element, scv, elemSol);
328
329 // set the dimensionless numbers and obtain respective quantities
330 const unsigned int vIdxGlobal = scv.dofIndex();
331 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
332 {
333 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
334 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
335 const auto density = fluidState.density(phaseIdx);
336 const auto kinematicViscosity = dynamicViscosity/density;
337
338 using FluidSystem = typename Traits::FluidSystem;
339 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
340 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
341 const auto porosity = this->porosity();
342
343 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_,kinematicViscosity);
344 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
345 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
346 prandtlNumber_[phaseIdx],
347 porosity,
348 ModelTraits::nusseltFormulation());
349 }
350 }
351
352
363 template<class ElemSol, class Problem, class Element, class Scv>
364 void updateInterfacialArea(const ElemSol& elemSol,
365 const FluidState& fluidState,
366 const ParameterCache& paramCache,
367 const Problem& problem,
368 const Element& element,
369 const Scv& scv)
370 {
371 using FluidSolidInterfacialAreaFormulation = typename Problem::SpatialParams::FluidSolidInterfacialAreaFormulation;
372 interfacialArea_ = FluidSolidInterfacialAreaFormulation::fluidSolidInterfacialArea(this->porosity(), characteristicLength());
373 }
374
376 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
377 { return reynoldsNumber_[phaseIdx]; }
378
380 const Scalar prandtlNumber(const unsigned int phaseIdx) const
381 { return prandtlNumber_[phaseIdx]; }
382
384 const Scalar nusseltNumber(const unsigned int phaseIdx) const
385 { return nusseltNumber_[phaseIdx]; }
386
388 const Scalar characteristicLength() const
389 { return characteristicLength_; }
390
392 const Scalar factorEnergyTransfer() const
393 { return factorEnergyTransfer_; }
394
395 const Scalar fluidSolidInterfacialArea() const
396 {return interfacialArea_;}
397
398private:
400 NumFluidPhasesArray reynoldsNumber_;
401 NumFluidPhasesArray prandtlNumber_;
402 NumFluidPhasesArray nusseltNumber_;
403
404 Scalar characteristicLength_;
405 Scalar factorEnergyTransfer_;
406 Scalar interfacialArea_ ;
407};
408
410// specialization for the case of (only) kinetic mass transfer. Be careful, we do not test this!
412template<class Traits, class EquilibriumVolumeVariables>
414 EquilibriumVolumeVariables,
415 true/*chemicalNonEquilibrium?*/,
416 false/*thermalNonEquilibrium?*/,
417 0>
418: public EquilibriumVolumeVariables
419{
420 using ParentType = EquilibriumVolumeVariables;
421 using FluidState = typename Traits::FluidState;
422 using FS = typename Traits::FluidSystem;
423 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
424 using Scalar = typename Traits::PrimaryVariables::value_type;
425
426 using ModelTraits = typename Traits::ModelTraits;
427
428 static constexpr auto phase0Idx = FS::phase0Idx;
429 static constexpr auto phase1Idx = FS::phase1Idx;
430 static constexpr auto wCompIdx = FS::comp0Idx;
431 static constexpr auto nCompIdx = FS::comp1Idx;
432
434
435 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
436
437public:
438 using Indices = typename ModelTraits::Indices;
448 template<class ElemSol, class Problem, class Element, class Scv>
449 void update(const ElemSol &elemSol,
450 const Problem &problem,
451 const Element &element,
452 const Scv& scv)
453 {
454 // Update parent type (also completes the fluid state)
455 ParentType::update(elemSol, problem, element, scv);
456
457 ParameterCache paramCache;
458 paramCache.updateAll(this->fluidState());
459 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
460 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
461 }
462
473 template<class ElemSol, class Problem, class Element, class Scv>
474 void updateDimLessNumbers(const ElemSol& elemSol,
475 const FluidState& fluidState,
476 const ParameterCache& paramCache,
477 const Problem& problem,
478 const Element& element,
479 const Scv& scv)
480 {
481 const auto& spatialParams = problem.spatialParams();
482 factorMassTransfer_ = spatialParams.factorMassTransfer(element, scv, elemSol);
483 characteristicLength_ = spatialParams.characteristicLength(element, scv, elemSol);
484
485 // set the dimensionless numbers and obtain respective quantities.
486 const unsigned int vIdxGlobal = scv.dofIndex();
487 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
488 {
489 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
490 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
491 const auto density = fluidState.density(phaseIdx);
492 const auto kinematicViscosity = dynamicViscosity/density;
493
494 // diffusion coefficient of nonwetting component in wetting phase
495 using FluidSystem = typename Traits::FluidSystem;
496 const auto diffCoeff = FluidSystem::binaryDiffusionCoefficient(fluidState,
497 paramCache,
498 phaseIdx,
499 wCompIdx,
500 nCompIdx);
501
502 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_, kinematicViscosity);
503 schmidtNumber_[phaseIdx] = DimLessNum::schmidtNumber(dynamicViscosity, density, diffCoeff);
504 sherwoodNumber_[phaseIdx] = DimLessNum::sherwoodNumber(reynoldsNumber_[phaseIdx],
505 schmidtNumber_[phaseIdx],
506 ModelTraits::sherwoodFormulation());
507 }
508 }
509
520 template<class ElemSol, class Problem, class Element, class Scv>
521 void updateInterfacialArea(const ElemSol& elemSol,
522 const FluidState& fluidState,
523 const ParameterCache& paramCache,
524 const Problem& problem,
525 const Element& element,
526 const Scv& scv)
527 {
528 const auto Sw = fluidState.saturation(phase0Idx) ;
529 const auto pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
530
531 // when we only consider chemical non-equilibrium there is only mass transfer between
532 // the fluid phases, so in 2p only interfacial area between wetting and non-wetting
533 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
534 interfacialArea_ = fluidMatrixInteraction.wettingNonwettingInterface().area(Sw, pc);
535 }
536
540 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
541 {
542 // so far there is only a model for kinetic mass transfer between fluid phases
543 assert( (phaseIIdx == phase1Idx && phaseJIdx == phase0Idx)
544 || (phaseIIdx == phase0Idx && phaseJIdx == phase1Idx) );
545 return interfacialArea_;
546 }
547
549 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
550 { return reynoldsNumber_[phaseIdx]; }
551
553 const Scalar schmidtNumber(const unsigned int phaseIdx) const
554 { return schmidtNumber_[phaseIdx]; }
555
557 const Scalar sherwoodNumber(const unsigned int phaseIdx) const
558 { return sherwoodNumber_[phaseIdx]; }
559
561 const Scalar characteristicLength() const
562 { return characteristicLength_; }
563
565 const Scalar factorMassTransfer() const
566 { return factorMassTransfer_; }
567
568private:
569 Scalar characteristicLength_;
570 Scalar factorMassTransfer_;
571 Scalar interfacialArea_ ;
572 NumFluidPhasesArray sherwoodNumber_;
573 NumFluidPhasesArray schmidtNumber_;
574 NumFluidPhasesArray reynoldsNumber_;
575};
576
577// Specialization where everything is in non-equilibrium.
578template<class Traits, class EquilibriumVolumeVariables>
580 EquilibriumVolumeVariables,
581 true/*chemicalNonEquilibrium?*/,
582 true/*thermalNonEquilibrium?*/,
583 2>
584: public EquilibriumVolumeVariables
585{
586 using ParentType = EquilibriumVolumeVariables;
587 using FluidState = typename Traits::FluidState;
588 using FS = typename Traits::FluidSystem;
589 using ParameterCache = typename Traits::FluidSystem::ParameterCache;
590 using Scalar = typename Traits::PrimaryVariables::value_type;
591
592 using ModelTraits = typename Traits::ModelTraits;
593 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
594 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
595
596 static constexpr auto phase0Idx = FS::phase0Idx;
597 static constexpr auto phase1Idx = FS::phase1Idx;
598 static constexpr auto sPhaseIdx = FS::numPhases;
599 static constexpr auto wCompIdx = FS::comp0Idx;
600 static constexpr auto nCompIdx = FS::comp1Idx;
601
603 static_assert((numEnergyEqFluid > 1), "This model only deals with energy transfer between two fluids and one solid phase");
604
605 using NumFluidPhasesArray = std::array<Scalar, ModelTraits::numFluidPhases()>;
606 using InterfacialAreasArray = std::array<std::array<Scalar, ModelTraits::numFluidPhases()+numEnergyEqSolid>,
607 ModelTraits::numFluidPhases()+numEnergyEqSolid>;
608
609public:
610 using Indices = typename ModelTraits::Indices;
611
621 template<class ElemSol, class Problem, class Element, class Scv>
622 void update(const ElemSol &elemSol,
623 const Problem &problem,
624 const Element &element,
625 const Scv& scv)
626 {
627 // Update parent type (also completes the fluid state)
628 ParentType::update(elemSol, problem, element, scv);
629
630 ParameterCache paramCache;
631 paramCache.updateAll(this->fluidState());
632 updateDimLessNumbers(elemSol, this->fluidState(), paramCache, problem, element, scv);
633 updateInterfacialArea(elemSol, this->fluidState(), paramCache, problem, element, scv);
634 }
635
646 template<class ElemSol, class Problem, class Element, class Scv>
647 void updateDimLessNumbers(const ElemSol& elemSol,
648 const FluidState& fluidState,
649 const ParameterCache& paramCache,
650 const Problem& problem,
651 const Element& element,
652 const Scv& scv)
653 {
654 const auto& spatialParams = problem.spatialParams();
655 factorMassTransfer_ = spatialParams.factorMassTransfer(element, scv, elemSol);
656 factorEnergyTransfer_ = spatialParams.factorEnergyTransfer(element, scv, elemSol);
657 characteristicLength_ = spatialParams.characteristicLength(element, scv, elemSol);
658
659 const auto vIdxGlobal = scv.dofIndex();
660 using FluidSystem = typename Traits::FluidSystem;
661 for (int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
662 {
663 const auto darcyMagVelocity = problem.gridVariables().volumeDarcyMagVelocity(phaseIdx, vIdxGlobal);
664 const auto dynamicViscosity = fluidState.viscosity(phaseIdx);
665 const auto density = fluidState.density(phaseIdx);
666 const auto kinematicViscosity = dynamicViscosity/density;
667 const auto heatCapacity = FluidSystem::heatCapacity(fluidState, paramCache, phaseIdx);
668 const auto thermalConductivity = FluidSystem::thermalConductivity(fluidState, paramCache, phaseIdx);
669
670 // diffusion coefficient of nonwetting component in wetting phase
671 const auto porosity = this->porosity();
672 const auto diffCoeff = FluidSystem::binaryDiffusionCoefficient(fluidState,
673 paramCache,
674 phaseIdx,
675 wCompIdx,
676 nCompIdx);
677
678 reynoldsNumber_[phaseIdx] = DimLessNum::reynoldsNumber(darcyMagVelocity, characteristicLength_, kinematicViscosity);
679 prandtlNumber_[phaseIdx] = DimLessNum::prandtlNumber(dynamicViscosity, heatCapacity, thermalConductivity);
680 schmidtNumber_[phaseIdx] = DimLessNum::schmidtNumber(dynamicViscosity, density, diffCoeff);
681 nusseltNumber_[phaseIdx] = DimLessNum::nusseltNumberForced(reynoldsNumber_[phaseIdx],
682 prandtlNumber_[phaseIdx],
683 porosity,
684 ModelTraits::nusseltFormulation());
685 // If Diffusion is not enabled, Sherwood is divided by zero
686 sherwoodNumber_[phaseIdx] = DimLessNum::sherwoodNumber(reynoldsNumber_[phaseIdx],
687 schmidtNumber_[phaseIdx],
688 ModelTraits::sherwoodFormulation());
689 }
690 }
691
702 template<class ElemSol, class Problem, class Element, class Scv>
703 void updateInterfacialArea(const ElemSol& elemSol,
704 const FluidState& fluidState,
705 const ParameterCache& paramCache,
706 const Problem& problem,
707 const Element& element,
708 const Scv& scv)
709 {
710 const Scalar pc = fluidState.pressure(phase1Idx) - fluidState.pressure(phase0Idx);
711 const Scalar Sw = fluidState.saturation(phase0Idx);
712
713 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
714
715 const auto awn = fluidMatrixInteraction.wettingNonwettingInterface().area(Sw, pc);
716 interfacialArea_[phase0Idx][phase1Idx] = awn;
717 interfacialArea_[phase1Idx][phase0Idx] = interfacialArea_[phase0Idx][phase1Idx];
718 interfacialArea_[phase0Idx][phase0Idx] = 0.;
719
720 const auto areaNS = fluidMatrixInteraction.nonwettingSolidInterface().area(Sw, pc);
721
722 // Switch for using a a_{wn} relations that has some "maximum capillary pressure" as parameter.
723 // That value is obtained by regularization of the pc(Sw) function.
724 static const bool computeAwsFromAnsAndPcMax = getParam<bool>("SpatialParams.ComputeAwsFromAnsAndPcMax", true);
725 if (computeAwsFromAnsAndPcMax)
726 {
727 // I know the solid surface from the pore network. But it is more consistent to use the fit value.
728 const Scalar pcMax = fluidMatrixInteraction.wettingNonwettingInterface().basicParams().pcMax();
729 const auto solidSurface = fluidMatrixInteraction.nonwettingSolidInterface().area(/*Sw=*/0., pcMax);
730 interfacialArea_[phase0Idx][sPhaseIdx] = solidSurface - areaNS;
731 }
732 else
733 interfacialArea_[phase0Idx][sPhaseIdx] = fluidMatrixInteraction.wettingSolidInterface().area(Sw, pc);
734
735 interfacialArea_[sPhaseIdx][phase0Idx] = interfacialArea_[phase0Idx][sPhaseIdx];
736 interfacialArea_[sPhaseIdx][sPhaseIdx] = 0.;
737
738 interfacialArea_[phase1Idx][sPhaseIdx] = areaNS;
739 interfacialArea_[sPhaseIdx][phase1Idx] = interfacialArea_[phase1Idx][sPhaseIdx];
740 interfacialArea_[phase1Idx][phase1Idx] = 0.;
741 }
742
747 const Scalar interfacialArea(const unsigned int phaseIIdx, const unsigned int phaseJIdx) const
748 {
749 // there is no interfacial area between a phase and itself
750 assert(phaseIIdx not_eq phaseJIdx);
751 return interfacialArea_[phaseIIdx][phaseJIdx];
752 }
753
755 const Scalar reynoldsNumber(const unsigned int phaseIdx) const
756 { return reynoldsNumber_[phaseIdx]; }
757
759 const Scalar prandtlNumber(const unsigned int phaseIdx) const
760 { return prandtlNumber_[phaseIdx]; }
761
763 const Scalar nusseltNumber(const unsigned int phaseIdx) const
764 { return nusseltNumber_[phaseIdx]; }
765
767 const Scalar schmidtNumber(const unsigned int phaseIdx) const
768 { return schmidtNumber_[phaseIdx]; }
769
771 const Scalar sherwoodNumber(const unsigned int phaseIdx) const
772 { return sherwoodNumber_[phaseIdx]; }
773
775 const Scalar characteristicLength() const
776 { return characteristicLength_; }
777
779 const Scalar factorEnergyTransfer() const
780 { return factorEnergyTransfer_; }
781
783 const Scalar factorMassTransfer() const
784 { return factorMassTransfer_; }
785
786private:
788 NumFluidPhasesArray reynoldsNumber_;
789 NumFluidPhasesArray prandtlNumber_;
790 NumFluidPhasesArray nusseltNumber_;
791 NumFluidPhasesArray schmidtNumber_;
792 NumFluidPhasesArray sherwoodNumber_;
793 Scalar characteristicLength_;
794 Scalar factorEnergyTransfer_;
795 Scalar factorMassTransfer_;
796 InterfacialAreasArray interfacialArea_;
797};
798
799} // namespace Dumux
800
801#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Collection of functions, calculating dimensionless numbers.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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:65
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:173
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:239
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:227
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:215
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:223
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:231
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:235
typename Traits::FluidState FluidState
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:91
typename ModelTraits::Indices Indices
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:278
const Scalar fluidSolidInterfacialArea() const
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:395
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:388
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:376
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:318
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:364
const Scalar nusseltNumber(const unsigned int phaseIdx) const
access function Nusselt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:384
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:380
typename Traits::FluidState FluidState
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:279
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:392
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:291
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:449
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:549
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:474
typename ModelTraits::Indices Indices
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:438
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:561
const Scalar sherwoodNumber(const unsigned int phaseIdx) const
access function Sherwood Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:557
const Scalar factorMassTransfer() const
access function pre factor mass transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:565
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:540
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:521
const Scalar schmidtNumber(const unsigned int phaseIdx) const
access function Schmidt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:553
const Scalar characteristicLength() const
access function characteristic length
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:775
const Scalar sherwoodNumber(const unsigned int phaseIdx) const
access function Sherwood Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:771
typename ModelTraits::Indices Indices
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:610
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:622
const Scalar factorMassTransfer() const
access function pre factor mass transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:783
const Scalar factorEnergyTransfer() const
access function pre factor energy transfer
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:779
const Scalar schmidtNumber(const unsigned int phaseIdx) const
access function Schmidt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:767
const Scalar nusseltNumber(const unsigned int phaseIdx) const
access function Nusselt Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:763
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:747
const Scalar prandtlNumber(const unsigned int phaseIdx) const
access function Prandtl Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:759
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:647
const Scalar reynoldsNumber(const unsigned int phaseIdx) const
access function Reynolds Number
Definition: porousmediumflow/nonequilibrium/volumevariables.hh:755
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:703