version 3.8
staggered/velocitygradients.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//
12#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
13#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
14
15#include <optional>
18
19namespace Dumux {
20
25template<class Scalar, class GridGeometry, class BoundaryTypes, class Indices>
26class StaggeredVelocityGradients
27{
28 using FVElementGeometry = typename GridGeometry::LocalView;
29 using GridView = typename GridGeometry::GridView;
30 using Element = typename GridView::template Codim<0>::Entity;
31 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
32 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
33
34public:
35
51 template<class FaceVariables>
52 static Scalar velocityGradII(const SubControlVolumeFace& scvf,
53 const FaceVariables& faceVars)
54 {
55 // The velocities of the dof at interest and the one of the opposite scvf.
56 const Scalar velocitySelf = faceVars.velocitySelf();
57 const Scalar velocityOpposite = faceVars.velocityOpposite();
58
59 return ((velocityOpposite - velocitySelf) / scvf.selfToOppositeDistance()) * scvf.directionSign();
60 }
61
105 template<class Problem, class FaceVariables>
106 static Scalar velocityGradIJ(const Problem& problem,
107 const Element& element,
108 const FVElementGeometry& fvGeometry,
109 const SubControlVolumeFace& scvf,
110 const FaceVariables& faceVars,
111 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
112 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
113 const std::size_t localSubFaceIdx)
114 {
115 const auto eIdx = scvf.insideScvIdx();
116 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
117
118 // For the velocityGrad_ij derivative, get the velocities at the current (own) scvf
119 // and at the parallel one at the neighboring scvf.
120 const Scalar innerParallelVelocity = faceVars.velocitySelf();
121
122 const auto outerParallelVelocity = [&]()
123 {
124 if (!lateralScvf.boundary())
125 return faceVars.velocityParallel(localSubFaceIdx, 0);
126 else if (lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())))
127 {
128 // Sample the value of the Dirichlet BC at the center of the staggered lateral face.
129 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
130 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralScvf, lateralBoundaryFacePos);
131 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())];
132 }
133 else if (lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
134 {
135 return beaversJosephVelocityAtLateralScvf(problem, element, fvGeometry, scvf, faceVars,
136 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
137 }
138 else
139 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary type at " << lateralScvf.center());
140 }();
141
142 // The velocity gradient already accounts for the orientation
143 // of the staggered face's outer normal vector. This also correctly accounts for the reduced
144 // distance used in the gradient if the lateral scvf lies on a boundary.
145 return (outerParallelVelocity - innerParallelVelocity)
146 / scvf.parallelDofsDistance(localSubFaceIdx, 0) * lateralScvf.directionSign();
147 }
148
180 template<class Problem, class FaceVariables>
181 static Scalar velocityGradJI(const Problem& problem,
182 const Element& element,
183 const FVElementGeometry& fvGeometry,
184 const SubControlVolumeFace& scvf,
185 const FaceVariables& faceVars,
186 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
187 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
188 const std::size_t localSubFaceIdx)
189 {
190 const auto eIdx = scvf.insideScvIdx();
191 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
192
193 // Assume a zero velocity gradient for pressure boundary conditions.
194 if (currentScvfBoundaryTypes && currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx))
195 return 0.0;
196
197 // For the velocityGrad_ji gradient, get the velocities perpendicular to the velocity at the current scvf.
198 // The inner one is located at staggered face within the own element,
199 // the outer one at the respective staggered face of the element on the other side of the
200 // current scvf.
201 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
202 const Scalar outerLateralVelocity = [&]()
203 {
204 if (!scvf.boundary())
205 return faceVars.velocityLateralOutside(localSubFaceIdx);
206 else if (currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralScvf.directionIndex())))
207 {
208 // Sample the value of the Dirichlet BC at the center of the lateral face intersecting with the boundary.
209 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
210 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralBoundaryFacePos);
211 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralScvf.directionIndex())];
212 }
213 else if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
214 {
215 return beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
216 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
217 }
218 else
219 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary types at " << lateralScvf.center());
220 }();
221
222 // Calculate the velocity gradient in positive coordinate direction.
223 const Scalar lateralDeltaV = scvf.normalInPosCoordDir()
224 ? (outerLateralVelocity - innerLateralVelocity)
225 : (innerLateralVelocity - outerLateralVelocity);
226
227 return lateralDeltaV / scvf.pairData(localSubFaceIdx).lateralDistance;
228 }
229
250 template<class Problem, class FaceVariables>
251 static Scalar beaversJosephVelocityAtCurrentScvf(const Problem& problem,
252 const Element& element,
253 const FVElementGeometry& fvGeometry,
254 const SubControlVolumeFace& scvf,
255 const FaceVariables& faceVars,
256 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
257 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
258 const std::size_t localSubFaceIdx)
259 {
260 const auto eIdx = scvf.insideScvIdx();
261 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
262 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
263
264 const auto tangentialVelocityGradient = [&]()
265 {
266 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
267 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
268 // (towards the current scvf).
269 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
270 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
271
272 if (unsymmetrizedGradientForBJ)
273 return 0.0;
274
275 if (lateralScvf.boundary())
276 {
277 if (lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
278 lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
279 return 0.0;
280 }
281
282 return velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
283 }();
284
285 return problem.beaversJosephVelocity(element,
286 fvGeometry.scv(scvf.insideScvIdx()),
287 lateralScvf,
288 scvf, /*on boundary*/
289 innerLateralVelocity,
290 tangentialVelocityGradient);
291 }
292
310 template<class Problem, class FaceVariables>
311 static Scalar beaversJosephVelocityAtLateralScvf(const Problem& problem,
312 const Element& element,
313 const FVElementGeometry& fvGeometry,
314 const SubControlVolumeFace& scvf,
315 const FaceVariables& faceVars,
316 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
317 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
318 const std::size_t localSubFaceIdx)
319 {
320 const auto eIdx = scvf.insideScvIdx();
321 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
322 const Scalar innerParallelVelocity = faceVars.velocitySelf();
323
324 const auto tangentialVelocityGradient = [&]()
325 {
326 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
327 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
328 // (towards the current scvf).
329 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
330 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
331
332 if (unsymmetrizedGradientForBJ)
333 return 0.0;
334
335 if (scvf.boundary())
336 {
337 if (currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
338 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
339 return 0.0;
340 }
341
342 return velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
343 }();
344
345 return problem.beaversJosephVelocity(element,
346 fvGeometry.scv(scvf.insideScvIdx()),
347 scvf,
348 lateralScvf, /*on boundary*/
349 innerParallelVelocity,
350 tangentialVelocityGradient);
351 }
352
353private:
354
369 static const GlobalPosition& lateralStaggeredFaceCenter_(const SubControlVolumeFace& scvf, const int localSubFaceIdx)
370 {
371 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
372 };
373};
374
375} // end namespace Dumux
376
377#endif
static Scalar velocityGradIJ(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
Returns the velocity gradient perpendicular to the orientation of our current scvf.
Definition: staggered/velocitygradients.hh:106
static Scalar beaversJosephVelocityAtLateralScvf(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
Returns the Beavers-Jospeh slip velocity for a lateral scvf which lies on the boundary.
Definition: staggered/velocitygradients.hh:311
static Scalar velocityGradII(const SubControlVolumeFace &scvf, const FaceVariables &faceVars)
Returns the in-axis velocity gradient.
Definition: staggered/velocitygradients.hh:52
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 Scalar velocityGradJI(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
Returns the velocity gradient in line with our current scvf.
Definition: staggered/velocitygradients.hh:181
static Scalar beaversJosephVelocityAtCurrentScvf(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
Returns the Beavers-Jospeh slip velocity for a scvf which lies on the boundary itself.
Definition: staggered/velocitygradients.hh:251
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
Some exceptions thrown in DuMux
SubControlVolumeFace makeStaggeredBoundaryFace(const SubControlVolumeFace &scvf, const typename SubControlVolumeFace::GlobalPosition &newCenter)
Helper function to turn a given cell scvface into a fake boundary face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:62
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.