24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
37template<
class Scalar,
class Gr
idGeometry,
class BoundaryTypes,
class Indices>
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;
63 template<
class FaceVariables>
65 const FaceVariables& faceVars)
68 const Scalar velocitySelf = faceVars.velocitySelf();
69 const Scalar velocityOpposite = faceVars.velocityOpposite();
71 return ((velocityOpposite - velocitySelf) / scvf.selfToOppositeDistance()) * scvf.directionSign();
97 template<
class Problem,
class FaceVariables>
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)
107 const auto eIdx = scvf.insideScvIdx();
108 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
112 const Scalar innerParallelVelocity = faceVars.velocitySelf();
114 const auto outerParallelVelocity = [&]()
116 if (!lateralScvf.boundary())
117 return faceVars.velocityParallel(localSubFaceIdx, 0);
118 else if (lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())))
121 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
122 return problem.dirichlet(element, lateralScvf.makeBoundaryFace(lateralBoundaryFacePos))[Indices::velocity(scvf.directionIndex())];
124 else if (lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
127 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
130 DUNE_THROW(Dune::InvalidStateException,
"Invalid lateral boundary type at " << lateralScvf.center());
136 return (outerParallelVelocity - innerParallelVelocity)
137 / scvf.parallelDofsDistance(localSubFaceIdx, 0) * lateralScvf.directionSign();
171 template<
class Problem,
class FaceVariables>
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)
181 const auto eIdx = scvf.insideScvIdx();
182 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
185 if (currentScvfBoundaryTypes && currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx))
192 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
193 const Scalar outerLateralVelocity = [&]()
195 if (!scvf.boundary())
196 return faceVars.velocityLateralOutside(localSubFaceIdx);
197 else if (currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralScvf.directionIndex())))
200 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
201 return problem.dirichlet(element, scvf.makeBoundaryFace(lateralBoundaryFacePos))[Indices::velocity(lateralScvf.directionIndex())];
203 else if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
206 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
209 DUNE_THROW(Dune::InvalidStateException,
"Invalid lateral boundary types at " << lateralScvf.center());
213 const Scalar lateralDeltaV = scvf.normalInPosCoordDir()
214 ? (outerLateralVelocity - innerLateralVelocity)
215 : (innerLateralVelocity - outerLateralVelocity);
217 return lateralDeltaV / scvf.pairData(localSubFaceIdx).lateralDistance;
240 template<
class Problem,
class FaceVariables>
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)
250 const auto eIdx = scvf.insideScvIdx();
251 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
252 const Scalar innerLateralVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
254 const auto tangentialVelocityGradient = [&]()
259 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
260 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph",
false);
262 if (unsymmetrizedGradientForBJ)
265 if (lateralScvf.boundary())
267 if (lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
268 lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())))
272 return velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
275 return problem.beaversJosephVelocity(element,
276 fvGeometry.scv(scvf.insideScvIdx()),
279 innerLateralVelocity,
280 tangentialVelocityGradient);
300 template<
class Problem,
class FaceVariables>
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)
310 const auto eIdx = scvf.insideScvIdx();
311 const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
312 const Scalar innerParallelVelocity = faceVars.velocitySelf();
314 const auto tangentialVelocityGradient = [&]()
319 static const bool unsymmetrizedGradientForBJ = getParamFromGroup<bool>(problem.paramGroup(),
320 "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph",
false);
322 if (unsymmetrizedGradientForBJ)
327 if (currentScvfBoundaryTypes->isDirichlet(Indices::pressureIdx) ||
328 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())))
332 return velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
335 return problem.beaversJosephVelocity(element,
336 fvGeometry.scv(scvf.insideScvIdx()),
339 innerParallelVelocity,
340 tangentialVelocityGradient);
359 static const GlobalPosition& lateralStaggeredFaceCenter_(
const SubControlVolumeFace& scvf,
const int localSubFaceIdx)
361 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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