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-FileCopyrightInfo: 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::NavierStokesFluxVariablesImpl | ||
11 | */ | ||
12 | #ifndef DUMUX_NAVIERSTOKES_MOMENTUM_FLUXVARIABLES_HH | ||
13 | #define DUMUX_NAVIERSTOKES_MOMENTUM_FLUXVARIABLES_HH | ||
14 | |||
15 | #include <array> | ||
16 | |||
17 | #include <dumux/common/numeqvector.hh> | ||
18 | #include <dumux/common/math.hh> | ||
19 | #include <dumux/common/exceptions.hh> | ||
20 | #include <dumux/common/parameters.hh> | ||
21 | #include <dumux/common/properties.hh> | ||
22 | |||
23 | #include <dumux/discretization/extrusion.hh> | ||
24 | #include <dumux/discretization/method.hh> | ||
25 | |||
26 | // #include "staggeredupwindhelper.hh" | ||
27 | #include "velocitygradients.hh" | ||
28 | |||
29 | namespace Dumux { | ||
30 | |||
31 | /*! | ||
32 | * \ingroup NavierStokesModel | ||
33 | * \brief The flux variables class for the Navier-Stokes model using the staggered grid discretization. | ||
34 | */ | ||
35 | template<class TypeTag> | ||
36 | class NavierStokesMomentumFluxVariables | ||
37 | { | ||
38 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
39 | |||
40 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
41 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
42 | using VolumeVariables = typename GridVolumeVariables::VolumeVariables; | ||
43 | |||
44 | using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache; | ||
45 | using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView; | ||
46 | using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache; | ||
47 | |||
48 | using GridGeometry = typename GridVariables::GridGeometry; | ||
49 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
50 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
51 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
52 | |||
53 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
54 | |||
55 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
56 | using GridView = typename GridGeometry::GridView; | ||
57 | using Element = typename GridView::template Codim<0>::Entity; | ||
58 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
59 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
60 | using Indices = typename ModelTraits::Indices; | ||
61 | using VelocityGradients = StaggeredVelocityGradients; | ||
62 | |||
63 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
64 | using Extrusion = Extrusion_t<GridGeometry>; | ||
65 | |||
66 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
67 | |||
68 | public: | ||
69 | |||
70 | 279435948 | NavierStokesMomentumFluxVariables(const Problem& problem, | |
71 | const Element& element, | ||
72 | const FVElementGeometry& fvGeometry, | ||
73 | const SubControlVolumeFace& scvFace, | ||
74 | const ElementVolumeVariables& elemVolVars, | ||
75 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
76 | const ElementBoundaryTypes& elemBcTypes) | ||
77 | : problemPtr_(&problem) | ||
78 | , elementPtr_(&element) | ||
79 | , fvGeometryPtr_(&fvGeometry) | ||
80 | , scvFacePtr_(&scvFace) | ||
81 | , elemVolVarsPtr_(&elemVolVars) | ||
82 | , elemFluxVarsCachePtr_(&elemFluxVarsCache) | ||
83 | 279435948 | , elemBcTypesPtr_(&elemBcTypes) | |
84 | { | ||
85 | static_assert( | ||
86 | std::decay_t<decltype(problem.dirichlet(element, scvFace))>::size() | ||
87 | == static_cast<std::size_t>(GridView::dimension), | ||
88 | "Expects problem.dirichlet to return an array with as many entries as dimensions." | ||
89 | ); | ||
90 | } | ||
91 | |||
92 | ✗ | const Problem& problem() const | |
93 | ✗ | { return *problemPtr_; } | |
94 | |||
95 | ✗ | const Element& element() const | |
96 | ✗ | { return *elementPtr_; } | |
97 | |||
98 | ✗ | const SubControlVolumeFace& scvFace() const | |
99 | ✗ | { return *scvFacePtr_; } | |
100 | |||
101 | ✗ | const FVElementGeometry& fvGeometry() const | |
102 | ✗ | { return *fvGeometryPtr_; } | |
103 | |||
104 | ✗ | const ElementVolumeVariables& elemVolVars() const | |
105 | ✗ | { return *elemVolVarsPtr_; } | |
106 | |||
107 | const ElementFluxVariablesCache& elemFluxVarsCache() const | ||
108 | { return *elemFluxVarsCachePtr_; } | ||
109 | |||
110 | ✗ | const ElementBoundaryTypes& elemBcTypes() const | |
111 | ✗ | { return *elemBcTypesPtr_; } | |
112 | |||
113 | /*! | ||
114 | * \brief Returns the diffusive momentum flux due to viscous forces | ||
115 | */ | ||
116 | 279435948 | NumEqVector advectiveMomentumFlux() const | |
117 | { | ||
118 |
2/2✓ Branch 0 taken 83467152 times.
✓ Branch 1 taken 195968796 times.
|
279435948 | if (!this->problem().enableInertiaTerms()) |
119 | 83467152 | return NumEqVector(0.0); | |
120 | |||
121 |
2/2✓ Branch 0 taken 65276866 times.
✓ Branch 1 taken 130691930 times.
|
195968796 | if (this->scvFace().isFrontal()) |
122 | 65276866 | return frontalAdvectiveMomentumFlux(); | |
123 | else | ||
124 | 130691930 | return lateralAdvectiveMomentumFlux(); | |
125 | } | ||
126 | |||
127 | /*! | ||
128 | * \brief Returns the diffusive momentum flux due to viscous forces | ||
129 | */ | ||
130 | 279435948 | NumEqVector diffusiveMomentumFlux() const | |
131 | { | ||
132 |
2/2✓ Branch 0 taken 93209932 times.
✓ Branch 1 taken 186226016 times.
|
279435948 | if (this->scvFace().isFrontal()) |
133 | 93209932 | return frontalDiffusiveMomentumFlux(); | |
134 | else | ||
135 | 186226016 | return lateralDiffusiveMomentumFlux(); | |
136 | } | ||
137 | |||
138 | /*! | ||
139 | * \brief Returns the frontal part of the momentum flux. | ||
140 | * This treats the flux over the staggered face at the center of an element, | ||
141 | * parallel to the current scvf where the velocity dof of interest lives. | ||
142 | * | ||
143 | * \verbatim | ||
144 | * scvf | ||
145 | * ---------======= == and # staggered half-control-volume | ||
146 | * | # | current scv | ||
147 | * | # | # staggered face over which fluxes are calculated | ||
148 | * vel.Opp |~~> O~~~> x~~~~> vel.Self | ||
149 | * | # | x dof position | ||
150 | * | # | | ||
151 | * --------======== -- element | ||
152 | * scvf | ||
153 | * O integration point | ||
154 | * \endverbatim | ||
155 | */ | ||
156 | 93209932 | NumEqVector frontalDiffusiveMomentumFlux() const | |
157 | { | ||
158 | 93209932 | const auto& scvf = this->scvFace(); | |
159 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 93209932 times.
|
93209932 | assert(scvf.isFrontal()); |
160 | |||
161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 93209932 times.
|
93209932 | NumEqVector result(0.0); |
162 | 93209932 | const auto& fvGeometry = this->fvGeometry(); | |
163 | 93209932 | const auto& elemVolVars = this->elemVolVars(); | |
164 | |||
165 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 93209932 times.
|
93209932 | if (scvf.boundary()) |
166 | ✗ | return result; | |
167 | |||
168 | // get the velocity gradient at the normal face's integration point | ||
169 | 93209932 | const auto velGradII = VelocityGradients::velocityGradII(fvGeometry, scvf, elemVolVars); | |
170 | |||
171 |
4/4✓ Branch 0 taken 83 times.
✓ Branch 1 taken 93209849 times.
✓ Branch 3 taken 60 times.
✓ Branch 4 taken 23 times.
|
93209992 | static const bool enableUnsymmetrizedVelocityGradient |
172 |
2/4✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 60 times.
✗ Branch 5 not taken.
|
180 | = getParamFromGroup<bool>(this->problem().paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false); |
173 |
4/6✓ Branch 0 taken 60 times.
✓ Branch 1 taken 93209872 times.
✓ Branch 3 taken 60 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 60 times.
✗ Branch 6 not taken.
|
93209992 | static const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0; |
174 | |||
175 | 93209932 | const auto mu = this->problem().effectiveViscosity(this->element(), this->fvGeometry(), this->scvFace()); | |
176 |
2/2✓ Branch 3 taken 104 times.
✓ Branch 4 taken 93209828 times.
|
279629796 | result -= factor * mu * velGradII * Extrusion::area(fvGeometry, scvf) * elemVolVars[scvf.insideScvIdx()].extrusionFactor() * scvf.directionSign(); |
177 | |||
178 |
6/8✓ Branch 0 taken 104 times.
✓ Branch 1 taken 93209828 times.
✓ Branch 3 taken 60 times.
✓ Branch 4 taken 44 times.
✓ Branch 6 taken 60 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 60 times.
✗ Branch 10 not taken.
|
93209932 | static const bool enableDilatationTerm = getParamFromGroup<bool>(this->problem().paramGroup(), "FreeFlow.EnableDilatationTerm", false); |
179 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 93209932 times.
|
93209932 | if (enableDilatationTerm) |
180 | { | ||
181 | ✗ | Scalar divergence = velGradII; | |
182 | ✗ | for (const auto& scv : scvs(fvGeometry)) | |
183 | { | ||
184 | ✗ | const auto otherFrontalScvf = *(scvfs(fvGeometry, scv).begin()); | |
185 | ✗ | assert(otherFrontalScvf.isFrontal() && !otherFrontalScvf.boundary()); | |
186 | ✗ | if (otherFrontalScvf.index() != scvf.index()) | |
187 | ✗ | divergence += VelocityGradients::velocityGradII(fvGeometry, otherFrontalScvf, elemVolVars); | |
188 | } | ||
189 | |||
190 | ✗ | result += 2.0/3.0 * mu * divergence * scvf.directionSign() * Extrusion::area(fvGeometry, scvf) * elemVolVars[scvf.insideScvIdx()].extrusionFactor(); | |
191 | } | ||
192 | |||
193 | 93209932 | return result; | |
194 | } | ||
195 | |||
196 | /*! | ||
197 | * \brief Returns the diffusive momentum flux over the staggered face perpendicular to the scvf | ||
198 | * where the velocity dof of interest lives (coinciding with the element's scvfs). | ||
199 | * | ||
200 | * \verbatim | ||
201 | * ---------------- | ||
202 | * | |vel. | ||
203 | * | in.norm. |Parallel | ||
204 | * | vel. |~~~~> | ||
205 | * | ^ | ^ out.norm.vel. | ||
206 | * | | | | | ||
207 | * scvf ---------#######::::::::: || and # staggered half-control-volume (own element) | ||
208 | * | || | curr. :: | ||
209 | * | || | scvf :: :: staggered half-control-volume (neighbor element) | ||
210 | * | || x~~~~> :: | ||
211 | * | || | vel. :: # lateral staggered faces over which fluxes are calculated | ||
212 | * scvf | || | Self :: | ||
213 | * ---------#######::::::::: x dof position | ||
214 | * scvf | ||
215 | * -- elements | ||
216 | * \endverbatim | ||
217 | */ | ||
218 | 186226016 | NumEqVector lateralDiffusiveMomentumFlux() const | |
219 | { | ||
220 | 186226016 | const auto& scvf = this->scvFace(); | |
221 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 186226016 times.
|
186226016 | assert(scvf.isLateral()); |
222 | |||
223 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 182152912 times.
|
186226016 | NumEqVector result(0.0); |
224 | 186226016 | const auto& fvGeometry = this->fvGeometry(); | |
225 | 186226016 | const auto& elemVolVars = this->elemVolVars(); | |
226 | 186226016 | const auto& problem = this->problem(); | |
227 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 182152912 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 182152912 times.
|
372452032 | const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); |
228 | |||
229 |
4/4✓ Branch 0 taken 81 times.
✓ Branch 1 taken 186225935 times.
✓ Branch 3 taken 59 times.
✓ Branch 4 taken 22 times.
|
186226075 | static const bool enableUnsymmetrizedVelocityGradient |
230 |
2/4✓ Branch 1 taken 59 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 59 times.
✗ Branch 5 not taken.
|
177 | = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false); |
231 | |||
232 | 186226016 | const auto mu = this->problem().effectiveViscosity(this->element(), this->fvGeometry(), this->scvFace()); | |
233 | |||
234 | // get the velocity gradient at the lateral face's integration point | ||
235 | 186226016 | const auto gradV = VelocityGradients::velocityGradient(fvGeometry, scvf, elemVolVars, this->elemBcTypes(), false); | |
236 | |||
237 | // Consider the shear stress caused by the gradient of the velocities parallel to our face of interest. | ||
238 | 186226016 | GlobalPosition gradVn(0.0); | |
239 | 371256136 | gradV.mv(scvf.unitOuterNormal(), gradVn); | |
240 |
1/2✓ Branch 0 taken 186226016 times.
✗ Branch 1 not taken.
|
186226016 | const Scalar velocityGrad_ij = gradVn[scv.dofAxis()]; |
241 |
1/2✓ Branch 0 taken 186226016 times.
✗ Branch 1 not taken.
|
186226016 | result -= mu * velocityGrad_ij; |
242 | |||
243 | // Consider the shear stress caused by the gradient of the velocities normal to our face of interest. | ||
244 |
1/2✓ Branch 0 taken 186226016 times.
✗ Branch 1 not taken.
|
186226016 | if (!enableUnsymmetrizedVelocityGradient) |
245 | { | ||
246 | 186226016 | GlobalPosition gradVTransposedN(0.0); | |
247 | 371256136 | gradV.mtv(scvf.unitOuterNormal(), gradVTransposedN); | |
248 | 186226016 | const Scalar velocityGrad_ji = gradVTransposedN[scv.dofAxis()]; | |
249 | 372452032 | result -= mu * velocityGrad_ji; | |
250 | } | ||
251 | |||
252 | // Account for the area of the staggered lateral face. | ||
253 | 372452032 | return result * Extrusion::area(fvGeometry, scvf) * elemVolVars[scvf.insideScvIdx()].extrusionFactor(); | |
254 | } | ||
255 | |||
256 | /*! | ||
257 | * \brief Returns the frontal pressure contribution. | ||
258 | * | ||
259 | * \verbatim | ||
260 | * | ||
261 | * ---------======= == and # staggered half-control-volume | ||
262 | * | # | current scv | ||
263 | * | # | # frontal staggered face for which pressure contribution is calculated | ||
264 | * | P x | ||
265 | * | # | x dof position | ||
266 | * | # | | ||
267 | * --------======== -- element | ||
268 | * | ||
269 | * P integration point, pressure DOF lives here | ||
270 | * \endverbatim | ||
271 | */ | ||
272 | 279435948 | NumEqVector pressureContribution() const | |
273 | { | ||
274 |
2/2✓ Branch 0 taken 93209932 times.
✓ Branch 1 taken 186226016 times.
|
279435948 | NumEqVector result(0.0); |
275 | 279435948 | const auto& scvf = this->scvFace(); | |
276 |
3/4✓ Branch 0 taken 93209932 times.
✓ Branch 1 taken 186226016 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 93209932 times.
|
279435948 | if (scvf.isLateral() || scvf.boundary()) |
277 | 186226016 | return result; | |
278 | |||
279 | // The pressure force needs to take the extruded scvf area into account. | ||
280 | 93209932 | const auto pressure = this->problem().pressure(this->element(), this->fvGeometry(), scvf); | |
281 | 279629796 | result = pressure*Extrusion::area(this->fvGeometry(), scvf)*this->elemVolVars()[scvf.insideScvIdx()].extrusionFactor(); | |
282 | |||
283 | // The pressure contribution calculated above might have a much larger numerical value compared to the viscous or inertial forces. | ||
284 | // This may lead to numerical inaccuracies due to loss of significance (cancellantion) for the final residual value. | ||
285 | // In the end, we are only interested in a pressure difference between the two relevant faces so we can | ||
286 | // subtract a reference value from the actual pressure contribution. Assuming an axisparallel cartesian grid, | ||
287 | // scvf.area() will have the same value at both opposing faces such that the reference pressure contribution | ||
288 | // cancels out in the final residual which combines the pressure contribution of two adjacent elements | ||
289 | // We explicitly do extrude the area here because that might yield different results in both elements. | ||
290 | // The multiplication by scvf.area() aims at having a reference value of the same order of magnitude as the actual pressure contribution. | ||
291 | 93209932 | const auto referencePressure = this->problem().referencePressure(this->element(), this->fvGeometry(), scvf); | |
292 | 93209932 | result -= referencePressure*scvf.area(); | |
293 | |||
294 | // Account for the orientation of the face. | ||
295 | 93209932 | result *= scvf.directionSign(); | |
296 | 93209932 | return result; | |
297 | } | ||
298 | |||
299 | /*! | ||
300 | * \brief Returns the frontal part of the momentum flux. | ||
301 | * This treats the flux over the staggered face at the center of an element, | ||
302 | * parallel to the current scvf where the velocity dof of interest lives. | ||
303 | * | ||
304 | * \verbatim | ||
305 | * scvf | ||
306 | * ---------======= == and # staggered half-control-volume | ||
307 | * | # | current scv | ||
308 | * | # | # staggered face over which fluxes are calculated | ||
309 | * vel.Opp |~~> O~~~> x~~~~> vel.Self | ||
310 | * | # | x dof position | ||
311 | * | # | | ||
312 | * --------======== -- element | ||
313 | * scvf | ||
314 | * O integration point | ||
315 | * \endverbatim | ||
316 | */ | ||
317 | 65276866 | NumEqVector frontalAdvectiveMomentumFlux() const | |
318 | { | ||
319 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 65276866 times.
|
65276866 | NumEqVector flux(0.0); |
320 | 65276866 | const auto& scvf = this->scvFace(); | |
321 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 65276866 times.
|
65276866 | assert(scvf.isFrontal()); |
322 | |||
323 | 65276866 | const auto& problem = this->problem(); | |
324 | 65276866 | const auto& elemVolVars = this->elemVolVars(); | |
325 | 130553732 | const auto velocitySelf = elemVolVars[scvf.insideScvIdx()].velocity(); | |
326 |
0/2✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
130553732 | const auto velocityOpposite = elemVolVars[scvf.outsideScvIdx()].velocity(); |
327 | |||
328 | // Get the average velocity at the center of the element (i.e. the location of the staggered face). | ||
329 | 65276866 | const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5; | |
330 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
65276866 | const Scalar density = this->problem().density(this->element(), this->fvGeometry(), scvf); |
331 |
2/2✓ Branch 0 taken 23 times.
✓ Branch 1 taken 65276843 times.
|
65276866 | const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity); |
332 | |||
333 | // TODO use higher order helper | ||
334 |
5/8✓ Branch 0 taken 23 times.
✓ Branch 1 taken 65276843 times.
✓ Branch 3 taken 23 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 23 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 23 times.
✗ Branch 10 not taken.
|
65276866 | static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight"); |
335 |
2/2✓ Branch 0 taken 29319344 times.
✓ Branch 1 taken 35957522 times.
|
65276866 | const Scalar transportedMomentum = selfIsUpstream ? (upwindWeight * velocitySelf + (1.0 - upwindWeight) * velocityOpposite) * density |
336 | 35957522 | : (upwindWeight * velocityOpposite + (1.0 - upwindWeight) * velocitySelf) * density; | |
337 | |||
338 | 130553732 | return transportingVelocity * transportedMomentum * scvf.directionSign() * Extrusion::area(this->fvGeometry(), scvf) * extrusionFactor_(elemVolVars, scvf); | |
339 | } | ||
340 | |||
341 | /*! | ||
342 | * \brief Returns the advective momentum flux over the staggered face perpendicular to the scvf | ||
343 | * where the velocity dof of interest lives (coinciding with the element's scvfs). | ||
344 | * | ||
345 | * \verbatim | ||
346 | * ---------------- | ||
347 | * | inner | | ||
348 | * | transp. | | ||
349 | * | vel. |~~~~> outer vel. | ||
350 | * | ^ | | ||
351 | * | | | | ||
352 | * ---------######O || and # staggered half-control-volume | ||
353 | * | || | scv | ||
354 | * | || | # lateral staggered faces over which fluxes are calculated | ||
355 | * | || x~~~~> inner vel. | ||
356 | * | || | x dof position | ||
357 | * | || | | ||
358 | * ---------####### -- elements | ||
359 | * | ||
360 | * O integration point | ||
361 | * \endverbatim | ||
362 | */ | ||
363 | 130691930 | NumEqVector lateralAdvectiveMomentumFlux() const | |
364 | { | ||
365 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 130691930 times.
|
130691930 | NumEqVector flux(0.0); |
366 | 130691930 | const auto& scvf = this->scvFace(); | |
367 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 130691930 times.
|
130691930 | assert(scvf.isLateral()); |
368 | |||
369 | 130691930 | const auto& problem = this->problem(); | |
370 | 130691930 | const auto& elemVolVars = this->elemVolVars(); | |
371 | 130691930 | const auto fvGeometry = this->fvGeometry(); | |
372 | |||
373 | // get the transporting velocity which is perpendicular to our own (inner) velocity | ||
374 | 261383860 | const Scalar transportingVelocity = [&]() | |
375 | { | ||
376 | 130691930 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf); | |
377 |
2/2✓ Branch 2 taken 22 times.
✓ Branch 3 taken 130691908 times.
|
261383860 | const Scalar innerTransportingVelocity = elemVolVars[orthogonalScvf.insideScvIdx()].velocity(); |
378 | |||
379 |
4/6✓ Branch 0 taken 22 times.
✓ Branch 1 taken 130691908 times.
✓ Branch 3 taken 22 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 22 times.
✗ Branch 7 not taken.
|
130691930 | static const bool useOldScheme = getParam<bool>("FreeFlow.UseOldTransportingVelocity", true); // TODO how to deprecate? |
380 | |||
381 |
1/2✓ Branch 0 taken 130691930 times.
✗ Branch 1 not taken.
|
130691930 | if (useOldScheme) |
382 | { | ||
383 |
5/8✓ Branch 0 taken 1619662 times.
✓ Branch 1 taken 129072268 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1544510 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1619662 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 1544510 times.
|
130691930 | if (scvf.boundary() && fvGeometry.scv(scvf.insideScvIdx()).boundary()) |
384 | { | ||
385 | ✗ | if (this->elemBcTypes()[scvf.localIndex()].isDirichlet(scvf.normalAxis())) | |
386 | ✗ | return problem.dirichlet(this->element(), scvf)[scvf.normalAxis()]; | |
387 | } | ||
388 | else | ||
389 | 3775024 | return innerTransportingVelocity; | |
390 | } | ||
391 | |||
392 | // use the Dirichlet velocity as transporting velocity if the lateral face is on a Dirichlet boundary | ||
393 | ✗ | if (scvf.boundary()) | |
394 | { | ||
395 | ✗ | if (this->elemBcTypes()[scvf.localIndex()].isDirichlet(scvf.normalAxis())) | |
396 | ✗ | return 0.5*(problem.dirichlet(this->element(), scvf)[scvf.normalAxis()] + innerTransportingVelocity); | |
397 | } | ||
398 | |||
399 | ✗ | if (orthogonalScvf.boundary()) | |
400 | { | ||
401 | ✗ | if (this->elemBcTypes()[orthogonalScvf.localIndex()].isDirichlet(scvf.normalAxis())) | |
402 | ✗ | return 0.5*(problem.dirichlet(this->element(), scvf)[scvf.normalAxis()] + innerTransportingVelocity); | |
403 | else | ||
404 | return innerTransportingVelocity; // fallback value, should actually never be called | ||
405 | } | ||
406 | |||
407 | // average the transporting velocity by weighting with the scv volumes | ||
408 | ✗ | const auto insideVolume = fvGeometry.scv(orthogonalScvf.insideScvIdx()).volume(); | |
409 | ✗ | const auto outsideVolume = fvGeometry.scv(orthogonalScvf.outsideScvIdx()).volume(); | |
410 | ✗ | const auto outerTransportingVelocity = elemVolVars[orthogonalScvf.outsideScvIdx()].velocity(); | |
411 | ✗ | return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume); | |
412 | 130691930 | }(); | |
413 | |||
414 | 130691930 | const Scalar transportedMomentum = [&]() | |
415 | { | ||
416 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 126916906 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 126916906 times.
|
261383860 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); |
417 | |||
418 | 269922118 | auto getDirichletMomentumFlux = [&]() | |
419 | { | ||
420 |
0/8✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
2342718 | return problem.dirichlet(this->element(), scvf)[insideScv.dofAxis()] * this->problem().density(this->element(), insideScv); |
421 | }; | ||
422 | |||
423 | // use the Dirichlet velocity as for transported momentum if the lateral face is on a Dirichlet boundary | ||
424 | 130691930 | if (scvf.boundary()) | |
425 | { | ||
426 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 1619662 times.
|
1619662 | if (!this->elemBcTypes()[scvf.localIndex()].isDirichlet(insideScv.dofAxis())) |
427 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Neither Dirichlet nor Neumann BC set at " << scvf.ipGlobal()); | |
428 | |||
429 | 1619662 | return getDirichletMomentumFlux(); | |
430 | } | ||
431 | else | ||
432 | { | ||
433 |
2/3✗ Branch 0 not taken.
✓ Branch 1 taken 127956616 times.
✓ Branch 2 taken 1115652 times.
|
129072268 | if (fvGeometry.scvfIntegrationPointInConcaveCorner(scvf)) |
434 | { | ||
435 | // TODO: we could put this into an assert, as the construction of outsideScvfWithSameIntegrationPoint is quite expensive | ||
436 | 15488 | const auto& outsideScvfWithSameIntegrationPoint = fvGeometry.outsideScvfWithSameIntegrationPoint(scvf); | |
437 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 15488 times.
|
15488 | if (!this->problem().boundaryTypes(this->element(), outsideScvfWithSameIntegrationPoint).isDirichlet(insideScv.dofAxis())) |
438 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Neither Dirichlet nor Neumann BC set at " << outsideScvfWithSameIntegrationPoint.ipGlobal()); | |
439 | |||
440 | 15488 | return getDirichletMomentumFlux(); | |
441 | } | ||
442 | } | ||
443 | |||
444 | 129056780 | const bool selfIsUpstream = scvf.directionSign() == sign(transportingVelocity); | |
445 | |||
446 | 258113560 | const auto innerVelocity = elemVolVars[scvf.insideScvIdx()].velocity(); | |
447 |
0/2✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
258113560 | const auto outerVelocity = elemVolVars[scvf.outsideScvIdx()].velocity(); |
448 | |||
449 | 387170340 | const auto rho = this->problem().insideAndOutsideDensity(this->element(), fvGeometry, scvf); | |
450 | |||
451 | 129056780 | const auto insideMomentum = innerVelocity * rho.first; | |
452 | 129056780 | const auto outsideMomentum = outerVelocity * rho.second; | |
453 | |||
454 | // TODO use higher order helper | ||
455 |
5/8✓ Branch 0 taken 22 times.
✓ Branch 1 taken 129056758 times.
✓ Branch 3 taken 22 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 22 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 22 times.
✗ Branch 10 not taken.
|
129056780 | static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight"); |
456 | |||
457 |
2/2✓ Branch 0 taken 57225834 times.
✓ Branch 1 taken 71830946 times.
|
129056780 | return selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum) |
458 | 71830946 | : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum); | |
459 |
2/5✗ Branch 0 not taken.
✓ Branch 1 taken 1619662 times.
✓ Branch 2 taken 129072268 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
390456128 | }(); |
460 | |||
461 | 261383860 | return transportingVelocity * transportedMomentum * scvf.directionSign() * Extrusion::area(fvGeometry, scvf) * extrusionFactor_(elemVolVars, scvf); | |
462 | } | ||
463 | |||
464 | private: | ||
465 | |||
466 | template<class ElementVolumeVariables, class SubControlVolumeFace> | ||
467 | ✗ | Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf) const | |
468 | { | ||
469 | ✗ | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; | |
470 | ✗ | const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()]; | |
471 | ✗ | return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor()); | |
472 | } | ||
473 | |||
474 | |||
475 | const Problem* problemPtr_; //!< Pointer to the problem | ||
476 | const Element* elementPtr_; //!< Pointer to the element at hand | ||
477 | const FVElementGeometry* fvGeometryPtr_; //!< Pointer to the current FVElementGeometry | ||
478 | const SubControlVolumeFace* scvFacePtr_; //!< Pointer to the sub control volume face for which the flux variables are created | ||
479 | const ElementVolumeVariables* elemVolVarsPtr_; //!< Pointer to the current element volume variables | ||
480 | const ElementFluxVariablesCache* elemFluxVarsCachePtr_; //!< Pointer to the current element flux variables cache | ||
481 | const ElementBoundaryTypes* elemBcTypesPtr_; //!< Pointer to element boundary types | ||
482 | }; | ||
483 | |||
484 | } // end namespace Dumux | ||
485 | |||
486 | #endif // DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH | ||
487 |