Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup NavierStokesModel | ||
10 | * \copydoc Dumux::StaggeredVelocityGradients | ||
11 | */ | ||
12 | #ifndef DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH | ||
13 | #define DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH | ||
14 | |||
15 | #include <optional> | ||
16 | #include <dumux/common/exceptions.hh> | ||
17 | #include <dumux/common/parameters.hh> | ||
18 | #include <dune/common/fmatrix.hh> | ||
19 | #include <dumux/common/typetraits/problem.hh> | ||
20 | #include <dumux/discretization/facecentered/staggered/elementboundarytypes.hh> | ||
21 | #include <dumux/io/grid/periodicgridtraits.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup NavierStokesModel | ||
27 | * \brief Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered grid discretization. | ||
28 | */ | ||
29 | class StaggeredVelocityGradients | ||
30 | { | ||
31 | public: | ||
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; | ||
41 | using ElementBoundaryTypes = FaceCenteredStaggeredElementBoundaryTypes<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 | 204950960 | static auto velocityGradient(const FVElementGeometry& fvGeometry, | |
49 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
50 | const ElemVolVars& elemVolVars, | ||
51 | const FaceCenteredStaggeredElementBoundaryTypes<BoundaryTypes>& elemBcTypes, | ||
52 | bool fullGradient = false) | ||
53 | { | ||
54 | using Scalar = typename FVElementGeometry::GridGeometry::GlobalCoordinate::value_type; | ||
55 | static constexpr auto dim = FVElementGeometry::GridGeometry::GlobalCoordinate::dimension; | ||
56 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 200877856 times.
|
409901920 | Dune::FieldMatrix<Scalar, dim, dim> gradient(0.0); |
57 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 200877856 times.
✓ Branch 2 taken 4073104 times.
✓ Branch 3 taken 200877856 times.
|
204950960 | const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); |
58 | |||
59 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 204950960 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
204950960 | 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 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 204950960 times.
|
204950960 | if (fullGradient) |
105 | ✗ | DUNE_THROW(Dune::NotImplemented, "Full gradient for lateral staggered faces not implemented"); | |
106 | |||
107 | 204950960 | gradient[scv.dofAxis()][scvf.normalAxis()] = velocityGradIJ(fvGeometry, scvf, elemVolVars); | |
108 | 204950960 | gradient[scvf.normalAxis()][scv.dofAxis()] = velocityGradJI(fvGeometry, scvf, elemVolVars); | |
109 | } | ||
110 | |||
111 | 204950960 | return gradient; | |
112 | } | ||
113 | |||
114 | /*! | ||
115 | * \brief Returns the in-axis velocity gradient. | ||
116 | * | ||
117 | * \verbatim | ||
118 | * ---------======= == and # staggered half-control-volume | ||
119 | * | # | current scv | ||
120 | * | # | # staggered face over which fluxes are calculated | ||
121 | * vel.Opp <~~| O~~> x~~~~> vel.Self | ||
122 | * | # | x dof position | ||
123 | * | # | | ||
124 | * --------======== -- element | ||
125 | * | ||
126 | * O position at which gradient is evaluated (integration point) | ||
127 | * \endverbatim | ||
128 | */ | ||
129 | template<class FVElementGeometry, class ElemVolVars> | ||
130 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 102701186 times.
|
102701186 | static auto velocityGradII(const FVElementGeometry& fvGeometry, |
131 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
132 | const ElemVolVars& elemVolVars) | ||
133 | { | ||
134 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 102701186 times.
|
102701186 | assert(scvf.isFrontal()); |
135 | // The velocities of the dof at interest and the one of the opposite scvf. | ||
136 | 102701186 | const auto velocitySelf = elemVolVars[scvf.insideScvIdx()].velocity(); | |
137 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 100648618 times.
|
102701186 | const auto velocityOpposite = elemVolVars[scvf.outsideScvIdx()].velocity(); |
138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 100648618 times.
|
205400822 | const auto distance = (fvGeometry.scv(scvf.outsideScvIdx()).dofPosition() - fvGeometry.scv(scvf.insideScvIdx()).dofPosition()).two_norm(); |
139 | |||
140 | 102701186 | return (velocityOpposite - velocitySelf) / distance * scvf.directionSign(); | |
141 | } | ||
142 | |||
143 | /*! | ||
144 | * \brief Returns the velocity gradient perpendicular to the orientation of our current scvf. | ||
145 | * | ||
146 | * \verbatim | ||
147 | * ---------------- | ||
148 | * | |outer || and # staggered half-control-volume (own element) | ||
149 | * | |vel. gradient | ||
150 | * | |~~~~> -------> :: staggered half-control-volume (neighbor element) | ||
151 | * | | ------> | ||
152 | * | lateral scvf | -----> x dof position | ||
153 | * ---------######O::::::::: ----> | ||
154 | * | || | :: ---> -- elements | ||
155 | * | || | :: --> | ||
156 | * | || scv x~~~~> :: -> O position at which gradient is evaluated (integration point) | ||
157 | * | || | inner :: | ||
158 | * | || | vel. :: | ||
159 | * ---------#######::::::::: | ||
160 | * | ||
161 | * | ||
162 | * \endverbatim | ||
163 | */ | ||
164 | template<class FVElementGeometry, class ElemVolVars> | ||
165 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 205639008 times.
|
205639008 | static auto velocityGradIJ(const FVElementGeometry& fvGeometry, |
166 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
167 | const ElemVolVars& elemVolVars) | ||
168 | { | ||
169 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 205639008 times.
|
205639008 | assert(scvf.isLateral()); |
170 | 205639008 | const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity(); | |
171 | 205639008 | const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity(); | |
172 | 205639008 | const auto distance = getDistanceIJ_(fvGeometry, scvf); | |
173 | 205639008 | return (outerVelocity - innerVelocity) / distance * scvf.directionSign(); | |
174 | } | ||
175 | |||
176 | /*! | ||
177 | * \brief Returns the velocity gradient in line with our current scvf. | ||
178 | * | ||
179 | * \verbatim | ||
180 | * ^ gradient | ||
181 | * | ^ | ||
182 | * | | ^ | ||
183 | * | | | ^ | ||
184 | * | | | | ^ | ||
185 | * | | | | | ^ || and # staggered half-control-volume (own element) | ||
186 | * | | | | | | | ||
187 | * :: staggered half-control-volume (neighbor element) | ||
188 | * ---------------- | ||
189 | * | inner | outer x dof position (of own scv) | ||
190 | * | vel. | vel. | ||
191 | * | ^ | ^ -- elements | ||
192 | * | | lat. | | | ||
193 | * | | scvf | | O position at which gradient is evaluated (integration point) | ||
194 | * ---------######O::::::::: | ||
195 | * | || ~ :: ~ ~ orthogonal scvf | ||
196 | * | || ~ :: | ||
197 | * | || scv x :: | ||
198 | * | || | :: | ||
199 | * | || | :: | ||
200 | * ---------#######::::::::: | ||
201 | * | ||
202 | * | ||
203 | * \endverbatim | ||
204 | */ | ||
205 | template<class FVElementGeometry, class ElemVolVars> | ||
206 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 205625362 times.
|
205625362 | static auto velocityGradJI(const FVElementGeometry& fvGeometry, |
207 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
208 | const ElemVolVars& elemVolVars) | ||
209 | { | ||
210 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 205625362 times.
|
205625362 | assert(scvf.isLateral()); |
211 | 205625362 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); | |
212 | 205625362 | const auto innerVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity(); | |
213 | 205625362 | const auto outerVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity(); | |
214 | 205625362 | const auto distance = getDistanceJI_(fvGeometry, scvf, orthogonalScvf); | |
215 | 205625362 | return (outerVelocity - innerVelocity) / distance * orthogonalScvf.directionSign(); | |
216 | } | ||
217 | |||
218 | private: | ||
219 | |||
220 | template<class FVElementGeometry> | ||
221 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 201548656 times.
|
205639008 | static auto getDistanceIJ_(const FVElementGeometry& fvGeometry, |
222 | const typename FVElementGeometry::SubControlVolumeFace& scvf) | ||
223 | { | ||
224 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 201633456 times.
✓ Branch 2 taken 7029506 times.
✓ Branch 3 taken 198524702 times.
|
205639008 | 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 |
2/2✓ Branch 0 taken 3108754 times.
✓ Branch 1 taken 202530254 times.
|
205639008 | if (scvf.boundary()) |
228 | 3108754 | return getDistanceToBoundary_(insideScv, scvf); | |
229 | |||
230 | 202530254 | 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 | 199031342 | return getStandardDistance_(insideScv, outsideScv); | |
235 | else | ||
236 | { | ||
237 | 3498912 | const auto& gridGeometry = fvGeometry.gridGeometry(); | |
238 |
1/2✓ Branch 1 taken 3498912 times.
✗ Branch 2 not taken.
|
3498912 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); |
239 |
1/2✓ Branch 0 taken 3498912 times.
✗ Branch 1 not taken.
|
3498912 | 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 |
3/4✓ Branch 0 taken 3498912 times.
✗ Branch 1 not taken.
✓ Branch 3 taken 3480124 times.
✓ Branch 4 taken 18788 times.
|
6997824 | if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(orthogonalInsideScv.dofIndex())) |
243 | 3480124 | 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 | 18788 | auto periodicFvGeometry = localView(gridGeometry); | |
252 | 18788 | const auto& periodicElement = gridGeometry.element(outsideScv.elementIndex()); | |
253 | 18788 | periodicFvGeometry.bindElement(periodicElement); | |
254 | |||
255 |
3/4✓ Branch 1 taken 28180 times.
✓ Branch 2 taken 18788 times.
✓ Branch 4 taken 46968 times.
✗ Branch 5 not taken.
|
103328 | for (const auto& outsideScvf : scvfs(periodicFvGeometry, outsideScv)) |
256 | { | ||
257 |
4/4✓ Branch 0 taken 28180 times.
✓ Branch 1 taken 18788 times.
✓ Branch 2 taken 9392 times.
✓ Branch 3 taken 18788 times.
|
46968 | if (outsideScvf.normalAxis() == scvf.normalAxis() && outsideScvf.directionSign() != scvf.directionSign()) |
258 | { | ||
259 | 56364 | const auto insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm(); | |
260 | 37576 | const auto outsideDistance = (outsideScv.dofPosition() - outsideScvf.ipGlobal()).two_norm(); | |
261 | 18788 | return insideDistance + outsideDistance; | |
262 | } | ||
263 | } | ||
264 |
0/20✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
3498912 | DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element"); |
265 | } | ||
266 | } | ||
267 | |||
268 | //! Get the distance between the DOFs at which the inner and outer velocities are defined. | ||
269 | template<class FVElementGeometry> | ||
270 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 201537474 times.
|
205625362 | static auto getDistanceJI_(const FVElementGeometry& fvGeometry, |
271 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
272 | const typename FVElementGeometry::SubControlVolumeFace& orthogonalScvf) | ||
273 | { | ||
274 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 201537474 times.
✓ Branch 2 taken 4091472 times.
✓ Branch 3 taken 201533890 times.
|
205625362 | 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 |
2/2✓ Branch 0 taken 3584 times.
✓ Branch 1 taken 205621778 times.
|
205625362 | if (orthogonalScvf.boundary()) |
279 | 3584 | return getDistanceToBoundary_(orthogonalInsideScv, orthogonalScvf); | |
280 | |||
281 |
1/2✓ Branch 0 taken 3680944 times.
✗ Branch 1 not taken.
|
205621778 | 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 | 201940834 | return getStandardDistance_(orthogonalInsideScv, orthogonalOutsideScv); | |
286 | else | ||
287 | { | ||
288 |
1/2✓ Branch 0 taken 3680944 times.
✗ Branch 1 not taken.
|
3680944 | const auto& gridGeometry = fvGeometry.gridGeometry(); |
289 | 3680944 | 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 |
3/4✓ Branch 0 taken 3680944 times.
✗ Branch 1 not taken.
✓ Branch 3 taken 3660592 times.
✓ Branch 4 taken 20352 times.
|
7361888 | if (!gridGeometry.isPeriodic() || !gridGeometry.dofOnPeriodicBoundary(insideScv.dofIndex())) |
293 | 3660592 | 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 | 20352 | auto periodicFvGeometry = localView(gridGeometry); | |
302 | 20352 | const auto& periodicElement = gridGeometry.element(orthogonalOutsideScv.elementIndex()); | |
303 | 20352 | periodicFvGeometry.bindElement(periodicElement); | |
304 | |||
305 |
3/4✓ Branch 1 taken 30528 times.
✓ Branch 2 taken 20352 times.
✓ Branch 4 taken 50880 times.
✗ Branch 5 not taken.
|
111936 | for (const auto& outsideOrthogonalScvf : scvfs(periodicFvGeometry, orthogonalOutsideScv)) |
306 | { | ||
307 |
4/4✓ Branch 0 taken 30528 times.
✓ Branch 1 taken 20352 times.
✓ Branch 2 taken 10176 times.
✓ Branch 3 taken 20352 times.
|
50880 | if (outsideOrthogonalScvf.normalAxis() == orthogonalScvf.normalAxis() && outsideOrthogonalScvf.directionSign() != orthogonalScvf.directionSign()) |
308 | { | ||
309 | 61056 | const auto insideDistance = (orthogonalInsideScv.dofPosition() - orthogonalScvf.ipGlobal()).two_norm(); | |
310 | 40704 | const auto outsideDistance = (orthogonalOutsideScv.dofPosition() - outsideOrthogonalScvf.ipGlobal()).two_norm(); | |
311 | 20352 | return insideDistance + outsideDistance; | |
312 | } | ||
313 | } | ||
314 |
0/20✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
3680944 | DUNE_THROW(Dune::InvalidStateException, "Could not find scvf in periodic element"); |
315 | } | ||
316 | } | ||
317 | |||
318 | template<class SubControlVolume> | ||
319 | 408112892 | static auto getStandardDistance_(const SubControlVolume& insideScv, | |
320 | const SubControlVolume& outsideScv) | ||
321 | { | ||
322 | 816225784 | return (insideScv.dofPosition() - outsideScv.dofPosition()).two_norm(); | |
323 | } | ||
324 | |||
325 | template<class SubControlVolume, class SubControlVolumeFace> | ||
326 | 3112338 | static auto getDistanceToBoundary_(const SubControlVolume& scv, | |
327 | const SubControlVolumeFace& scvf) | ||
328 | { | ||
329 | 6224676 | return (scv.dofPosition() - scvf.ipGlobal()).two_norm(); | |
330 | } | ||
331 | }; | ||
332 | |||
333 | } // end namespace Dumux | ||
334 | |||
335 | #endif // DUMUX_NAVIERSTOKES_STAGGERED_VELOCITYGRADIENTS_HH | ||
336 |