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