3.2-git
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
97 template<class Problem, class FaceVariables>
98 static Scalar velocityGradIJ(const Problem& problem,
99 const Element& element,
100 const FVElementGeometry& fvGeometry,
101 const SubControlVolumeFace& scvf,
102 const FaceVariables& faceVars,
103 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
104 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
105 const std::size_t localSubFaceIdx)
106 {
107 const auto eIdx = scvf.insideScvIdx();
108 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
109
110 // For the velocityGrad_ij derivative, get the velocities at the current (own) scvf
111 // and at the parallel one at the neighboring scvf.
112 const Scalar innerParallelVelocity = faceVars.velocitySelf();
113
114 const auto outerParallelVelocity = [&]()
115 {
116 if (!lateralScvf.boundary())
117 return faceVars.velocityParallel(localSubFaceIdx, 0);
118 else if (lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())))
119 {
120 // Sample the value of the Dirichlet BC at the center of the staggered lateral face.
121 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
122 return problem.dirichlet(element, lateralScvf.makeBoundaryFace(lateralBoundaryFacePos))[Indices::velocity(scvf.directionIndex())];
123 }
124 else if (lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
125 {
126 return beaversJosephVelocityAtLateralScvf(problem, element, fvGeometry, scvf, faceVars,
127 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
128 }
129 else
130 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary type at " << lateralScvf.center());
131 }();
132
133 // The velocity gradient already accounts for the orientation
134 // of the staggered face's outer normal vector. This also correctly accounts for the reduced
135 // distance used in the gradient if the lateral scvf lies on a boundary.
136 return (outerParallelVelocity - innerParallelVelocity)
137 / scvf.parallelDofsDistance(localSubFaceIdx, 0) * lateralScvf.directionSign();
138 }
139
171 template<class Problem, class FaceVariables>
172 static Scalar velocityGradJI(const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const SubControlVolumeFace& scvf,
176 const FaceVariables& faceVars,
177 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
178 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
179 const std::size_t localSubFaceIdx)
180 {
181 const auto eIdx = scvf.insideScvIdx();
182 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
183
184 // Assume a zero velocity gradient for pressure boundary conditions.
185 if (currentScvfBoundaryTypes && currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx))
186 return 0.0;
187
188 // For the velocityGrad_ji gradient, get the velocities perpendicular to the velocity at the current scvf.
189 // The inner one is located at staggered face within the own element,
190 // the outer one at the respective staggered face of the element on the other side of the
191 // current scvf.
192 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
193 const Scalar outerLateralVelocity = [&]()
194 {
195 if (!scvf.boundary())
196 return faceVars.velocityLateralOutside(localSubFaceIdx);
197 else if (currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralScvf.directionIndex())))
198 {
199 // Sample the value of the Dirichlet BC at the center of the lateral face intersecting with the boundary.
200 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
201 return problem.dirichlet(element, scvf.makeBoundaryFace(lateralBoundaryFacePos))[Indices::velocity(lateralScvf.directionIndex())];
202 }
203 else if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
204 {
205 return beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
206 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
207 }
208 else
209 DUNE_THROW(Dune::InvalidStateException, "Invalid lateral boundary types at " << lateralScvf.center());
210 }();
211
212 // Calculate the velocity gradient in positive coordinate direction.
213 const Scalar lateralDeltaV = scvf.normalInPosCoordDir()
214 ? (outerLateralVelocity - innerLateralVelocity)
215 : (innerLateralVelocity - outerLateralVelocity);
216
217 return lateralDeltaV / scvf.pairData(localSubFaceIdx).lateralDistance;
218 }
219
240 template<class Problem, class FaceVariables>
241 static Scalar beaversJosephVelocityAtCurrentScvf(const Problem& problem,
242 const Element& element,
243 const FVElementGeometry& fvGeometry,
244 const SubControlVolumeFace& scvf,
245 const FaceVariables& faceVars,
246 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
247 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
248 const std::size_t localSubFaceIdx)
249 {
250 const auto eIdx = scvf.insideScvIdx();
251 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
252 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
253
254 const auto tangentialVelocityGradient = [&]()
255 {
256 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
257 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
258 // (towards the current scvf).
259 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
260 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
261
262 if (unsymmetrizedGradientForBJ)
263 return 0.0;
264
265 if (lateralScvf.boundary())
266 {
267 if (lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
268 lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
269 return 0.0;
270 }
271
272 return velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
273 }();
274
275 return problem.beaversJosephVelocity(element,
276 fvGeometry.scv(scvf.insideScvIdx()),
277 lateralScvf,
278 scvf, /*on boundary*/
279 innerLateralVelocity,
280 tangentialVelocityGradient);
281 }
282
300 template<class Problem, class FaceVariables>
301 static Scalar beaversJosephVelocityAtLateralScvf(const Problem& problem,
302 const Element& element,
303 const FVElementGeometry& fvGeometry,
304 const SubControlVolumeFace& scvf,
305 const FaceVariables& faceVars,
306 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
307 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
308 const std::size_t localSubFaceIdx)
309 {
310 const auto eIdx = scvf.insideScvIdx();
311 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
312 const Scalar innerParallelVelocity = faceVars.velocitySelf();
313
314 const auto tangentialVelocityGradient = [&]()
315 {
316 // If the current scvf is on a boundary and if a Dirichlet BC for the pressure or a BJ condition for
317 // the slip velocity is set there, assume a tangential velocity gradient of zero along the lateral face
318 // (towards the current scvf).
319 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
320 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false);
321
322 if (unsymmetrizedGradientForBJ)
323 return 0.0;
324
325 if (scvf.boundary())
326 {
327 if (currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
328 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
329 return 0.0;
330 }
331
332 return velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
333 }();
334
335 return problem.beaversJosephVelocity(element,
336 fvGeometry.scv(scvf.insideScvIdx()),
337 scvf,
338 lateralScvf, /*on boundary*/
339 innerParallelVelocity,
340 tangentialVelocityGradient);
341 }
342
343private:
344
359 static const GlobalPosition& lateralStaggeredFaceCenter_(const SubControlVolumeFace& scvf, const int localSubFaceIdx)
360 {
361 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
362 };
363};
364
365} // end namespace Dumux
366
367#endif // DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:98
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:301
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:172
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:241