version 3.9
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
7
13#ifndef DUMUX_FREEFLOW_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
14#define DUMUX_FREEFLOW_NAVIERSTOKES_MOMENTUM_BOUNDARY_FLUXHELPER_HH
15
16#include <dumux/common/math.hh>
21
22namespace Dumux {
23
35template<class DiscretizationMethod, class SlipVelocityPolicy = void>
37
48template<class SlipVelocityPolicy>
49struct NavierStokesMomentumBoundaryFlux<DiscretizationMethods::FCStaggered, SlipVelocityPolicy>
50{
65 template<class Problem, class FVElementGeometry, class ElementVolumeVariables, class ElementFluxVariablesCache, class Scalar>
66 static auto fixedPressureMomentumFlux(const Problem& problem,
67 const FVElementGeometry& fvGeometry,
68 const typename FVElementGeometry::SubControlVolumeFace& scvf,
69 const ElementVolumeVariables& elemVolVars,
70 const ElementFluxVariablesCache& elemFluxVarsCache,
71 const Scalar pressure,
72 const bool zeroNormalVelocityGradient = true)
73 {
74 // TODO density upwinding?
75 using MomentumFlux = typename Problem::MomentumFluxType;
76 constexpr std::size_t dim = static_cast<std::size_t>(FVElementGeometry::GridGeometry::GridView::dimension);
77 static_assert(
78 MomentumFlux::size() == dim,
79 "Expects momentum flux to have as many entries as dimensions."
80 );
81
82 const auto& element = fvGeometry.element();
83
84 static_assert(
85 std::decay_t<decltype(problem.dirichlet(element, scvf))>::size() == dim,
86 "Expects problem.dirichlet to return an array with as many entries as dimensions."
87 );
88
89 if (scvf.isFrontal())
90 {
91 MomentumFlux flux(0.0);
92
93 // pressure contribution
94 flux[scvf.normalAxis()] = (pressure - problem.referencePressure(element, fvGeometry, scvf)) * scvf.directionSign();
95
96 // advective terms
97 if (problem.enableInertiaTerms())
98 {
99 const auto v = elemVolVars[scvf.insideScvIdx()].velocity();
100 flux[scvf.normalAxis()] += v*v * problem.density(element, fvGeometry, scvf) * scvf.directionSign();
101 }
102
103 return flux;
104 }
105 else if (scvf.isLateral())
106 {
107 // if no zero velocity gradient is desired the lateral flux is zero
108 if (!zeroNormalVelocityGradient)
109 return MomentumFlux(0.0);
110
111 // otherwise, make sure the flow does not diverge by accounting for the off-diagonal entries of the stress tensor
112
113 // if the lateral scvf is adjacent to a slip boundary, call the helper function for this special case
114 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
115 if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
116 return slipVelocityMomentumFlux(problem, fvGeometry, scvf, elemVolVars, elemFluxVarsCache);
117
118 // viscous terms
119 MomentumFlux flux(0.0);
120 const Scalar mu = problem.effectiveViscosity(element, fvGeometry, scvf);
121
122 // lateral face normal to boundary (integration point touches boundary)
123 if (scv.boundary())
124 flux[scv.dofAxis()] -= mu * StaggeredVelocityGradients::velocityGradIJ(fvGeometry, scvf, elemVolVars)
125 * scvf.directionSign();
126 // lateral face coinciding with boundary
127 else if (scvf.boundary())
128 flux[scv.dofAxis()] -= mu * StaggeredVelocityGradients::velocityGradJI(fvGeometry, scvf, elemVolVars)
129 * scvf.directionSign();
130
131 // advective terms
132 if (problem.enableInertiaTerms()) // TODO revise advection terms! Take care with upwinding!
133 {
134 const auto transportingVelocity = [&]()
135 {
136 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
137 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
138
139 static const bool useOldScheme = getParam<bool>("FreeFlow.UseOldTransportingVelocity", true); // TODO how to deprecate?
140
141 if (useOldScheme)
142 {
143 if (scvf.boundary())
144 {
145 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
146 return problem.dirichlet(element, scvf)[scvf.normalAxis()];
147 }
148 else
149 return innerTransportingVelocity;
150 }
151
152 // use the Dirichlet velocity as transporting velocity if the lateral face is on a Dirichlet boundary
153 if (scvf.boundary())
154 {
155 if (problem.boundaryTypes(element, scvf).isDirichlet(scvf.normalAxis()))
156 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
157 }
158
159 if (orthogonalScvf.boundary())
160 {
161 if (problem.boundaryTypes(element, orthogonalScvf).isDirichlet(scvf.normalAxis()))
162 return 0.5*(problem.dirichlet(element, scvf)[scvf.normalAxis()] + innerTransportingVelocity);
163 else
164 return innerTransportingVelocity; // fallback value, should actually never be called
165 }
166
167 // average the transporting velocity by weighting with the scv volumes
168 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
169 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
170 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
171 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
172 }();
173
174 // lateral face normal to boundary (integration point touches boundary)
175 if (scv.boundary())
176 {
177 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
178 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
179
180 const auto rho = problem.insideAndOutsideDensity(element, fvGeometry, scvf);
181 const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity);
182
183 const auto insideMomentum = innerVelocity * rho.first;
184 const auto outsideMomentum = outerVelocity * rho.second;
185
186 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
187
188 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
189 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
190
191 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
192 }
193
194 // lateral face coinciding with boundary
195 else if (scvf.boundary())
196 {
197 const auto insideDensity = problem.density(element, fvGeometry.scv(scvf.insideScvIdx()));
198 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
199 flux[scv.dofAxis()] += innerVelocity * transportingVelocity * insideDensity * scvf.directionSign();
200 }
201 }
202
203 return flux;
204 }
205
206 DUNE_THROW(Dune::InvalidStateException, "Face is neither lateral nor frontal.");
207 }
208
219 template<class Problem, class FVElementGeometry, class ElementVolumeVariables, class ElementFluxVariablesCache, class TangentialVelocityGradient = double>
220 static auto slipVelocityMomentumFlux(const Problem& problem,
221 const FVElementGeometry& fvGeometry,
222 const typename FVElementGeometry::SubControlVolumeFace& scvf,
223 const ElementVolumeVariables& elemVolVars,
224 const ElementFluxVariablesCache& elemFluxVarsCache,
225 const TangentialVelocityGradient& tangentialVelocityGradient = TangentialVelocityGradient(0.0))
226 {
227 using MomentumFlux = typename Problem::MomentumFluxType;
228 constexpr std::size_t dim = static_cast<std::size_t>(FVElementGeometry::GridGeometry::GridView::dimension);
229 static_assert(
230 MomentumFlux::size() == dim,
231 "Expects momentum flux to have as many entries as dimensions."
232 );
233
234 MomentumFlux flux(0.0);
235
236 // slip velocity momentum contribution only makes sense for lateral scvfs
237 if (scvf.isFrontal())
238 return flux;
239
240 using Scalar = std::decay_t<decltype(flux[0])>;
241 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
242 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
243 const auto& orthogonalScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
244
245 if (scvf.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(orthogonalScv)))
246 {
247 /*
248 * ---------#######::::::::: x dof position ***** porous boundary at bottom
249 * slip | || | ::
250 * gradient | || | velI :: -- element
251 * -------> | || scv x~~~~> ::
252 * ------> | || | :: O position at which gradient is evaluated (integration point)
253 * -----> | velJ ^ | :^
254 * ----> -------|-######O::::::::| <----This velocity dof (outer velJ) does not exist if the scv itself lies on a
255 * ***************~~~>****** non-Dirichlet boundary. In that case, use the given tangentialVelocityGradient.
256 * frontal scvf v_slip
257 * on porous || and # staggered half-control-volume (own element)
258 * boundary
259 * :: staggered half-control-volume (neighbor element)
260 *
261 */
262
263 if constexpr (!std::is_void_v<SlipVelocityPolicy>)
264 {
265 const Scalar velI = elemVolVars[scvf.insideScvIdx()].velocity();
266
267 // viscous terms
268 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
269 const Scalar distance = (scv.dofPosition()- scvf.ipGlobal()).two_norm();
270 const Scalar velGradJI = [&]
271 {
272 if (elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
273 return StaggeredVelocityGradients::velocityGradJI(fvGeometry, scvf, elemVolVars);
274 else
275 return tangentialVelocityGradient;
276 }();
277
278 const Scalar slipVelocity = SlipVelocityPolicy::velocity(problem, fvGeometry, scvf, elemVolVars, velGradJI)[scv.dofAxis()];
279 const Scalar velGradIJ = (slipVelocity - velI) / distance * scvf.directionSign();
280
281 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
282
283 // advective terms
284 if (problem.enableInertiaTerms())
285 {
286 // transporting velocity corresponds to velJ
287 const auto transportingVelocity = [&]()
288 {
289 const auto innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
290
291 if (!elemVolVars.hasVolVars(orthogonalScvf.outsideScvIdx()))
292 return innerTransportingVelocity; // fallback
293
294 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
295
296 // if the orthogonal scvf lies on a boundary and if there are outside volvars, we assume that these come from a Dirichlet condition
297 if (orthogonalScvf.boundary())
298 return outerTransportingVelocity;
299
300 static const bool useOldScheme = getParam<bool>("FreeFlow.UseOldTransportingVelocity", true); // TODO how to deprecate?
301 if (useOldScheme)
302 return innerTransportingVelocity;
303 else
304 {
305 // average the transporting velocity by weighting with the scv volumes
306 const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume();
307 const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume();
308 const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
309 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
310 }
311 }();
312
313 // Do not use upwinding here but directly take the slip velocity located on the boundary. Upwinding with a weight of 0.5
314 // would actually prevent second order grid convergence.
315 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
316 const auto transportedMomentum = slipVelocity * rho.second;
317
318 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
319 }
320 }
321 else
322 DUNE_THROW(Dune::InvalidStateException, "SlipVelocityPolicy needs to be specified as template argument");
323 }
324 else if (scv.boundary() && problem.onSlipBoundary(fvGeometry, fvGeometry.frontalScvfOnBoundary(scv)))
325 {
326 /* *
327 * ---------#######* x dof position ***** porous boundary at right side
328 * | || |*
329 * | || |* velI -- element
330 * | || scv x~~~~>
331 * | || |* O position at which gradient is evaluated (integration point)
332 * | velJ ^ |*^
333 * -------|-######O*| v_slip || and # staggered half-control-volume (own element)
334 * | |*
335 * | neighbor |*
336 * | element ~~~~> velI outside (Does not exist if lower later lateral scvf lies on non-Dirichlet
337 * | |* boundary. In that case, use the given tangentialVelocityGradient.)
338 * | |*
339 * ----------------* :: staggered half-control-volume (neighbor element)
340 * ^
341 * | ^
342 * slip gradient | | ^
343 * | | |
344 *
345 */
346
347 if constexpr (!std::is_void_v<SlipVelocityPolicy>)
348 {
349 const Scalar velJ = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
350
351 // viscous terms
352 const Scalar mu = problem.effectiveViscosity(fvGeometry.element(), fvGeometry, scvf);
353 const Scalar distance = (fvGeometry.scv(orthogonalScvf.insideScvIdx()).dofPosition()- scvf.ipGlobal()).two_norm();
354
355 const Scalar velGradIJ = [&]
356 {
357 if (elemVolVars.hasVolVars(scvf.outsideScvIdx()))
358 return StaggeredVelocityGradients::velocityGradIJ(fvGeometry, scvf, elemVolVars);
359 else
360 return tangentialVelocityGradient;
361 }();
362
363 const Scalar slipVelocity = SlipVelocityPolicy::velocity(problem, fvGeometry, orthogonalScvf, elemVolVars, velGradIJ)[scvf.normalAxis()];
364 const Scalar velGradJI = (slipVelocity - velJ) / distance * orthogonalScvf.directionSign();
365
366 flux[scv.dofAxis()] -= (mu * (velGradIJ + velGradJI))*scvf.directionSign();
367
368 // advective terms
369 if (problem.enableInertiaTerms())
370 {
371 // transporting velocity corresponds to velJ
372 const auto transportingVelocity = slipVelocity;
373
374 // if the scvf lies on a boundary and if there are outside volvars, we assume that these come from a Dirichlet condition
375 if (scvf.boundary() && elemVolVars.hasVolVars(scvf.outsideScvIdx()))
376 {
377 flux[scv.dofAxis()] += problem.density(fvGeometry.element(), scv)
378 * elemVolVars[scvf.outsideScvIdx()].velocity() * transportingVelocity * scvf.directionSign(); // TODO revise density
379 return flux;
380 }
381
382 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
383 const auto outerVelocity = [&]
384 {
385 if (!elemVolVars.hasVolVars(scvf.outsideScvIdx()))
386 return innerVelocity; // fallback
387 else
388 return elemVolVars[scvf.outsideScvIdx()].velocity();
389 }();
390
391 const auto rho = problem.insideAndOutsideDensity(fvGeometry.element(), fvGeometry, scvf);
392 const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity);
393
394 const auto insideMomentum = innerVelocity * rho.first;
395 const auto outsideMomentum = outerVelocity * rho.second;
396
397 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
398
399 const auto transportedMomentum = selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
400 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
401
402 flux[scv.dofAxis()] += transportingVelocity * transportedMomentum * scvf.directionSign();
403 }
404 }
405 else
406 DUNE_THROW(Dune::InvalidStateException, "SlipVelocityPolicy needs to be specified as template argument");
407 }
408
409 return flux;
410 }
411};
412
414 [[deprecated("Replace with implementation class `NavierStokesMomentumBoundaryFlux`with template arguments `DiscretizationMethod` and `SlipVelocityPolicy`. This will be removed after 3.9.")]]
417
418} // end namespace Dumux
419
420#endif
Navier Stokes slip velocity policy.
Definition: slipcondition.hh:59
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:222
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:181
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:658
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:282
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
Define some often used mathematical functions.
The available discretization methods in Dumux.
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Navier Stokes slip condition.
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:66
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:220
Class to compute the boundary flux for the momentum balance of the Navier-Stokes model This helper cl...
Definition: fluxhelper.hh:36