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::NavierStokesResidualImpl | ||
11 | */ | ||
12 | #ifndef DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH | ||
13 | #define DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH | ||
14 | |||
15 | #include <dumux/common/numeqvector.hh> | ||
16 | #include <dumux/common/properties.hh> | ||
17 | #include <dumux/discretization/extrusion.hh> | ||
18 | #include <dumux/discretization/method.hh> | ||
19 | #include <dumux/assembly/fclocalresidual.hh> | ||
20 | #include <dune/common/hybridutilities.hh> | ||
21 | |||
22 | namespace Dumux { | ||
23 | |||
24 | /*! | ||
25 | * \ingroup NavierStokesModel | ||
26 | * \brief Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization | ||
27 | */ | ||
28 | template<class TypeTag> | ||
29 | class NavierStokesMomentumResidual | ||
30 | : public FaceCenteredLocalResidual<TypeTag> | ||
31 | { | ||
32 | using ParentType = FaceCenteredLocalResidual<TypeTag>; | ||
33 | friend class FaceCenteredLocalResidual<TypeTag>; | ||
34 | |||
35 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
36 | |||
37 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
38 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
39 | using VolumeVariables = typename GridVolumeVariables::VolumeVariables; | ||
40 | |||
41 | using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache; | ||
42 | using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView; | ||
43 | |||
44 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
45 | using Implementation = GetPropType<TypeTag, Properties::LocalResidual>; | ||
46 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
47 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
48 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
49 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
50 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
51 | using GridView = typename GridGeometry::GridView; | ||
52 | using Element = typename GridView::template Codim<0>::Entity; | ||
53 | using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; | ||
54 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
55 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
56 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
57 | |||
58 | using Extrusion = Extrusion_t<GridGeometry>; | ||
59 | |||
60 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
61 | |||
62 | static constexpr auto dim = GridView::dimension; | ||
63 | |||
64 | public: | ||
65 | //! Use the parent type's constructor | ||
66 |
2/4✓ Branch 1 taken 2574802 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2574802 times.
✗ Branch 5 not taken.
|
5211676 | using ParentType::ParentType; |
67 | |||
68 | /*! | ||
69 | * \brief Calculate the source term of the equation | ||
70 | * | ||
71 | * \param problem The problem to solve | ||
72 | * \param scv The sub-control volume over which we integrate the storage term | ||
73 | * \param volVars The volume variables associated with the scv | ||
74 | * \param isPreviousStorage Bool transferring the information if the storage term is computed at the current or previous time step | ||
75 | * \note has to be implemented by the model specific residual class | ||
76 | * | ||
77 | */ | ||
78 | 104823040 | NumEqVector computeStorage(const Problem& problem, | |
79 | const SubControlVolume& scv, | ||
80 | const VolumeVariables& volVars, | ||
81 | const bool isPreviousStorage = false) const | ||
82 | { | ||
83 | 209646080 | const auto& element = problem.gridGeometry().element(scv.elementIndex()); | |
84 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
104823040 | return problem.density(element, scv, isPreviousStorage) * volVars.velocity(); |
85 | } | ||
86 | |||
87 | /*! | ||
88 | * \brief Calculate the source term of the equation | ||
89 | * | ||
90 | * \param problem The problem to solve | ||
91 | * \param element The DUNE Codim<0> entity for which the residual | ||
92 | * ought to be calculated | ||
93 | * \param fvGeometry The finite-volume geometry of the element | ||
94 | * \param elemVolVars The volume variables associated with the element stencil | ||
95 | * \param scv The sub-control volume over which we integrate the source term | ||
96 | * \note This is the default implementation for all models as sources are computed | ||
97 | * in the user interface of the problem | ||
98 | * | ||
99 | */ | ||
100 | 43559688 | NumEqVector computeSource(const Problem& problem, | |
101 | const Element& element, | ||
102 | const FVElementGeometry& fvGeometry, | ||
103 | const ElementVolumeVariables& elemVolVars, | ||
104 | const SubControlVolume& scv) const | ||
105 | { | ||
106 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
43790472 | NumEqVector source = ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv); |
107 |
2/7✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1376760 times.
✓ Branch 4 taken 1376760 times.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
130679064 | source += problem.gravity()[scv.dofAxis()] * problem.density(element, scv); |
108 | |||
109 | // Axisymmetric problems in 2D feature an extra source terms arising from the transformation to cylindrical coordinates. | ||
110 | // See Ferziger/Peric: Computational methods for fluid dynamics chapter 8. | ||
111 | // https://doi.org/10.1007/978-3-540-68228-8 (page 301) | ||
112 | if constexpr (dim == 2 && isRotationalExtrusion<Extrusion>) | ||
113 | { | ||
114 |
2/2✓ Branch 0 taken 1376760 times.
✓ Branch 1 taken 1376760 times.
|
2753520 | if (scv.dofAxis() == Extrusion::radialAxis) |
115 | { | ||
116 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 1376760 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1376760 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1376760 times.
|
4130280 | const auto r = scv.center()[scv.dofAxis()] - fvGeometry.gridGeometry().bBoxMin()[scv.dofAxis()]; |
117 | 1376760 | const auto& scvf = (*scvfs(fvGeometry, scv).begin()); // the frontal scvf belonging to the scv | |
118 | |||
119 | // Velocity term | ||
120 | 1376760 | source -= -2.0*problem.effectiveViscosity(element, fvGeometry, scvf) * elemVolVars[scv].velocity() / (r*r); | |
121 | |||
122 | // Pressure term (needed because we incorporate pressure in terms of a surface integral). | ||
123 | // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term | ||
124 | // is new with respect to Cartesian coordinates and handled below as a source term. | ||
125 | 1376760 | source += problem.pressure(element, fvGeometry, scvf)/r; | |
126 | } | ||
127 | } | ||
128 | |||
129 | 43559688 | return source; | |
130 | } | ||
131 | |||
132 | /*! | ||
133 | * \brief Evaluates the momentum flux over a face of a sub control volume. | ||
134 | * | ||
135 | * \param problem The problem | ||
136 | * \param element The element | ||
137 | * \param fvGeometry The finite volume geometry context | ||
138 | * \param elemVolVars The volume variables for all flux stencil elements | ||
139 | * \param scvf The sub control volume face to compute the flux on | ||
140 | * \param elemFluxVarsCache The cache related to flux computation | ||
141 | * \param elemBcTypes The element boundary condition types | ||
142 | */ | ||
143 | 279435948 | NumEqVector computeFlux(const Problem& problem, | |
144 | const Element& element, | ||
145 | const FVElementGeometry& fvGeometry, | ||
146 | const ElementVolumeVariables& elemVolVars, | ||
147 | const SubControlVolumeFace& scvf, | ||
148 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
149 | const ElementBoundaryTypes& elemBcTypes) const | ||
150 | { | ||
151 | 279435948 | FluxVariables fluxVars(problem, element, fvGeometry, scvf, elemVolVars, elemFluxVarsCache, elemBcTypes); | |
152 | |||
153 | 279435948 | NumEqVector flux(0.0); | |
154 | 279435948 | flux += fluxVars.advectiveMomentumFlux(); | |
155 | 279435948 | flux += fluxVars.diffusiveMomentumFlux(); | |
156 | 279435948 | flux += fluxVars.pressureContribution(); | |
157 | 279435948 | return flux; | |
158 | } | ||
159 | |||
160 | //! Evaluate flux residuals for one sub control volume face when the element features at least one Dirichlet boundary condition | ||
161 | //! Skip the flux calculation if the DOF of the associated sub control volume features an appropriate Dirichlet condition. | ||
162 | 17156648 | NumEqVector maybeHandleDirichletBoundary(const Problem& problem, | |
163 | const Element& element, | ||
164 | const FVElementGeometry& fvGeometry, | ||
165 | const ElementVolumeVariables& elemVolVars, | ||
166 | const ElementBoundaryTypes& elemBcTypes, | ||
167 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
168 | const SubControlVolumeFace& scvf) const | ||
169 | { | ||
170 |
5/6✗ Branch 0 not taken.
✓ Branch 1 taken 16553428 times.
✓ Branch 2 taken 193984 times.
✓ Branch 3 taken 16962664 times.
✓ Branch 4 taken 4997468 times.
✓ Branch 5 taken 11555960 times.
|
34313296 | if (const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); scv.boundary()) |
171 | { | ||
172 | 5191452 | const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv); | |
173 | 5191452 | const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()]; | |
174 |
2/2✓ Branch 0 taken 110048 times.
✓ Branch 1 taken 5081404 times.
|
5191452 | if (bcTypes.isDirichlet(scv.dofAxis())) |
175 | 5081404 | return NumEqVector(0.0); // skip calculation as Dirichlet BC will be incorporated explicitly by assembler. | |
176 | } | ||
177 | |||
178 | // Default: The scv does not lie on a boundary with a Dirichlet condition for the velocity in normal direction. | ||
179 | // Check for a Neumann condition or calculate the flux as normal. | ||
180 |
2/2✓ Branch 0 taken 667494 times.
✓ Branch 1 taken 11407750 times.
|
12075244 | if (elemBcTypes.hasNeumann()) |
181 | 667494 | return this->asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf); | |
182 | else | ||
183 | 11407750 | return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes); | |
184 | } | ||
185 | |||
186 | //! Evaluate flux residuals for one sub control volume face when the element features at least one Neumann boundary condition | ||
187 | //! This requires special care. | ||
188 | 3478292 | NumEqVector maybeHandleNeumannBoundary(const Problem& problem, | |
189 | const Element& element, | ||
190 | const FVElementGeometry& fvGeometry, | ||
191 | const ElementVolumeVariables& elemVolVars, | ||
192 | const ElementBoundaryTypes& elemBcTypes, | ||
193 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
194 | const SubControlVolumeFace& scvf) const | ||
195 | { | ||
196 | static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered); // TODO overload this method for different discretizations | ||
197 | static_assert( | ||
198 | std::decay_t<decltype( | ||
199 | problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf) | ||
200 | )>::size() == ModelTraits::dim(), | ||
201 | "The momentum model expects problem.neumann to return a vector of size dim. " | ||
202 | "When in doubt you should be able to use 'using NumEqVector = typename ParentType::NumEqVector;'." | ||
203 | ); | ||
204 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 3478292 times.
|
3478292 | assert(elemBcTypes.hasNeumann()); |
205 | |||
206 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 3370272 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 3370272 times.
|
6956584 | const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); |
207 | |||
208 |
2/2✓ Branch 0 taken 819120 times.
✓ Branch 1 taken 2659172 times.
|
3478292 | if (scv.boundary()) |
209 | { | ||
210 | // This is the most simple case for a frontal scvf lying directly on the boundary. | ||
211 | // Just evaluate the Neumann flux for the normal velocity component. | ||
212 | // | ||
213 | // ---------------- * || frontal face of staggered half-control-volume | ||
214 | // | || # * | ||
215 | // | || # * D dof position (coincides with integration point where | ||
216 | // | || scv D * the Neumann flux is evaluated, here for x component) | ||
217 | // | || # * | ||
218 | // | || # * -- element boundaries | ||
219 | // ---------------- * | ||
220 | // y * domain boundary | ||
221 | // ^ | ||
222 | // | # current scvf over which Neumann flux is evaluated | ||
223 | // ----> x | ||
224 | // | ||
225 |
4/4✓ Branch 0 taken 208894 times.
✓ Branch 1 taken 610226 times.
✓ Branch 2 taken 186614 times.
✓ Branch 3 taken 22280 times.
|
819120 | if (scvf.boundary() && scvf.isFrontal()) |
226 | { | ||
227 | 186614 | const auto& bcTypes = elemBcTypes[scvf.localIndex()]; | |
228 |
3/6✓ Branch 0 taken 186614 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 186614 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 186614 times.
|
373228 | if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis())) |
229 | { | ||
230 | 186614 | const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf); | |
231 | 568466 | return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor(); | |
232 | } | ||
233 | } | ||
234 |
2/2✓ Branch 0 taken 445892 times.
✓ Branch 1 taken 186614 times.
|
632506 | else if (scvf.isLateral()) |
235 | { | ||
236 | // If a sub control volume lies on a boundary, Neumann fluxes also need to be considered for its lateral faces, | ||
237 | // even though those do not necessarily lie on a boundary (only their edges / integration points are touching the boundary). | ||
238 | // We cannot simply calculate momentum fluxes across the lateral boundary in this situation because we might not | ||
239 | // have enough information at the boundary to do so. | ||
240 | // We first need to find the scv's frontal face on the boundary to check if there actually is a Neumann BC. | ||
241 | // Afterwards, we use the actual (lateral) scvf to retrieve the Neumann flux at its integration point intersecting with the boundary. | ||
242 | // | ||
243 | // | ||
244 | // ---------######O * || frontal face of staggered half-control-volume | ||
245 | // | || : * | ||
246 | // | || : * D dof position | ||
247 | // | || scv D * | ||
248 | // | || : * -- element boundaries | ||
249 | // | || : * | ||
250 | // ---------------- * * domain boundary | ||
251 | // | ||
252 | // y # current scvf over which Neumann flux is evaluated | ||
253 | // ^ | ||
254 | // | : frontal scvf on boundary | ||
255 | // ----> x | ||
256 | // O integration point at which Neumann flux is evaluated (here for y component) | ||
257 | // | ||
258 | 445892 | const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv); | |
259 |
2/4✓ Branch 0 taken 445892 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 445892 times.
|
445892 | assert(frontalScvfOnBoundary.isFrontal() && frontalScvfOnBoundary.boundary()); |
260 | 445892 | const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()]; | |
261 |
4/6✓ Branch 0 taken 445892 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 445892 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 3584 times.
✓ Branch 5 taken 442308 times.
|
891784 | if (bcTypes.hasNeumann() && bcTypes.isNeumann(scvf.normalAxis())) |
262 | { | ||
263 | 442308 | const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf); | |
264 | 1326924 | return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor(); | |
265 | } | ||
266 | } | ||
267 | } | ||
268 |
2/2✓ Branch 0 taken 564220 times.
✓ Branch 1 taken 2094952 times.
|
2659172 | else if (scvf.boundary()) |
269 | { | ||
270 | // Here, the lateral face does lie on a boundary. Retrieve the Neumann flux component for | ||
271 | // the corresponding scvs' DOF orientation. | ||
272 | // | ||
273 | // | ||
274 | // ***************** | ||
275 | // ---------######O || frontal face of staggered half-control-volume | ||
276 | // | || : | ||
277 | // | || : D dof position | ||
278 | // | || scv D | ||
279 | // | || : -- element boundaries | ||
280 | // | || : | ||
281 | // ---------------- * domain boundary | ||
282 | // | ||
283 | // y # current scvf over which Neumann flux is evaluated | ||
284 | // ^ | ||
285 | // | : frontal scvf on boundary | ||
286 | // ----> x | ||
287 | // O integration point at which Neumann flux is evaluated (here for x component) | ||
288 | // | ||
289 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 564220 times.
|
564220 | assert(scvf.isLateral()); |
290 | 564220 | const auto& bcTypes = elemBcTypes[scvf.localIndex()]; | |
291 |
6/6✓ Branch 0 taken 503864 times.
✓ Branch 1 taken 60356 times.
✓ Branch 2 taken 503864 times.
✓ Branch 3 taken 60356 times.
✓ Branch 4 taken 3920 times.
✓ Branch 5 taken 499944 times.
|
1128440 | if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis())) |
292 | { | ||
293 | 499944 | const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf); | |
294 | 1499832 | return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor(); | |
295 | } | ||
296 | } | ||
297 | |||
298 | // Default: Neither the scvf itself nor a lateral one considers a Neumann flux. We just calculate the flux as normal. | ||
299 | 2349426 | return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes); | |
300 | } | ||
301 | |||
302 | //! Returns the implementation of the problem (i.e. static polymorphism) | ||
303 | Implementation &asImp_() | ||
304 | { return *static_cast<Implementation *>(this); } | ||
305 | |||
306 | //! \copydoc asImp_() | ||
307 | const Implementation &asImp_() const | ||
308 | { return *static_cast<const Implementation *>(this); } | ||
309 | }; | ||
310 | |||
311 | } // end namespace Dumux | ||
312 | |||
313 | #endif | ||
314 |