3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fluxhelper.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 *****************************************************************************/
19
25#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
26#define DUMUX_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
27
28#include <dumux/common/math.hh>
32
33namespace Dumux {
34
40{
55 template<class Problem, class FVElementGeometry, class ElementVolumeVariables, class ElementFluxVariablesCache, class Scalar>
56 static auto fixedPressureMomentumFlux(const Problem& problem,
57 const FVElementGeometry& fvGeometry,
58 const typename FVElementGeometry::SubControlVolumeFace& scvf,
59 const ElementVolumeVariables& elemVolVars,
60 const ElementFluxVariablesCache& elemFluxVarsCache,
61 const Scalar pressure,
62 const bool zeroNormalVelocityGradient = true)
63 {
64 // TODO density upwinding?
65 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered);
66 using MomentumFlux = typename Problem::MomentumFluxType;
67 constexpr std::size_t dim = static_cast<std::size_t>(FVElementGeometry::GridGeometry::GridView::dimension);
68 static_assert(
69 MomentumFlux::size() == dim,
70 "Expects momentum flux to have as many entries as dimensions."
71 );
72
73 const auto& element = fvGeometry.element();
74
75 static_assert(
76 std::decay_t<decltype(problem.dirichlet(element, scvf))>::size() == dim,
77 "Expects problem.dirichlet to return an array with as many entries as dimensions."
78 );
79
80 if (scvf.isFrontal())
81 {
82 MomentumFlux flux(0.0);
83
84 // pressure contribution
85 flux[scvf.normalAxis()] = (pressure - problem.referencePressure(element, fvGeometry, scvf)) * scvf.directionSign();
86
87 // advective terms
88 if (problem.enableInertiaTerms())
89 {
90 const auto v = elemVolVars[scvf.insideScvIdx()].velocity();
91 flux[scvf.normalAxis()] += v*v * problem.density(element, fvGeometry, scvf) * scvf.directionSign();
92 }
93
94 return flux;
95 }
96 else if (scvf.isLateral())
97 {
98 // if no zero velocity gradient is desired the lateral flux is zero
99 if (!zeroNormalVelocityGradient)
100 return MomentumFlux(0.0);
101
102 // otherwise, make sure the flow does not diverge by accounting for the off-diagonal entries of the stress tensor
103
104 // if the lateral scvf is adjacent to a slip boundary, call the helper function for this special case
105 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
106 if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
107 return slipVelocityMomentumFlux(problem, fvGeometry, scvf, elemVolVars, elemFluxVarsCache);
108
109 // viscous terms
110 MomentumFlux flux(0.0);
111 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, scvf);
112
113 // lateral face normal to boundary (integration point touches boundary)
114 if (scv.boundary())
115 flux[scv.dofAxis()] -= mu * StaggeredVelocityGradients::velocityGradIJ(fvGeometry, scvf, elemVolVars)
116 * scvf.directionSign();
117 // lateral face coinciding with boundary
118 else if (scvf.boundary())
119 flux[scv.dofAxis()] -= mu * StaggeredVelocityGradients::velocityGradJI(fvGeometry, scvf, elemVolVars)
120 * scvf.directionSign();
121
122 // advective terms
123 if (problem.enableInertiaTerms()) // TODO revise advection terms! Take care with upwinding!
124 {
125 const auto transportingVelocity = [&]()
126 {
127 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
128 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
129
130 static const bool useOldScheme = getParam<bool>("FreeFlow.UseOldTransportingVelocity", true); // TODO how to deprecate?
131
132 if (useOldScheme)
133 {
134 if (scvf.boundary())
135 {
136 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
137 return problem.dirichlet(element, scvf)[scvf.normalAxis()];
138 }
139 else
140 return innerTransportingVelocity;
141 }
142
143 // use the Dirichlet velocity as transporting velocity if the lateral face is on a Dirichlet boundary
144 if (scvf.boundary())
145 {
146 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
147 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
148 }
149
150 if (orthogonalScvf.boundary())
151 {
152 if (problem.boundaryTypes(element, orthogonalScvf).isDirichlet(scvf.normalAxis()))
153 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
154 else
155 return innerTransportingVelocity; // fallback value, should actually never be called
156 }
157
158 // average the transporting velocity by weighting with the scv volumes
159 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
160 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
161 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
162 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
163 }();
164
165 // lateral face normal to boundary (integration point touches boundary)
166 if (scv.boundary())
167 {
168 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
169 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
170
171 const auto rho = problem.insideAndOutsideDensity(element, fvGeometry, scvf);
172 const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity);
173
174 const auto insideMomentum = innerVelocity * rho.first;
175 const auto outsideMomentum = outerVelocity * rho.second;
176
177 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
178
179 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
180 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
181
182 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
183 }
184
185 // lateral face coinciding with boundary
186 else if (scvf.boundary())
187 {
188 const auto insideDensity = problem.density(element, fvGeometry.scv(scvf.insideScvIdx()));
189 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
190 flux[scv.dofAxis()] += innerVelocity * transportingVelocity * insideDensity * scvf.directionSign();
191 }
192 }
193
194 return flux;
195 }
196
197 DUNE_THROW(Dune::InvalidStateException, "Face is neither lateral nor frontal.");
198 }
199
210 template<class Problem, class FVElementGeometry, class ElementVolumeVariables, class ElementFluxVariablesCache, class TangentialVelocityGradient = double>
211 static auto slipVelocityMomentumFlux(const Problem& problem,
212 const FVElementGeometry& fvGeometry,
213 const typename FVElementGeometry::SubControlVolumeFace& scvf,
214 const ElementVolumeVariables& elemVolVars,
215 const ElementFluxVariablesCache& elemFluxVarsCache,
216 const TangentialVelocityGradient& tangentialVelocityGradient = TangentialVelocityGradient(0.0))
217 {
218 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered);
219 using MomentumFlux = typename Problem::MomentumFluxType;
220 constexpr std::size_t dim = static_cast<std::size_t>(FVElementGeometry::GridGeometry::GridView::dimension);
221 static_assert(
222 MomentumFlux::size() == dim,
223 "Expects momentum flux to have as many entries as dimensions."
224 );
225
226 MomentumFlux flux(0.0);
227
228 // slip velocity momentum contribution only makes sense for lateral scvfs
229 if (scvf.isFrontal())
230 return flux;
231
232 using Scalar = std::decay_t<decltype(flux[0])>;
233 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
234 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
235 const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
236
237 if (scvf.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(orthogonalScv)))
238 {
239 /*
240 * ---------#######::::::::: x dof position ***** porous boundary at bottom
241 * slip | || | ::
242 * gradient | || | velI :: -- element
243 * -------> | || scv x~~~~> ::
244 * ------> | || | :: O position at which gradient is evaluated (integration point)
245 * -----> | velJ ^ | :^
246 * ----> -------|-######O::::::::| <----This velocity dof (outer velJ) does not exist if the scv itself lies on a
247 * ***************~~~>****** non-Dirichlet boundary. In that case, use the given tangentialVelocityGradient.
248 * frontal scvf v_slip
249 * on porous || and # staggered half-control-volume (own element)
250 * boundary
251 * :: staggered half-control-volume (neighbor element)
252 *
253 */
254 const Scalar velI = elemVolVars[scvf.insideScvIdx()].velocity();
255
256 // viscous terms
257 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
258 const Scalar distance = (scv.dofPosition()- scvf.ipGlobal()).two_norm();
259 const Scalar velGradJI = [&]
260 {
261 if (elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
262 return StaggeredVelocityGradients::velocityGradJI(fvGeometry, scvf, elemVolVars);
263 else
264 return tangentialVelocityGradient;
265 }();
266
267 const Scalar slipVelocity = problem.beaversJosephVelocity(fvGeometry, scvf, elemVolVars, velGradJI)[scv.dofAxis()]; // TODO rename to slipVelocity
268 const Scalar velGradIJ = (slipVelocity - velI) / distance * scvf.directionSign();
269
270 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
271
272 // advective terms
273 if (problem.enableInertiaTerms())
274 {
275 // transporting velocity corresponds to velJ
276 const auto transportingVelocity = [&]()
277 {
278 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
279
280 if (!elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
281 return innerTransportingVelocity; // fallback
282
283 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
284
285 // if the orthogonal scvf lies on a boundary and if there are outside volvars, we assume that these come from a Dirichlet condition
286 if (orthogonalScvf.boundary())
287 return outerTransportingVelocity;
288
289 static const bool useOldScheme = getParam<bool>("FreeFlow.UseOldTransportingVelocity", true); // TODO how to deprecate?
290 if (useOldScheme)
291 return innerTransportingVelocity;
292 else
293 {
294 // average the transporting velocity by weighting with the scv volumes
295 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
296 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
297 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
298 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
299 }
300 }();
301
302 // Do not use upwinding here but directly take the slip velocity located on the boundary. Upwinding with a weight of 0.5
303 // would actually prevent second order grid convergence.
304 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
305 const auto transportedMomentum = slipVelocity * rho.second;
306
307 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
308 }
309 }
310 else if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
311 {
312 /* *
313 * ---------#######* x dof position ***** porous boundary at right side
314 * | || |*
315 * | || |* velI -- element
316 * | || scv x~~~~>
317 * | || |* O position at which gradient is evaluated (integration point)
318 * | velJ ^ |*^
319 * -------|-######O*| v_slip || and # staggered half-control-volume (own element)
320 * | |*
321 * | neighbor |*
322 * | element ~~~~> velI outside (Does not exist if lower later lateral scvf lies on non-Dirichlet
323 * | |* boundary. In that case, use the given tangentialVelocityGradient.)
324 * | |*
325 * ----------------* :: staggered half-control-volume (neighbor element)
326 * ^
327 * | ^
328 * slip gradient | | ^
329 * | | |
330 *
331 */
332 const Scalar velJ = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
333
334 // viscous terms
335 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
336 const Scalar distance = (fvGeometry.scv(orthogonalScvf.insideScvIdx()).dofPosition()- scvf.ipGlobal()).two_norm();
337
338 const Scalar velGradIJ = [&]
339 {
340 if (elemVolVars.hasVolVars(scvf.outsideScvIdx()))
341 return StaggeredVelocityGradients::velocityGradIJ(fvGeometry, scvf, elemVolVars);
342 else
343 return tangentialVelocityGradient;
344 }();
345
346 const Scalar slipVelocity = problem.beaversJosephVelocity(fvGeometry, orthogonalScvf, elemVolVars, velGradIJ)[scvf.normalAxis()]; // TODO rename to slipVelocity
347 const Scalar velGradJI = (slipVelocity - velJ) / distance * orthogonalScvf.directionSign();
348
349 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
350
351 // advective terms
352 if (problem.enableInertiaTerms())
353 {
354 // transporting velocity corresponds to velJ
355 const auto transportingVelocity = slipVelocity;
356
357 // if the scvf lies on a boundary and if there are outside volvars, we assume that these come from a Dirichlet condition
358 if (scvf.boundary() && elemVolVars.hasVolVars(scvf.outsideScvIdx()))
359 {
360 flux[scv.dofAxis()] += problem.density(fvGeometry.element(), scv)
361 * elemVolVars[scvf.outsideScvIdx()].velocity() * transportingVelocity * scvf.directionSign(); // TODO revise density
362 return flux;
363 }
364
365 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
366 const auto outerVelocity = [&]
367 {
368 if (!elemVolVars.hasVolVars(scvf.outsideScvIdx()))
369 return innerVelocity; // fallback
370 else
371 return elemVolVars[scvf.outsideScvIdx()].velocity();
372 }();
373
374 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
375 const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity);
376
377 const auto insideMomentum = innerVelocity * rho.first;
378 const auto outsideMomentum = outerVelocity * rho.second;
379
380 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
381
382 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
383 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
384
385 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
386 }
387 }
388
389 return flux;
390 }
391};
392
393} // end namespace Dumux
394
395#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
The available discretization methods in Dumux.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr FCStaggered fcstaggered
Definition: method.hh:140
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Struct containing flux helper functions to be used in the momentum problem's Neumann function.
Definition: fluxhelper.hh:40
static auto slipVelocityMomentumFlux(const Problem &problem, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const TangentialVelocityGradient &tangentialVelocityGradient=TangentialVelocityGradient(0.0))
Returns the momentum flux for a boundary with a slip condition.
Definition: fluxhelper.hh:211
static auto fixedPressureMomentumFlux(const Problem &problem, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const Scalar pressure, const bool zeroNormalVelocityGradient=true)
Returns the momentum flux a fixed-pressure boundary.
Definition: fluxhelper.hh:56
static auto velocityGradJI(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the velocity gradient in line with our current scvf.
Definition: momentum/velocitygradients.hh:234
static auto velocityGradIJ(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the velocity gradient perpendicular to the orientation of our current scvf.
Definition: momentum/velocitygradients.hh:193