version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
momentum/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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
13#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
14
15#include <optional>
18#include <dune/common/fmatrix.hh>
22
23namespace Dumux {
24
30{
31public:
32
33 template<class FVElementGeometry, class ElemVolVars>
34 static auto velocityGradient(const FVElementGeometry& fvGeometry,
35 const typename FVElementGeometry::SubControlVolumeFace& scvf,
36 const ElemVolVars& elemVolVars,
37 bool fullGradient = false)
38 {
39 const auto& problem = elemVolVars.gridVolVars().problem();
40 using BoundaryTypes = typename ProblemTraits<std::decay_t<decltype(problem)>>::BoundaryTypes;
42 ElementBoundaryTypes elemBcTypes;
43 elemBcTypes.update(problem, fvGeometry.element(), fvGeometry);
44 return velocityGradient(fvGeometry, scvf, elemVolVars, elemBcTypes, fullGradient);
45 }
46
47 template<class FVElementGeometry, class ElemVolVars, class BoundaryTypes>
48 static auto velocityGradient(const FVElementGeometry& fvGeometry,
49 const typename FVElementGeometry::SubControlVolumeFace& scvf,
50 const ElemVolVars& elemVolVars,
52 bool fullGradient = false)
53 {
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());
58
59 if (scvf.isFrontal() && !scvf.boundary())
60 {
61 gradient[scv.dofAxis()][scv.dofAxis()] = velocityGradII(fvGeometry, scvf, elemVolVars);
62
63 if (fullGradient)
64 {
65 Dune::FieldMatrix<std::size_t, dim, dim> counter(0.0);
66
67 for (const auto& otherScvf : scvfs(fvGeometry))
68 {
69 if (otherScvf.index() == scvf.index())
70 continue;
71
72 const auto& otherScv = fvGeometry.scv(otherScvf.insideScvIdx());
73
74 if (otherScvf.isFrontal() && !otherScvf.boundary())
75 gradient[otherScv.dofAxis()][otherScv.dofAxis()] = velocityGradII(fvGeometry, otherScvf, elemVolVars);
76 else if (otherScvf.isLateral())
77 {
78 assert(otherScv.dofAxis() != otherScvf.normalAxis());
79 const auto i = otherScv.dofAxis();
80 const auto j = otherScvf.normalAxis();
81
82 if (!fvGeometry.hasBoundaryScvf() || !elemBcTypes[otherScvf.localIndex()].isNeumann(i))
83 {
84 gradient[i][j] += velocityGradIJ(fvGeometry, otherScvf, elemVolVars);
85 ++counter[i][j];
86 }
87 if (!fvGeometry.hasBoundaryScvf() || !otherScv.boundary() ||
88 !elemBcTypes[fvGeometry.frontalScvfOnBoundary(otherScv).localIndex()].isNeumann(j))
89 {
90 gradient[j][i] += velocityGradJI(fvGeometry, otherScvf, elemVolVars);
91 ++counter[j][i];
92 }
93 }
94 }
95
96 for (int i = 0; i < dim; ++i)
97 for (int j = 0; j < dim; ++j)
98 if (i != j)
99 gradient[i][j] /= counter[i][j];
100 }
101 }
102 else
103 {
104 if (fullGradient)
105 DUNE_THROW(Dune::NotImplemented, "Full gradient for lateral staggered faces not implemented");
106
107 gradient[scv.dofAxis()][scvf.normalAxis()] = velocityGradIJ(fvGeometry, scvf, elemVolVars);
108 gradient[scvf.normalAxis()][scv.dofAxis()] = velocityGradJI(fvGeometry, scvf, elemVolVars);
109 }
110
111 return gradient;
112 }
113
129 template<class FVElementGeometry, class ElemVolVars>
130 static auto velocityGradII(const FVElementGeometry& fvGeometry,
131 const typename FVElementGeometry::SubControlVolumeFace& scvf,
132 const ElemVolVars& elemVolVars)
133 {
134 assert(scvf.isFrontal());
135 // The velocities of the dof at interest and the one of the opposite scvf.
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();
139
140 return (velocityOpposite - velocitySelf) / distance * scvf.directionSign();
141 }
142
164 template<class FVElementGeometry, class ElemVolVars>
165 static auto velocityGradIJ(const FVElementGeometry& fvGeometry,
166 const typename FVElementGeometry::SubControlVolumeFace& scvf,
167 const ElemVolVars& elemVolVars)
168 {
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();
174 }
175
205 template<class FVElementGeometry, class ElemVolVars>
206 static auto velocityGradJI(const FVElementGeometry& fvGeometry,
207 const typename FVElementGeometry::SubControlVolumeFace& scvf,
208 const ElemVolVars& elemVolVars)
209 {
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();
216 }
217
218private:
219
220 template<class FVElementGeometry>
221 static auto getDistanceIJ_(const FVElementGeometry& fvGeometry,
222 const typename FVElementGeometry::SubControlVolumeFace& scvf)
223 {
224 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
225 // The scvf is on a boundary, hence there is no outer DOF.
226 // We take the distance to the boundary instead.
227 if (scvf.boundary())
228 return getDistanceToBoundary_(insideScv, scvf);
229
230 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
231
232 // The standard case: Our grid does not support periodic boundaries.
233 if constexpr (!supportsPeriodicity<typename FVElementGeometry::GridGeometry>())
234 return getStandardDistance_(insideScv, outsideScv);
235 else
236 {
237 const auto& gridGeometry = fvGeometry.gridGeometry();
238 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
239 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
240
241 // The standard case: Our grid is not periodic or the lateral scvf does not lie on a periodic boundary.
242 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex()))
243 return getStandardDistance_(insideScv, outsideScv);
244
245 // Treating periodic boundaries is more involved:
246 // 1. Consider the outside scv within the element adjacent to the other periodic boundary.
247 // 2. Iterate over this scv's faces until you find the face parallel to our own scvf.
248 // This face would lie directly next to our own scvf if the grid was not periodic.
249 // 3. Calculate the total distance by summing up the distances between the DOFs and the respective integration points
250 // corresponding to the inside and outside scvfs.
251 auto periodicFvGeometry = localView(gridGeometry);
252 const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex());
253 periodicFvGeometry.bindElement(periodicElement);
254
255 for (const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv))
256 {
257 if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign())
258 {
259 const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
260 const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm();
261 return insideDistance + outsideDistance;
262 }
263 }
264 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element");
265 }
266 }
267
269 template<class FVElementGeometry>
270 static auto getDistanceJI_(const FVElementGeometry& fvGeometry,
271 const typename FVElementGeometry::SubControlVolumeFace& scvf,
272 const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf)
273 {
274 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
275
276 // The orthogonal scvf is on a boundary, hence there is no outer DOF.
277 // We take the distance to the boundary instead.
278 if (orthogonalScvf.boundary())
279 return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf);
280
281 const auto& orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx());
282
283 // The standard case: Our grid does not support periodic boundaries.
284 if constexpr (!supportsPeriodicity<typename FVElementGeometry::GridGeometry>())
285 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
286 else
287 {
288 const auto& gridGeometry = fvGeometry.gridGeometry();
289 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
290
291 // The standard case: Our grid is not periodic or our own DOF does not lie on a periodic boundary.
292 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex()))
293 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
294
295 // Treating periodic boundaries is more involved:
296 // 1. Consider the orthogonal outside scv within the element adjacent to the other periodic boundary.
297 // 2. Iterate over this scv's faces until you find the orthogonal face parallel to our own orthogonal scvf.
298 // This face would lie directly next to our own orthogonal scvf if the grid was not periodic.
299 // 3. Calculate the total distance by summing up the distances between the DOFs and the respective integration points
300 // corresponding to the orthogonal scvfs.
301 auto periodicFvGeometry = localView(gridGeometry);
302 const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex());
303 periodicFvGeometry.bindElement(periodicElement);
304
305 for (const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv))
306 {
307 if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign())
308 {
309 const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm();
310 const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm();
311 return insideDistance + outsideDistance;
312 }
313 }
314 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element");
315 }
316 }
317
318 template<class SubControlVolume>
319 static auto getStandardDistance_(const SubControlVolume& insideScv,
320 const SubControlVolume& outsideScv)
321 {
322 return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm();
323 }
324
325 template<class SubControlVolume, class SubControlVolumeFace>
326 static auto getDistanceToBoundary_(const SubControlVolume& scv,
327 const SubControlVolumeFace& scvf)
328 {
329 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
330 }
331};
332
333} // end namespace Dumux
334
335#endif // DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
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
Definition: adapt.hh:17
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