version 3.8
advection.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_FLUX_PNM_ADVECTION_HH
14#define DUMUX_FLUX_PNM_ADVECTION_HH
15
16#include <array>
18
20
21template<class... TransmissibilityLawTypes>
22struct Transmissibility : public TransmissibilityLawTypes... {};
23
24} // end namespace Dumux::PoreNetwork::Detail
25
26namespace Dumux::PoreNetwork {
27
32template<class ScalarT, class... TransmissibilityLawTypes>
34{
35
36public:
38 using Scalar = ScalarT;
39
41 using Transmissibility = Detail::Transmissibility<TransmissibilityLawTypes...>;
42
50 template<class Problem, class Element, class FVElementGeometry,
51 class ElementVolumeVariables, class SubControlVolumeFace, class ElemFluxVarsCache>
52 static Scalar flux(const Problem& problem,
53 const Element& element,
54 const FVElementGeometry& fvGeometry,
55 const ElementVolumeVariables& elemVolVars,
56 const SubControlVolumeFace& scvf,
57 const int phaseIdx,
58 const ElemFluxVarsCache& elemFluxVarsCache)
59 {
60 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
61
62 // Get the inside and outside volume variables
63 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
64 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
65 const auto& insideVolVars = elemVolVars[insideScv];
66 const auto& outsideVolVars = elemVolVars[outsideScv];
67
68 // calculate the pressure difference
69 const Scalar deltaP = insideVolVars.pressure(phaseIdx) - outsideVolVars.pressure(phaseIdx);
70 const Scalar transmissibility = fluxVarsCache.transmissibility(phaseIdx);
71 using std::isfinite;
72 assert(isfinite(transmissibility));
73
74 Scalar volumeFlow = transmissibility*deltaP;
75
76 // add gravity term
77 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
78 if (enableGravity)
79 {
80 const Scalar rho = 0.5*insideVolVars.density(phaseIdx) + 0.5*outsideVolVars.density(phaseIdx);
81 const Scalar g = problem.spatialParams().gravity(scvf.center()) * scvf.unitOuterNormal();
82
83 // The transmissibility is with respect to the effective throat length (potentially dropping the pore body radii).
84 // For gravity, we need to consider the total throat length (i.e., the cell-center to cell-center distance).
85 // This might cause some inconsistencies TODO: is there a better way?
86 volumeFlow += transmissibility * fluxVarsCache.poreToPoreDistance() * rho * g;
87 }
88
89 return volumeFlow;
90 }
91
95 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
96 static Scalar calculateTransmissibility(const Problem& problem,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const typename FVElementGeometry::SubControlVolumeFace& scvf,
100 const ElementVolumeVariables& elemVolVars,
101 const FluxVariablesCache& fluxVarsCache,
102 const int phaseIdx)
103 {
104 if constexpr (ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1)
105 return Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
106 else
107 {
108 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 2);
109
110 const auto& spatialParams = problem.spatialParams();
111 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
112 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, elemVolVars);
113 const bool invaded = fluxVarsCache.invaded();
114
115 if (phaseIdx == wPhaseIdx)
116 {
117 return invaded ? Transmissibility::wettingLayerTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
118 : Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
119 }
120 else // non-wetting phase
121 {
122 return invaded ? Transmissibility::nonWettingPhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
123 : 0.0;
124 }
125 }
126 }
127
128 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
129 static std::array<Scalar, 2> calculateTransmissibilities(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const typename FVElementGeometry::SubControlVolumeFace& scvf,
134 const FluxVariablesCache& fluxVarsCache)
135 {
136 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
137 const Scalar t = calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, 0);
138 return {t, -t};
139 }
140};
141
148template <class ScalarT, class... TransmissibilityLawTypes>
150{
151public:
153 using Scalar = ScalarT;
154
156 using TransmissibilityCreepingFlow = Detail::Transmissibility<TransmissibilityLawTypes...>;
157
160 {
161 public:
162 class SinglePhaseCache : public TransmissibilityCreepingFlow::SinglePhaseCache
163 {
164 using ParentType = typename TransmissibilityCreepingFlow::SinglePhaseCache;
165 public:
166 using Scalar = ScalarT;
167
168 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
169 void fill(const Problem& problem,
170 const Element& element,
171 const FVElementGeometry& fvGeometry,
172 const typename FVElementGeometry::SubControlVolumeFace& scvf,
173 const ElementVolumeVariables& elemVolVars,
174 const FluxVariablesCache& fluxVarsCache,
175 const int phaseIdx)
176 {
177 ParentType::fill(problem, element, fvGeometry,scvf, elemVolVars, fluxVarsCache, phaseIdx);
178 const auto elemSol = elementSolution(element, elemVolVars, fvGeometry);
179
180 for (const auto& scv : scvs(fvGeometry))
181 {
182 const auto localIdx = scv.indexInElement();
183 const Scalar throatToPoreAreaRatio = fluxVarsCache.throatCrossSectionalArea() / problem.spatialParams().poreCrossSectionalArea(element, scv, elemSol);
184
185 // dimensionless momentum coefficient
186 const Scalar momentumCoefficient = problem.spatialParams().momentumCoefficient(element, scv, elemSol);
187
188 // dimensionless kinetik-energy coefficient
189 const Scalar kineticEnergyCoefficient = problem.spatialParams().kineticEnergyCoefficient(element, scv, elemSol);
190
191 expansionCoefficient_[localIdx] = getExpansionCoefficient_(throatToPoreAreaRatio, momentumCoefficient, kineticEnergyCoefficient);
192 contractionCoefficient_[localIdx] = getContractionCoefficient_(throatToPoreAreaRatio, momentumCoefficient, kineticEnergyCoefficient);
193 }
194 }
195
196 Scalar expansionCoefficient(int downstreamIdx) const
197 { return expansionCoefficient_[downstreamIdx]; }
198
199 Scalar contractionCoefficient(int upstreamIdx) const
200 { return contractionCoefficient_[upstreamIdx]; }
201
202 private:
203
204 Scalar getExpansionCoefficient_(const Scalar throatToPoreAreaRatio, const Scalar momentumCoefficient, const Scalar kineticEnergyCoefficient) const
205 {
206 Scalar expansionCoefficient = throatToPoreAreaRatio * throatToPoreAreaRatio * (2 * momentumCoefficient - kineticEnergyCoefficient)
207 + kineticEnergyCoefficient - 2 * momentumCoefficient * throatToPoreAreaRatio
208 - (1 - throatToPoreAreaRatio * throatToPoreAreaRatio);
209
211 }
212
213 Scalar getContractionCoefficient_(const Scalar throatToPoreAreaRatio, const Scalar momentumCoefficient, const Scalar kineticEnergyCoefficient) const
214 {
215 const Scalar contractionAreaRatio = getContractionAreaRatio_(throatToPoreAreaRatio);
216 Scalar contractionCoefficient = (1 - (throatToPoreAreaRatio * throatToPoreAreaRatio * kineticEnergyCoefficient - 2 * momentumCoefficient /*+1-throatToPoreAreaRatio*throatToPoreAreaRatio*/)
217 * contractionAreaRatio * contractionAreaRatio - 2 * contractionAreaRatio) / (contractionAreaRatio * contractionAreaRatio);
218
220 }
221
222 Scalar getContractionAreaRatio_(const Scalar throatToPoreAreaRatio) const
223 {
224 return 1.0-(1.0-throatToPoreAreaRatio)/(2.08*(1.0-throatToPoreAreaRatio)+0.5371);
225 }
226
227 std::array<Scalar, 2> expansionCoefficient_;
228 std::array<Scalar, 2> contractionCoefficient_;
229 };
230 };
231
239 template<class Problem, class Element, class FVElementGeometry,
240 class ElementVolumeVariables, class SubControlVolumeFace, class ElemFluxVarsCache>
241 static Scalar flux(const Problem& problem,
242 const Element& element,
243 const FVElementGeometry& fvGeometry,
244 const ElementVolumeVariables& elemVolVars,
245 const SubControlVolumeFace& scvf,
246 const int phaseIdx,
247 const ElemFluxVarsCache& elemFluxVarsCache)
248 {
249 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
250
251 // Get the inside and outside volume variables
252 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
253 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
254 const auto& insideVolVars = elemVolVars[insideScv];
255 const auto& outsideVolVars = elemVolVars[outsideScv];
256
257 // calculate the pressure difference
258 Scalar deltaP = insideVolVars.pressure(phaseIdx) - outsideVolVars.pressure(phaseIdx);
259
260 // add gravity term
261 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
262 if (enableGravity)
263 {
264 const Scalar rho = 0.5*insideVolVars.density(phaseIdx) + 0.5*outsideVolVars.density(phaseIdx);
265
266 // projected distance between pores in gravity field direction
267 const auto poreToPoreVector = fluxVarsCache.poreToPoreDistance() * scvf.unitOuterNormal();
268 const auto& gravityVector = problem.spatialParams().gravity(scvf.center());
269
270 deltaP += rho * (gravityVector * poreToPoreVector);
271 }
272 const Scalar creepingFlowTransmissibility = fluxVarsCache.transmissibility(phaseIdx);
273 const Scalar throatCrossSectionalArea = fluxVarsCache.throatCrossSectionalArea();
274 using std::isfinite;
275 assert(isfinite(creepingFlowTransmissibility));
276
277 assert(scvf.insideScvIdx() == 0);
278 assert(scvf.outsideScvIdx() == 1);
279
280 // determine the flow direction to predict contraction and expansion
281 const auto [upstreamIdx, downstreamIdx] = deltaP > 0 ? std::pair(0, 1) : std::pair(1, 0);
282
283 const Scalar contractionCoefficient = fluxVarsCache.singlePhaseFlowVariables().contractionCoefficient(upstreamIdx);
284 const Scalar expansionCoefficient = fluxVarsCache.singlePhaseFlowVariables().expansionCoefficient(downstreamIdx);
285 const Scalar mu = elemVolVars[upstreamIdx].viscosity();
286
296 const Scalar A = (contractionCoefficient * elemVolVars[upstreamIdx].density() + expansionCoefficient * elemVolVars[downstreamIdx].density())
297 / (2.0 * mu * mu * throatCrossSectionalArea * throatCrossSectionalArea);
298 const Scalar B = 1.0/creepingFlowTransmissibility;
299 using std::abs;
300 const Scalar C = -abs(deltaP);
301
302 using std::sqrt;
303 const auto discriminant = B*B - 4*A*C;
304 const auto q = (-B + sqrt(discriminant)) / (2*A);
305
307 if (upstreamIdx == 0)
308 return q;
309 else
310 return -q;
311
312 }
313
314 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
315 static Scalar calculateTransmissibility(const Problem& problem,
316 const Element& element,
317 const FVElementGeometry& fvGeometry,
318 const typename FVElementGeometry::SubControlVolumeFace& scvf,
319 const ElementVolumeVariables& elemVolVars,
320 const FluxVariablesCache& fluxVarsCache,
321 const int phaseIdx)
322 {
323 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
324 return Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
325 }
326
327 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
328 static std::array<Scalar, 2> calculateTransmissibilities(const Problem& problem,
329 const Element& element,
330 const FVElementGeometry& fvGeometry,
331 const ElementVolumeVariables& elemVolVars,
332 const typename FVElementGeometry::SubControlVolumeFace& scvf,
333 const FluxVariablesCache& fluxVarsCache)
334 {
335 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
336 const Scalar t = calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, 0);
337 return {t, -t};
338 }
339
340};
341
342} // end namespace Dumux
343
344#endif
Hagen–Poiseuille-type flux law to describe the advective flux for pore-network models.
Definition: advection.hh:34
static std::array< Scalar, 2 > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: advection.hh:129
ScalarT Scalar
Export the Scalar type.
Definition: advection.hh:38
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Returns the throat conductivity.
Definition: advection.hh:96
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElemFluxVarsCache &elemFluxVarsCache)
Returns the advective flux of a fluid phase across the given sub-control volume face (corresponding t...
Definition: advection.hh:52
void fill(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: advection.hh:169
Scalar contractionCoefficient(int upstreamIdx) const
Definition: advection.hh:199
Scalar expansionCoefficient(int downstreamIdx) const
Definition: advection.hh:196
Inherit transmissibility from creeping flow transmissibility law to cache non-creeping flow-related p...
Definition: advection.hh:160
Non-creeping flow flux law to describe the advective flux for pore-network models based on El-Zehairy...
Definition: advection.hh:150
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElemFluxVarsCache &elemFluxVarsCache)
Returns the advective flux of a fluid phase across the given sub-control volume face (corresponding t...
Definition: advection.hh:241
ScalarT Scalar
Export the Scalar type.
Definition: advection.hh:153
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: advection.hh:315
static std::array< Scalar, 2 > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: advection.hh:328
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
Definition: advection.hh:19
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.