3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
26
27#include <optional>
30
31namespace Dumux {
32
37template<class Scalar, class GridGeometry, class BoundaryTypes, class Indices>
39{
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using GridView = typename GridGeometry::GridView;
42 using Element = typename GridView::template Codim<0>::Entity;
43 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
44 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
45
46public:
47
63 template<class FaceVariables>
64 static Scalar velocityGradII(const SubControlVolumeFace& scvf,
65 const FaceVariables& faceVars)
66 {
67 // The velocities of the dof at interest and the one of the opposite scvf.
68 const Scalar velocitySelf = faceVars.velocitySelf();
69 const Scalar velocityOpposite = faceVars.velocityOpposite();
70
71 return ((velocityOpposite - velocitySelf) / scvf.selfToOppositeDistance()) * scvf.directionSign();
72 }
73
117 template<class Problem, class FaceVariables>
118 static Scalar velocityGradIJ(const Problem& problem,
119 const Element& element,
120 const FVElementGeometry& fvGeometry,
121 const SubControlVolumeFace& scvf,
122 const FaceVariables& faceVars,
123 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
124 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
125 const std::size_t localSubFaceIdx)
126 {
127 const auto eIdx = scvf.insideScvIdx();
128 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
129
130 // For the velocityGrad_ij derivative, get the velocities at the current (own) scvf
131 // and at the parallel one at the neighboring scvf.
132 const Scalar innerParallelVelocity = faceVars.velocitySelf();
133
134 const auto outerParallelVelocity = [&]()
135 {
136 if (!lateralScvf.boundary())
137 return faceVars.velocityParallel(localSubFaceIdx, 0);
138 else if (lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())))
139 {
140 // Sample the value of the Dirichlet BC at the center of the staggered lateral face.
141 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
142 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralScvf, lateralBoundaryFacePos);
143 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())];
144 }
145 else if (lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
146 {
147 return beaversJosephVelocityAtLateralScvf(problem, element, fvGeometry, scvf, faceVars,
148 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
149 }
150 else
151 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary type at " << lateralScvf.center());
152 }();
153
154 // The velocity gradient already accounts for the orientation
155 // of the staggered face's outer normal vector. This also correctly accounts for the reduced
156 // distance used in the gradient if the lateral scvf lies on a boundary.
157 return (outerParallelVelocity - innerParallelVelocity)
158 / scvf.parallelDofsDistance(localSubFaceIdx, 0) * lateralScvf.directionSign();
159 }
160
192 template<class Problem, class FaceVariables>
193 static Scalar velocityGradJI(const Problem& problem,
194 const Element& element,
195 const FVElementGeometry& fvGeometry,
196 const SubControlVolumeFace& scvf,
197 const FaceVariables& faceVars,
198 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
199 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
200 const std::size_t localSubFaceIdx)
201 {
202 const auto eIdx = scvf.insideScvIdx();
203 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
204
205 // Assume a zero velocity gradient for pressure boundary conditions.
206 if (currentScvfBoundaryTypes && currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx))
207 return 0.0;
208
209 // For the velocityGrad_ji gradient, get the velocities perpendicular to the velocity at the current scvf.
210 // The inner one is located at staggered face within the own element,
211 // the outer one at the respective staggered face of the element on the other side of the
212 // current scvf.
213 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
214 const Scalar outerLateralVelocity = [&]()
215 {
216 if (!scvf.boundary())
217 return faceVars.velocityLateralOutside(localSubFaceIdx);
218 else if (currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralScvf.directionIndex())))
219 {
220 // Sample the value of the Dirichlet BC at the center of the lateral face intersecting with the boundary.
221 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
222 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralBoundaryFacePos);
223 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralScvf.directionIndex())];
224 }
225 else if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
226 {
227 return beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
228 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
229 }
230 else
231 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary types at " << lateralScvf.center());
232 }();
233
234 // Calculate the velocity gradient in positive coordinate direction.
235 const Scalar lateralDeltaV = scvf.normalInPosCoordDir()
236 ? (outerLateralVelocity - innerLateralVelocity)
237 : (innerLateralVelocity - outerLateralVelocity);
238
239 return lateralDeltaV / scvf.pairData(localSubFaceIdx).lateralDistance;
240 }
241
262 template<class Problem, class FaceVariables>
263 static Scalar beaversJosephVelocityAtCurrentScvf(const Problem& problem,
264 const Element& element,
265 const FVElementGeometry& fvGeometry,
266 const SubControlVolumeFace& scvf,
267 const FaceVariables& faceVars,
268 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
269 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
270 const std::size_t localSubFaceIdx)
271 {
272 const auto eIdx = scvf.insideScvIdx();
273 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
274 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
275
276 const auto tangentialVelocityGradient = [&]()
277 {
278 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
279 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
280 // (towards the current scvf).
281 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
282 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
283
284 if (unsymmetrizedGradientForBJ)
285 return 0.0;
286
287 if (lateralScvf.boundary())
288 {
289 if (lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
290 lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
291 return 0.0;
292 }
293
294 return velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
295 }();
296
297 return problem.beaversJosephVelocity(element,
298 fvGeometry.scv(scvf.insideScvIdx()),
299 lateralScvf,
300 scvf, /*on boundary*/
301 innerLateralVelocity,
302 tangentialVelocityGradient);
303 }
304
322 template<class Problem, class FaceVariables>
323 static Scalar beaversJosephVelocityAtLateralScvf(const Problem& problem,
324 const Element& element,
325 const FVElementGeometry& fvGeometry,
326 const SubControlVolumeFace& scvf,
327 const FaceVariables& faceVars,
328 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
329 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
330 const std::size_t localSubFaceIdx)
331 {
332 const auto eIdx = scvf.insideScvIdx();
333 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
334 const Scalar innerParallelVelocity = faceVars.velocitySelf();
335
336 const auto tangentialVelocityGradient = [&]()
337 {
338 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
339 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
340 // (towards the current scvf).
341 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
342 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
343
344 if (unsymmetrizedGradientForBJ)
345 return 0.0;
346
347 if (scvf.boundary())
348 {
349 if (currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
350 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
351 return 0.0;
352 }
353
354 return velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
355 }();
356
357 return problem.beaversJosephVelocity(element,
358 fvGeometry.scv(scvf.insideScvIdx()),
359 scvf,
360 lateralScvf, /*on boundary*/
361 innerParallelVelocity,
362 tangentialVelocityGradient);
363 }
364
365private:
366
381 static const GlobalPosition& lateralStaggeredFaceCenter_(const SubControlVolumeFace& scvf, const int localSubFaceIdx)
382 {
383 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
384 };
385};
386
387} // end namespace Dumux
388
389#endif
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:104
Definition: adapt.hh:29
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: velocitygradients.hh:39
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: velocitygradients.hh:118
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: velocitygradients.hh:323
static Scalar velocityGradII(const SubControlVolumeFace &scvf, const FaceVariables &faceVars)
Returns the in-axis velocity gradient.
Definition: velocitygradients.hh:64
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: velocitygradients.hh:193
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: velocitygradients.hh:263