12#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
13#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
18#include <dune/common/fmatrix.hh>
33 template<
class FVElementGeometry,
class ElemVolVars>
35 const typename FVElementGeometry::SubControlVolumeFace& scvf,
36 const ElemVolVars& elemVolVars,
37 bool fullGradient =
false)
39 const auto& problem = elemVolVars.gridVolVars().problem();
42 ElementBoundaryTypes elemBcTypes;
43 elemBcTypes.
update(problem, fvGeometry.element(), fvGeometry);
44 return velocityGradient(fvGeometry, scvf, elemVolVars, elemBcTypes, fullGradient);
47 template<
class FVElementGeometry,
class ElemVolVars,
class BoundaryTypes>
49 const typename FVElementGeometry::SubControlVolumeFace& scvf,
50 const ElemVolVars& elemVolVars,
52 bool fullGradient =
false)
54 using Scalar =
typename FVElementGeometry::GridGeometry::GlobalCoordinate::value_type;
55 static constexpr auto dim = FVElementGeometry::GridGeometry::GlobalCoordinate::dimension;
56 Dune::FieldMatrix<Scalar, dim, dim> gradient(0.0);
57 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
59 if (scvf.isFrontal() && !scvf.boundary())
61 gradient[scv.dofAxis()][scv.dofAxis()] =
velocityGradII(fvGeometry, scvf, elemVolVars);
65 Dune::FieldMatrix<std::size_t, dim, dim> counter(0.0);
67 for (
const auto& otherScvf : scvfs(fvGeometry))
69 if (otherScvf.index() == scvf.index())
72 const auto& otherScv = fvGeometry.scv(otherScvf.insideScvIdx());
74 if (otherScvf.isFrontal() && !otherScvf.boundary())
75 gradient[otherScv.dofAxis()][otherScv.dofAxis()] =
velocityGradII(fvGeometry, otherScvf, elemVolVars);
76 else if (otherScvf.isLateral())
78 assert(otherScv.dofAxis() != otherScvf.normalAxis());
79 const auto i = otherScv.dofAxis();
80 const auto j = otherScvf.normalAxis();
82 if (!fvGeometry.hasBoundaryScvf() || !elemBcTypes[otherScvf.localIndex()].isNeumann(i))
84 gradient[i][j] +=
velocityGradIJ(fvGeometry, otherScvf, elemVolVars);
87 if (!fvGeometry.hasBoundaryScvf() || !otherScv.boundary() ||
88 !elemBcTypes[fvGeometry.frontalScvfOnBoundary(otherScv).localIndex()].isNeumann(j))
90 gradient[j][i] +=
velocityGradJI(fvGeometry, otherScvf, elemVolVars);
96 for (
int i = 0; i < dim; ++i)
97 for (
int j = 0; j < dim; ++j)
99 gradient[i][j] /= counter[i][j];
105 DUNE_THROW(Dune::NotImplemented,
"Full gradient for lateral staggered faces not implemented");
107 gradient[scv.dofAxis()][scvf.normalAxis()] =
velocityGradIJ(fvGeometry, scvf, elemVolVars);
108 gradient[scvf.normalAxis()][scv.dofAxis()] =
velocityGradJI(fvGeometry, scvf, elemVolVars);
129 template<
class FVElementGeometry,
class ElemVolVars>
131 const typename FVElementGeometry::SubControlVolumeFace& scvf,
132 const ElemVolVars& elemVolVars)
134 assert(scvf.isFrontal());
136 const auto velocitySelf = elemVolVars[scvf.insideScvIdx()].velocity();
137 const auto velocityOpposite = elemVolVars[scvf.outsideScvIdx()].velocity();
138 const auto distance = (fvGeometry.scv(scvf.outsideScvIdx()).dofPosition() - fvGeometry.scv(scvf.insideScvIdx()).dofPosition()).two_norm();
140 return (velocityOpposite - velocitySelf) /
distance * scvf.directionSign();
164 template<
class FVElementGeometry,
class ElemVolVars>
166 const typename FVElementGeometry::SubControlVolumeFace& scvf,
167 const ElemVolVars& elemVolVars)
169 assert(scvf.isLateral());
170 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
171 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
172 const auto distance = getDistanceIJ_(fvGeometry, scvf);
173 return (outerVelocity - innerVelocity) /
distance * scvf.directionSign();
205 template<
class FVElementGeometry,
class ElemVolVars>
207 const typename FVElementGeometry::SubControlVolumeFace& scvf,
208 const ElemVolVars& elemVolVars)
210 assert(scvf.isLateral());
211 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
212 const auto innerVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
213 const auto outerVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
214 const auto distance = getDistanceJI_(fvGeometry, scvf, orthogonalScvf);
215 return (outerVelocity - innerVelocity) /
distance * orthogonalScvf.directionSign();
220 template<
class FVElementGeometry>
221 static auto getDistanceIJ_(
const FVElementGeometry& fvGeometry,
222 const typename FVElementGeometry::SubControlVolumeFace& scvf)
224 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
228 return getDistanceToBoundary_(insideScv, scvf);
230 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
233 if constexpr (!supportsPeriodicity<typename FVElementGeometry::GridGeometry>())
234 return getStandardDistance_(insideScv, outsideScv);
237 const auto& gridGeometry = fvGeometry.gridGeometry();
238 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
239 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
242 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex()))
243 return getStandardDistance_(insideScv, outsideScv);
251 auto periodicFvGeometry =
localView(gridGeometry);
252 const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex());
253 periodicFvGeometry.bindElement(periodicElement);
255 for (
const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv))
257 if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign())
259 const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
260 const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm();
261 return insideDistance + outsideDistance;
264 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
269 template<
class FVElementGeometry>
270 static auto getDistanceJI_(
const FVElementGeometry& fvGeometry,
271 const typename FVElementGeometry::SubControlVolumeFace& scvf,
272 const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf)
274 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
278 if (orthogonalScvf.boundary())
279 return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf);
281 const auto& orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx());
284 if constexpr (!supportsPeriodicity<typename FVElementGeometry::GridGeometry>())
285 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
288 const auto& gridGeometry = fvGeometry.gridGeometry();
289 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
292 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex()))
293 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
301 auto periodicFvGeometry =
localView(gridGeometry);
302 const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex());
303 periodicFvGeometry.bindElement(periodicElement);
305 for (
const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv))
307 if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign())
309 const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm();
310 const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm();
311 return insideDistance + outsideDistance;
314 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
318 template<
class SubControlVolume>
319 static auto getStandardDistance_(
const SubControlVolume& insideScv,
320 const SubControlVolume& outsideScv)
322 return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm();
325 template<
class SubControlVolume,
class SubControlVolumeFace>
326 static auto getDistanceToBoundary_(
const SubControlVolume& scv,
327 const SubControlVolumeFace& scvf)
329 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:26
This class stores an array of BoundaryTypes objects.
Definition: facecentered/staggered/elementboundarytypes.hh:25
void update(const Problem &problem, const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: facecentered/staggered/elementboundarytypes.hh:38
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:30
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:206
static auto velocityGradient(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars, bool fullGradient=false)
Definition: momentum/velocitygradients.hh:34
static auto velocityGradII(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the in-axis velocity gradient.
Definition: momentum/velocitygradients.hh:130
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:165
static auto velocityGradient(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars, const FaceCenteredStaggeredElementBoundaryTypes< BoundaryTypes > &elemBcTypes, bool fullGradient=false)
Definition: momentum/velocitygradients.hh:48
Type traits for problem classes.
Some exceptions thrown in DuMux
Boundary types gathered on an element.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Grid properties related to periodicity.
Type traits for problem classes.
Definition: common/typetraits/problem.hh:32