24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
30#include <dune/common/fmatrix.hh>
37template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
49template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
61 template<
class FVElementGeometry,
class ElemVolVars>
63 const typename FVElementGeometry::SubControlVolumeFace& scvf,
64 const ElemVolVars& elemVolVars,
65 bool fullGradient =
false)
67 const auto& problem = elemVolVars.gridVolVars().problem();
70 ElementBoundaryTypes elemBcTypes;
71 elemBcTypes.
update(problem, fvGeometry.element(), fvGeometry);
72 return velocityGradient(fvGeometry, scvf, elemVolVars, elemBcTypes, fullGradient);
75 template<
class FVElementGeometry,
class ElemVolVars,
class BoundaryTypes>
77 const typename FVElementGeometry::SubControlVolumeFace& scvf,
78 const ElemVolVars& elemVolVars,
80 bool fullGradient =
false)
82 using Scalar =
typename FVElementGeometry::GridGeometry::GlobalCoordinate::value_type;
83 static constexpr auto dim = FVElementGeometry::GridGeometry::GlobalCoordinate::dimension;
84 Dune::FieldMatrix<Scalar, dim, dim> gradient(0.0);
85 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
87 if (scvf.isFrontal() && !scvf.boundary())
89 gradient[scv.dofAxis()][scv.dofAxis()] =
velocityGradII(fvGeometry, scvf, elemVolVars);
93 Dune::FieldMatrix<std::size_t, dim, dim> counter(0.0);
95 for (
const auto& otherScvf : scvfs(fvGeometry))
97 if (otherScvf.index() == scvf.index())
100 const auto& otherScv = fvGeometry.scv(otherScvf.insideScvIdx());
102 if (otherScvf.isFrontal() && !otherScvf.boundary())
103 gradient[otherScv.dofAxis()][otherScv.dofAxis()] =
velocityGradII(fvGeometry, otherScvf, elemVolVars);
104 else if (otherScvf.isLateral())
106 assert(otherScv.dofAxis() != otherScvf.normalAxis());
107 const auto i = otherScv.dofAxis();
108 const auto j = otherScvf.normalAxis();
110 if (!fvGeometry.hasBoundaryScvf() || !elemBcTypes[otherScvf.localIndex()].isNeumann(i))
112 gradient[i][j] +=
velocityGradIJ(fvGeometry, otherScvf, elemVolVars);
115 if (!fvGeometry.hasBoundaryScvf() || !otherScv.boundary() ||
116 !elemBcTypes[fvGeometry.frontalScvfOnBoundary(otherScv).localIndex()].isNeumann(j))
118 gradient[j][i] +=
velocityGradJI(fvGeometry, otherScvf, elemVolVars);
124 for (
int i = 0; i < dim; ++i)
125 for (
int j = 0; j < dim; ++j)
127 gradient[i][j] /= counter[i][j];
133 DUNE_THROW(Dune::NotImplemented,
"Full gradient for lateral staggered faces not implemented");
135 gradient[scv.dofAxis()][scvf.normalAxis()] =
velocityGradIJ(fvGeometry, scvf, elemVolVars);
136 gradient[scvf.normalAxis()][scv.dofAxis()] =
velocityGradJI(fvGeometry, scvf, elemVolVars);
157 template<
class FVElementGeometry,
class ElemVolVars>
159 const typename FVElementGeometry::SubControlVolumeFace& scvf,
160 const ElemVolVars& elemVolVars)
162 assert(scvf.isFrontal());
164 const auto velocitySelf = elemVolVars[scvf.insideScvIdx()].velocity();
165 const auto velocityOpposite = elemVolVars[scvf.outsideScvIdx()].velocity();
166 const auto distance = (fvGeometry.scv(scvf.outsideScvIdx()).dofPosition() - fvGeometry.scv(scvf.insideScvIdx()).dofPosition()).two_norm();
168 return (velocityOpposite - velocitySelf) /
distance * scvf.directionSign();
192 template<
class FVElementGeometry,
class ElemVolVars>
194 const typename FVElementGeometry::SubControlVolumeFace& scvf,
195 const ElemVolVars& elemVolVars)
197 assert(scvf.isLateral());
198 const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity();
199 const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity();
200 const auto distance = getDistanceIJ_(fvGeometry, scvf);
201 return (outerVelocity - innerVelocity) /
distance * scvf.directionSign();
233 template<
class FVElementGeometry,
class ElemVolVars>
235 const typename FVElementGeometry::SubControlVolumeFace& scvf,
236 const ElemVolVars& elemVolVars)
238 assert(scvf.isLateral());
239 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
240 const auto innerVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
241 const auto outerVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
242 const auto distance = getDistanceJI_(fvGeometry, scvf, orthogonalScvf);
243 return (outerVelocity - innerVelocity) /
distance * orthogonalScvf.directionSign();
248 template<
class FVElementGeometry>
249 static auto getDistanceIJ_(
const FVElementGeometry& fvGeometry,
250 const typename FVElementGeometry::SubControlVolumeFace& scvf)
252 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
256 return getDistanceToBoundary_(insideScv, scvf);
258 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
262 return getStandardDistance_(insideScv, outsideScv);
265 const auto& gridGeometry = fvGeometry.gridGeometry();
266 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
267 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
270 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex()))
271 return getStandardDistance_(insideScv, outsideScv);
279 auto periodicFvGeometry =
localView(gridGeometry);
280 const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex());
281 periodicFvGeometry.bindElement(periodicElement);
283 for (
const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv))
285 if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign())
287 const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
288 const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm();
289 return insideDistance + outsideDistance;
292 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
297 template<
class FVElementGeometry>
298 static auto getDistanceJI_(
const FVElementGeometry& fvGeometry,
299 const typename FVElementGeometry::SubControlVolumeFace& scvf,
300 const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf)
302 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
306 if (orthogonalScvf.boundary())
307 return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf);
309 const auto& orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx());
312 if constexpr (!Detail::SupportsPeriodicity<typename FVElementGeometry::GridGeometry::Grid>())
313 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
316 const auto& gridGeometry = fvGeometry.gridGeometry();
317 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
320 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex()))
321 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
329 auto periodicFvGeometry =
localView(gridGeometry);
330 const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex());
331 periodicFvGeometry.bindElement(periodicElement);
333 for (
const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv))
335 if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign())
337 const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm();
338 const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm();
339 return insideDistance + outsideDistance;
342 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf in periodic element");
346 template<
class SubControlVolume>
347 static auto getStandardDistance_(
const SubControlVolume& insideScv,
348 const SubControlVolume& outsideScv)
350 return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm();
353 template<
class SubControlVolume,
class SubControlVolumeFace>
354 static auto getDistanceToBoundary_(
const SubControlVolume& scv,
355 const SubControlVolumeFace& scvf)
357 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:294
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
Type traits for problem classes.
Definition: common/typetraits/problem.hh:44
This class stores an array of BoundaryTypes objects.
Definition: facecentered/staggered/elementboundarytypes.hh:37
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:50
Definition: momentum/velocitygradients.hh:38
Definition: momentum/velocitygradients.hh:47
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:58
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:234
static auto velocityGradient(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars, bool fullGradient=false)
Definition: momentum/velocitygradients.hh:62
static auto velocityGradII(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the in-axis velocity gradient.
Definition: momentum/velocitygradients.hh:158
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:193
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:76
Type traits for problem classes.
Boundary types gathered on an element.