12#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
13#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
18#include <dune/common/fmatrix.hh>
25template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
37template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
49 template<
class FVElementGeometry,
class ElemVolVars>
51 const typename FVElementGeometry::SubControlVolumeFace& scvf,
52 const ElemVolVars& elemVolVars,
53 bool fullGradient =
false)
55 const auto& problem = elemVolVars.gridVolVars().problem();
58 ElementBoundaryTypes elemBcTypes;
59 elemBcTypes.
update(problem, fvGeometry.element(), fvGeometry);
60 return velocityGradient(fvGeometry, scvf, elemVolVars, elemBcTypes, fullGradient);
63 template<
class FVElementGeometry,
class ElemVolVars,
class BoundaryTypes>
65 const typename FVElementGeometry::SubControlVolumeFace& scvf,
66 const ElemVolVars& elemVolVars,
68 bool fullGradient =
false)
70 using Scalar =
typename FVElementGeometry::GridGeometry::GlobalCoordinate::value_type;
71 static constexpr auto dim = FVElementGeometry::GridGeometry::GlobalCoordinate::dimension;
72 Dune::FieldMatrix<Scalar, dim, dim> gradient(0.0);
73 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
75 if (scvf.isFrontal() && !scvf.boundary())
77 gradient[scv.dofAxis()][scv.dofAxis()] =
velocityGradII(fvGeometry, scvf, elemVolVars);
81 Dune::FieldMatrix<std::size_t, dim, dim> counter(0.0);
83 for (
const auto& otherScvf : scvfs(fvGeometry))
85 if (otherScvf.index() == scvf.index())
88 const auto& otherScv = fvGeometry.scv(otherScvf.insideScvIdx());
90 if (otherScvf.isFrontal() && !otherScvf.boundary())
91 gradient[otherScv.dofAxis()][otherScv.dofAxis()] =
velocityGradII(fvGeometry, otherScvf, elemVolVars);
92 else if (otherScvf.isLateral())
94 assert(otherScv.dofAxis() != otherScvf.normalAxis());
95 const auto i = otherScv.dofAxis();
96 const auto j = otherScvf.normalAxis();
98 if (!fvGeometry.hasBoundaryScvf() || !elemBcTypes[otherScvf.localIndex()].isNeumann(i))
100 gradient[i][j] +=
velocityGradIJ(fvGeometry, otherScvf, elemVolVars);
103 if (!fvGeometry.hasBoundaryScvf() || !otherScv.boundary() ||
104 !elemBcTypes[fvGeometry.frontalScvfOnBoundary(otherScv).localIndex()].isNeumann(j))
106 gradient[j][i] +=
velocityGradJI(fvGeometry, otherScvf, elemVolVars);
112 for (
int i = 0; i < dim; ++i)
113 for (
int j = 0; j < dim; ++j)
115 gradient[i][j] /= counter[i][j];
121 DUNE_THROW(Dune::NotImplemented,
"Full gradient for lateral staggered faces not implemented");
123 gradient[scv.dofAxis()][scvf.normalAxis()] =
velocityGradIJ(fvGeometry, scvf, elemVolVars);
124 gradient[scvf.normalAxis()][scv.dofAxis()] =
velocityGradJI(fvGeometry, scvf, elemVolVars);
145 template<
class FVElementGeometry,
class ElemVolVars>
147 const typename FVElementGeometry::SubControlVolumeFace& scvf,
148 const ElemVolVars& elemVolVars)
150 assert(scvf.isFrontal());
152 const auto velocitySelf = elemVolVars[scvf.insideScvIdx()].velocity();
153 const auto velocityOpposite = elemVolVars[scvf.outsideScvIdx()].velocity();
154 const auto distance = (fvGeometry.scv(scvf.outsideScvIdx()).dofPosition() - fvGeometry.scv(scvf.insideScvIdx()).dofPosition()).two_norm();
156 return (velocityOpposite - velocitySelf) /
distance * scvf.directionSign();
180 template<
class FVElementGeometry,
class ElemVolVars>
182 const typename FVElementGeometry::SubControlVolumeFace& scvf,
183 const ElemVolVars& elemVolVars)
185 assert(scvf.isLateral());
186 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
187 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
188 const auto distance = getDistanceIJ_(fvGeometry, scvf);
189 return (outerVelocity - innerVelocity) /
distance * scvf.directionSign();
221 template<
class FVElementGeometry,
class ElemVolVars>
223 const typename FVElementGeometry::SubControlVolumeFace& scvf,
224 const ElemVolVars& elemVolVars)
226 assert(scvf.isLateral());
227 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
228 const auto innerVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
229 const auto outerVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
230 const auto distance = getDistanceJI_(fvGeometry, scvf, orthogonalScvf);
231 return (outerVelocity - innerVelocity) /
distance * orthogonalScvf.directionSign();
236 template<
class FVElementGeometry>
237 static auto getDistanceIJ_(
const FVElementGeometry& fvGeometry,
238 const typename FVElementGeometry::SubControlVolumeFace& scvf)
240 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
244 return getDistanceToBoundary_(insideScv, scvf);
246 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
250 return getStandardDistance_(insideScv, outsideScv);
253 const auto& gridGeometry = fvGeometry.gridGeometry();
254 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
255 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
258 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex()))
259 return getStandardDistance_(insideScv, outsideScv);
267 auto periodicFvGeometry =
localView(gridGeometry);
268 const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex());
269 periodicFvGeometry.bindElement(periodicElement);
271 for (
const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv))
273 if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign())
275 const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
276 const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm();
277 return insideDistance + outsideDistance;
280 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
285 template<
class FVElementGeometry>
286 static auto getDistanceJI_(
const FVElementGeometry& fvGeometry,
287 const typename FVElementGeometry::SubControlVolumeFace& scvf,
288 const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf)
290 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
294 if (orthogonalScvf.boundary())
295 return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf);
297 const auto& orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx());
300 if constexpr (!Detail::SupportsPeriodicity<typename FVElementGeometry::GridGeometry::Grid>())
301 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
304 const auto& gridGeometry = fvGeometry.gridGeometry();
305 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
308 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex()))
309 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
317 auto periodicFvGeometry =
localView(gridGeometry);
318 const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex());
319 periodicFvGeometry.bindElement(periodicElement);
321 for (
const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv))
323 if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign())
325 const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm();
326 const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm();
327 return insideDistance + outsideDistance;
330 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
334 template<
class SubControlVolume>
335 static auto getStandardDistance_(
const SubControlVolume& insideScv,
336 const SubControlVolume& outsideScv)
338 return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm();
341 template<
class SubControlVolume,
class SubControlVolumeFace>
342 static auto getDistanceToBoundary_(
const SubControlVolume& scv,
343 const SubControlVolumeFace& scvf)
345 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:46
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 auto velocityGradient(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars, bool fullGradient=false)
Definition: momentum/velocitygradients.hh:50
static auto velocityGradII(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the in-axis velocity gradient.
Definition: momentum/velocitygradients.hh:146
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
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:64
Definition: momentum/velocitygradients.hh:26
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
Definition: common/pdesolver.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: momentum/velocitygradients.hh:35
Type traits for problem classes.
Definition: common/typetraits/problem.hh:32