3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
26
27#include <optional>
30#include <dune/common/fmatrix.hh>
33
34// forward declare
35namespace Dune {
36
37template<class ct, int dim, template< int > class Ref, class Comm>
38class SPGrid;
39}
40
41namespace Dumux {
42
43namespace Detail
44{
45
46template<class Grid>
47struct SupportsPeriodicity : public std::false_type {};
48
49template<class ct, int dim, template< int > class Ref, class Comm>
50struct SupportsPeriodicity<Dune::SPGrid<ct, dim, Ref, Comm>> : public std::true_type {};
51}
52
58{
59public:
60
61 template<class FVElementGeometry, class ElemVolVars>
62 static auto velocityGradient(const FVElementGeometry& fvGeometry,
63 const typename FVElementGeometry::SubControlVolumeFace& scvf,
64 const ElemVolVars& elemVolVars,
65 bool fullGradient = false)
66 {
67 const auto& problem = elemVolVars.gridVolVars().problem();
68 using BoundaryTypes = typename ProblemTraits<std::decay_t<decltype(problem)>>::BoundaryTypes;
70 ElementBoundaryTypes elemBcTypes;
71 elemBcTypes.update(problem, fvGeometry.element(), fvGeometry);
72 return velocityGradient(fvGeometry, scvf, elemVolVars, elemBcTypes, fullGradient);
73 }
74
75 template<class FVElementGeometry, class ElemVolVars, class BoundaryTypes>
76 static auto velocityGradient(const FVElementGeometry& fvGeometry,
77 const typename FVElementGeometry::SubControlVolumeFace& scvf,
78 const ElemVolVars& elemVolVars,
80 bool fullGradient = false)
81 {
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());
86
87 if (scvf.isFrontal() && !scvf.boundary())
88 {
89 gradient[scv.dofAxis()][scv.dofAxis()] = velocityGradII(fvGeometry, scvf, elemVolVars);
90
91 if (fullGradient)
92 {
93 Dune::FieldMatrix<std::size_t, dim, dim> counter(0.0);
94
95 for (const auto& otherScvf : scvfs(fvGeometry))
96 {
97 if (otherScvf.index() == scvf.index())
98 continue;
99
100 const auto& otherScv = fvGeometry.scv(otherScvf.insideScvIdx());
101
102 if (otherScvf.isFrontal() && !otherScvf.boundary())
103 gradient[otherScv.dofAxis()][otherScv.dofAxis()] = velocityGradII(fvGeometry, otherScvf, elemVolVars);
104 else if (otherScvf.isLateral())
105 {
106 assert(otherScv.dofAxis() != otherScvf.normalAxis());
107 const auto i = otherScv.dofAxis();
108 const auto j = otherScvf.normalAxis();
109
110 if (!fvGeometry.hasBoundaryScvf() || !elemBcTypes[otherScvf.localIndex()].isNeumann(i))
111 {
112 gradient[i][j] += velocityGradIJ(fvGeometry, otherScvf, elemVolVars);
113 ++counter[i][j];
114 }
115 if (!fvGeometry.hasBoundaryScvf() || !otherScv.boundary() ||
116 !elemBcTypes[fvGeometry.frontalScvfOnBoundary(otherScv).localIndex()].isNeumann(j))
117 {
118 gradient[j][i] += velocityGradJI(fvGeometry, otherScvf, elemVolVars);
119 ++counter[j][i];
120 }
121 }
122 }
123
124 for (int i = 0; i < dim; ++i)
125 for (int j = 0; j < dim; ++j)
126 if (i != j)
127 gradient[i][j] /= counter[i][j];
128 }
129 }
130 else
131 {
132 if (fullGradient)
133 DUNE_THROW(Dune::NotImplemented, "Full gradient for lateral staggered faces not implemented");
134
135 gradient[scv.dofAxis()][scvf.normalAxis()] = velocityGradIJ(fvGeometry, scvf, elemVolVars);
136 gradient[scvf.normalAxis()][scv.dofAxis()] = velocityGradJI(fvGeometry, scvf, elemVolVars);
137 }
138
139 return gradient;
140 }
141
157 template<class FVElementGeometry, class ElemVolVars>
158 static auto velocityGradII(const FVElementGeometry& fvGeometry,
159 const typename FVElementGeometry::SubControlVolumeFace& scvf,
160 const ElemVolVars& elemVolVars)
161 {
162 assert(scvf.isFrontal());
163 // The velocities of the dof at interest and the one of the opposite scvf.
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();
167
168 return (velocityOpposite - velocitySelf) / distance * scvf.directionSign();
169 }
170
192 template<class FVElementGeometry, class ElemVolVars>
193 static auto velocityGradIJ(const FVElementGeometry& fvGeometry,
194 const typename FVElementGeometry::SubControlVolumeFace& scvf,
195 const ElemVolVars& elemVolVars)
196 {
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();
202 }
203
233 template<class FVElementGeometry, class ElemVolVars>
234 static auto velocityGradJI(const FVElementGeometry& fvGeometry,
235 const typename FVElementGeometry::SubControlVolumeFace& scvf,
236 const ElemVolVars& elemVolVars)
237 {
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();
244 }
245
246private:
247
248 template<class FVElementGeometry>
249 static auto getDistanceIJ_(const FVElementGeometry& fvGeometry,
250 const typename FVElementGeometry::SubControlVolumeFace& scvf)
251 {
252 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
253 // The scvf is on a boundary, hence there is no outer DOF.
254 // We take the distance to the boundary instead.
255 if (scvf.boundary())
256 return getDistanceToBoundary_(insideScv, scvf);
257
258 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
259
260 // The standard case: Our grid does not support periodic boundaries.
262 return getStandardDistance_(insideScv, outsideScv);
263 else
264 {
265 const auto& gridGeometry = fvGeometry.gridGeometry();
266 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
267 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
268
269 // The standard case: Our grid is not periodic or the lateral scvf does not lie on a periodic boundary.
270 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex()))
271 return getStandardDistance_(insideScv, outsideScv);
272
273 // Treating periodic boundaries is more involved:
274 // 1. Consider the outside scv within the element adjacent to the other periodic boundary.
275 // 2. Iterate over this scv's faces until you find the face parallel to our own scvf.
276 // This face would lie directly next to our own scvf if the grid was not periodic.
277 // 3. Calculate the total distance by summing up the distances between the DOFs and the respective integration points
278 // corresponding to the inside and outside scvfs.
279 auto periodicFvGeometry = localView(gridGeometry);
280 const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex());
281 periodicFvGeometry.bindElement(periodicElement);
282
283 for (const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv))
284 {
285 if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign())
286 {
287 const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
288 const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm();
289 return insideDistance + outsideDistance;
290 }
291 }
292 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element");
293 }
294 }
295
297 template<class FVElementGeometry>
298 static auto getDistanceJI_(const FVElementGeometry& fvGeometry,
299 const typename FVElementGeometry::SubControlVolumeFace& scvf,
300 const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf)
301 {
302 const auto& orthogonalInsideScv = fvGeometry.scv(orthogonalScvf.insideScvIdx());
303
304 // The orthogonal scvf is on a boundary, hence there is no outer DOF.
305 // We take the distance to the boundary instead.
306 if (orthogonalScvf.boundary())
307 return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf);
308
309 const auto& orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx());
310
311 // The standard case: Our grid does not support periodic boundaries.
312 if constexpr (!Detail::SupportsPeriodicity<typename FVElementGeometry::GridGeometry::Grid>())
313 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
314 else
315 {
316 const auto& gridGeometry = fvGeometry.gridGeometry();
317 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
318
319 // The standard case: Our grid is not periodic or our own DOF does not lie on a periodic boundary.
320 if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex()))
321 return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv);
322
323 // Treating periodic boundaries is more involved:
324 // 1. Consider the orthogonal outside scv within the element adjacent to the other periodic boundary.
325 // 2. Iterate over this scv's faces until you find the orthogonal face parallel to our own orthogonal scvf.
326 // This face would lie directly next to our own orthogonal scvf if the grid was not periodic.
327 // 3. Calculate the total distance by summing up the distances between the DOFs and the respective integration points
328 // corresponding to the orthogonal scvfs.
329 auto periodicFvGeometry = localView(gridGeometry);
330 const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex());
331 periodicFvGeometry.bindElement(periodicElement);
332
333 for (const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv))
334 {
335 if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign())
336 {
337 const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm();
338 const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm();
339 return insideDistance + outsideDistance;
340 }
341 }
342 DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element");
343 }
344 }
345
346 template<class SubControlVolume>
347 static auto getStandardDistance_(const SubControlVolume& insideScv,
348 const SubControlVolume& outsideScv)
349 {
350 return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm();
351 }
352
353 template<class SubControlVolume, class SubControlVolumeFace>
354 static auto getDistanceToBoundary_(const SubControlVolume& scv,
355 const SubControlVolumeFace& scvf)
356 {
357 return (scv.dofPosition() - scvf.ipGlobal()).two_norm();
358 }
359};
360
361} // end namespace Dumux
362
363#endif // DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH
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.