3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
25#ifndef DUMUX_FLUX_PNM_ADVECTION_HH
26#define DUMUX_FLUX_PNM_ADVECTION_HH
27
28#include <array>
30
32
33template<class... TransmissibilityLawTypes>
34struct Transmissibility : public TransmissibilityLawTypes... {};
35
36} // end namespace Dumux::PoreNetwork::Detail
37
38namespace Dumux::PoreNetwork {
39
44template<class ScalarT, class... TransmissibilityLawTypes>
46{
47
48public:
50 using Scalar = ScalarT;
51
53 using Transmissibility = Detail::Transmissibility<TransmissibilityLawTypes...>;
54
62 template<class Problem, class Element, class FVElementGeometry,
63 class ElementVolumeVariables, class SubControlVolumeFace, class ElemFluxVarsCache>
64 static Scalar flux(const Problem& problem,
65 const Element& element,
66 const FVElementGeometry& fvGeometry,
67 const ElementVolumeVariables& elemVolVars,
68 const SubControlVolumeFace& scvf,
69 const int phaseIdx,
70 const ElemFluxVarsCache& elemFluxVarsCache)
71 {
72 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
73
74 // Get the inside and outside volume variables
75 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
76 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
77 const auto& insideVolVars = elemVolVars[insideScv];
78 const auto& outsideVolVars = elemVolVars[outsideScv];
79
80 // calculate the pressure difference
81 const Scalar deltaP = insideVolVars.pressure(phaseIdx) - outsideVolVars.pressure(phaseIdx);
82 const Scalar transmissibility = fluxVarsCache.transmissibility(phaseIdx);
83 using std::isfinite;
84 assert(isfinite(transmissibility));
85
86 Scalar volumeFlow = transmissibility*deltaP;
87
88 // add gravity term
89 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
90 if (enableGravity)
91 {
92 const Scalar rho = 0.5*insideVolVars.density(phaseIdx) + 0.5*outsideVolVars.density(phaseIdx);
93 const Scalar g = problem.spatialParams().gravity(scvf.center()) * scvf.unitOuterNormal();
94
95 // The transmissibility is with respect to the effective throat length (potentially dropping the pore body radii).
96 // For gravity, we need to consider the total throat length (i.e., the cell-center to cell-center distance).
97 // This might cause some inconsistencies TODO: is there a better way?
98 volumeFlow += transmissibility * fluxVarsCache.poreToPoreDistance() * rho * g;
99 }
100
101 return volumeFlow;
102 }
103
107 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
108 static Scalar calculateTransmissibility(const Problem& problem,
109 const Element& element,
110 const FVElementGeometry& fvGeometry,
111 const typename FVElementGeometry::SubControlVolumeFace& scvf,
112 const ElementVolumeVariables& elemVolVars,
113 const FluxVariablesCache& fluxVarsCache,
114 const int phaseIdx)
115 {
116 if constexpr (ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1)
117 return Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
118 else
119 {
120 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 2);
121
122 const auto& spatialParams = problem.spatialParams();
123 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
124 const int wPhaseIdx = spatialParams.template wettingPhase<FluidSystem>(element, elemVolVars);
125 const bool invaded = fluxVarsCache.invaded();
126
127 if (phaseIdx == wPhaseIdx)
128 {
129 return invaded ? Transmissibility::wettingLayerTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
130 : Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
131 }
132 else // non-wetting phase
133 {
134 return invaded ? Transmissibility::nonWettingPhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
135 : 0.0;
136 }
137 }
138 }
139
140 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
141 static std::array<Scalar, 2> calculateTransmissibilities(const Problem& problem,
142 const Element& element,
143 const FVElementGeometry& fvGeometry,
144 const ElementVolumeVariables& elemVolVars,
145 const typename FVElementGeometry::SubControlVolumeFace& scvf,
146 const FluxVariablesCache& fluxVarsCache)
147 {
148 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
149 const Scalar t = calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, 0);
150 return {t, -t};
151 }
152};
153
160template <class ScalarT, class... TransmissibilityLawTypes>
162{
163public:
165 using Scalar = ScalarT;
166
168 using TransmissibilityCreepingFlow = Detail::Transmissibility<TransmissibilityLawTypes...>;
169
172 {
173 public:
174 class SinglePhaseCache : public TransmissibilityCreepingFlow::SinglePhaseCache
175 {
176 using ParentType = typename TransmissibilityCreepingFlow::SinglePhaseCache;
177 public:
178 using Scalar = ScalarT;
179
180 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
181 void fill(const Problem& problem,
182 const Element& element,
183 const FVElementGeometry& fvGeometry,
184 const typename FVElementGeometry::SubControlVolumeFace& scvf,
185 const ElementVolumeVariables& elemVolVars,
186 const FluxVariablesCache& fluxVarsCache,
187 const int phaseIdx)
188 {
189 ParentType::fill(problem, element, fvGeometry,scvf, elemVolVars, fluxVarsCache, phaseIdx);
190 const auto elemSol = elementSolution(element, elemVolVars, fvGeometry);
191
192 for (const auto& scv : scvs(fvGeometry))
193 {
194 const auto localIdx = scv.indexInElement();
195 const Scalar throatToPoreAreaRatio = fluxVarsCache.throatCrossSectionalArea() / problem.spatialParams().poreCrossSectionalArea(element, scv, elemSol);
196
197 // dimensionless momentum coefficient
198 const Scalar momentumCoefficient = problem.spatialParams().momentumCoefficient(element, scv, elemSol);
199
200 // dimensionless kinetik-energy coefficient
201 const Scalar kineticEnergyCoefficient = problem.spatialParams().kineticEnergyCoefficient(element, scv, elemSol);
202
203 expansionCoefficient_[localIdx] = getExpansionCoefficient_(throatToPoreAreaRatio, momentumCoefficient, kineticEnergyCoefficient);
204 contractionCoefficient_[localIdx] = getContractionCoefficient_(throatToPoreAreaRatio, momentumCoefficient, kineticEnergyCoefficient);
205 }
206 }
207
208 Scalar expansionCoefficient(int downstreamIdx) const
209 { return expansionCoefficient_[downstreamIdx]; }
210
211 Scalar contractionCoefficient(int upstreamIdx) const
212 { return contractionCoefficient_[upstreamIdx]; }
213
214 private:
215
216 Scalar getExpansionCoefficient_(const Scalar throatToPoreAreaRatio, const Scalar momentumCoefficient, const Scalar kineticEnergyCoefficient) const
217 {
218 Scalar expansionCoefficient = throatToPoreAreaRatio * throatToPoreAreaRatio * (2 * momentumCoefficient - kineticEnergyCoefficient)
219 + kineticEnergyCoefficient - 2 * momentumCoefficient * throatToPoreAreaRatio
220 - (1 - throatToPoreAreaRatio * throatToPoreAreaRatio);
221
223 }
224
225 Scalar getContractionCoefficient_(const Scalar throatToPoreAreaRatio, const Scalar momentumCoefficient, const Scalar kineticEnergyCoefficient) const
226 {
227 const Scalar contractionAreaRatio = getContractionAreaRatio_(throatToPoreAreaRatio);
228 Scalar contractionCoefficient = (1 - (throatToPoreAreaRatio * throatToPoreAreaRatio * kineticEnergyCoefficient - 2 * momentumCoefficient /*+1-throatToPoreAreaRatio*throatToPoreAreaRatio*/)
229 * contractionAreaRatio * contractionAreaRatio - 2 * contractionAreaRatio) / (contractionAreaRatio * contractionAreaRatio);
230
232 }
233
234 Scalar getContractionAreaRatio_(const Scalar throatToPoreAreaRatio) const
235 {
236 return 1.0-(1.0-throatToPoreAreaRatio)/(2.08*(1.0-throatToPoreAreaRatio)+0.5371);
237 }
238
239 std::array<Scalar, 2> expansionCoefficient_;
240 std::array<Scalar, 2> contractionCoefficient_;
241 };
242 };
243
251 template<class Problem, class Element, class FVElementGeometry,
252 class ElementVolumeVariables, class SubControlVolumeFace, class ElemFluxVarsCache>
253 static Scalar flux(const Problem& problem,
254 const Element& element,
255 const FVElementGeometry& fvGeometry,
256 const ElementVolumeVariables& elemVolVars,
257 const SubControlVolumeFace& scvf,
258 const int phaseIdx,
259 const ElemFluxVarsCache& elemFluxVarsCache)
260 {
261 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
262
263 // Get the inside and outside volume variables
264 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
265 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
266 const auto& insideVolVars = elemVolVars[insideScv];
267 const auto& outsideVolVars = elemVolVars[outsideScv];
268
269 // calculate the pressure difference
270 Scalar deltaP = insideVolVars.pressure(phaseIdx) - outsideVolVars.pressure(phaseIdx);
271
272 // add gravity term
273 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
274 if (enableGravity)
275 {
276 const Scalar rho = 0.5*insideVolVars.density(phaseIdx) + 0.5*outsideVolVars.density(phaseIdx);
277
278 // projected distance between pores in gravity field direction
279 const auto poreToPoreVector = fluxVarsCache.poreToPoreDistance() * scvf.unitOuterNormal();
280 const auto& gravityVector = problem.spatialParams().gravity(scvf.center());
281
282 deltaP += rho * (gravityVector * poreToPoreVector);
283 }
284 const Scalar creepingFlowTransmissibility = fluxVarsCache.transmissibility(phaseIdx);
285 const Scalar throatCrossSectionalArea = fluxVarsCache.throatCrossSectionalArea();
286 using std::isfinite;
287 assert(isfinite(creepingFlowTransmissibility));
288
289 assert(scvf.insideScvIdx() == 0);
290 assert(scvf.outsideScvIdx() == 1);
291
292 // determine the flow direction to predict contraction and expansion
293 const auto [upstreamIdx, downstreamIdx] = deltaP > 0 ? std::pair(0, 1) : std::pair(1, 0);
294
295 const Scalar contractionCoefficient = fluxVarsCache.singlePhaseFlowVariables().contractionCoefficient(upstreamIdx);
296 const Scalar expansionCoefficient = fluxVarsCache.singlePhaseFlowVariables().expansionCoefficient(downstreamIdx);
297 const Scalar mu = elemVolVars[upstreamIdx].viscosity();
298
308 const Scalar A = (contractionCoefficient * elemVolVars[upstreamIdx].density() + expansionCoefficient * elemVolVars[downstreamIdx].density())
309 / (2.0 * mu * mu * throatCrossSectionalArea * throatCrossSectionalArea);
310 const Scalar B = 1.0/creepingFlowTransmissibility;
311 using std::abs;
312 const Scalar C = -abs(deltaP);
313
314 using std::sqrt;
315 const auto discriminant = B*B - 4*A*C;
316 const auto q = (-B + sqrt(discriminant)) / (2*A);
317
319 if (upstreamIdx == 0)
320 return q;
321 else
322 return -q;
323
324 }
325
326 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
327 static Scalar calculateTransmissibility(const Problem& problem,
328 const Element& element,
329 const FVElementGeometry& fvGeometry,
330 const typename FVElementGeometry::SubControlVolumeFace& scvf,
331 const ElementVolumeVariables& elemVolVars,
332 const FluxVariablesCache& fluxVarsCache,
333 const int phaseIdx)
334 {
335 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
336 return Transmissibility::singlePhaseTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, phaseIdx);
337 }
338
339 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
340 static std::array<Scalar, 2> calculateTransmissibilities(const Problem& problem,
341 const Element& element,
342 const FVElementGeometry& fvGeometry,
343 const ElementVolumeVariables& elemVolVars,
344 const typename FVElementGeometry::SubControlVolumeFace& scvf,
345 const FluxVariablesCache& fluxVarsCache)
346 {
347 static_assert(ElementVolumeVariables::VolumeVariables::numFluidPhases() == 1);
348 const Scalar t = calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, fluxVarsCache, 0);
349 return {t, -t};
350 }
351
352};
353
354} // end namespace Dumux
355
356#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Definition: discretization/porenetwork/fvelementgeometry.hh:34
Definition: advection.hh:31
Hagen–Poiseuille-type flux law to describe the advective flux for pore-network models.
Definition: advection.hh:46
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:141
ScalarT Scalar
Export the Scalar type.
Definition: advection.hh:50
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:108
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:64
Non-creeping flow flux law to describe the advective flux for pore-network models based on El-Zehairy...
Definition: advection.hh:162
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:253
ScalarT Scalar
Export the Scalar type.
Definition: advection.hh:165
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:327
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:340
Inherit transmissibility from creeping flow transmissibility law to cache non-creeping flow-related p...
Definition: advection.hh:172
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:181
Scalar contractionCoefficient(int upstreamIdx) const
Definition: advection.hh:211
Scalar expansionCoefficient(int downstreamIdx) const
Definition: advection.hh:208