version 3.8
freeflow/navierstokes/staggered/localresidual.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_STAGGERED_NAVIERSTOKES_LOCAL_RESIDUAL_HH
13#define DUMUX_STAGGERED_NAVIERSTOKES_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
16
22
23namespace Dumux {
24
25namespace Impl {
26template<class T>
27static constexpr bool isRotationalExtrusion = false;
28
29template<int radialAxis>
30static constexpr bool isRotationalExtrusion<RotationalExtrusion<radialAxis>> = true;
31} // end namespace Impl
32
38// forward declaration
39template<class TypeTag, class DiscretizationMethod>
41
42template<class TypeTag>
43class NavierStokesResidualImpl<TypeTag, DiscretizationMethods::Staggered>
44: public StaggeredLocalResidual<TypeTag>
45{
47 friend class StaggeredLocalResidual<TypeTag>;
48
50
51 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
52 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
53 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
54
55 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
56 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
57
58 using GridFaceVariables = typename GridVariables::GridFaceVariables;
59 using ElementFaceVariables = typename GridFaceVariables::LocalView;
60
65 using FVElementGeometry = typename GridGeometry::LocalView;
66 using GridView = typename GridGeometry::GridView;
67 using Element = typename GridView::template Codim<0>::Entity;
68 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
69 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
70 using Extrusion = Extrusion_t<GridGeometry>;
72 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
76
77 static constexpr bool normalizePressure = getPropValue<TypeTag, Properties::NormalizePressure>();
78
79 using CellCenterResidual = CellCenterPrimaryVariables;
80 using FaceResidual = FacePrimaryVariables;
81
83
84public:
85 using EnergyLocalResidual = FreeFlowEnergyLocalResidual<GridGeometry, FluxVariables, ModelTraits::enableEnergyBalance(), (ModelTraits::numFluidComponents() > 1)>;
86
87 // account for the offset of the cell center privars within the PrimaryVariables container
88 static constexpr auto cellCenterOffset = ModelTraits::numEq() - CellCenterPrimaryVariables::dimension;
89 static_assert(cellCenterOffset == ModelTraits::dim(), "cellCenterOffset must equal dim for staggered NavierStokes");
90
92 using ParentType::ParentType;
93
95 CellCenterPrimaryVariables computeFluxForCellCenter(const Problem& problem,
96 const Element &element,
97 const FVElementGeometry& fvGeometry,
98 const ElementVolumeVariables& elemVolVars,
99 const ElementFaceVariables& elemFaceVars,
100 const SubControlVolumeFace &scvf,
101 const ElementFluxVariablesCache& elemFluxVarsCache) const
102 {
103 FluxVariables fluxVars;
104 CellCenterPrimaryVariables flux = fluxVars.computeMassFlux(problem, element, fvGeometry, elemVolVars,
105 elemFaceVars, scvf, elemFluxVarsCache[scvf]);
106
107 EnergyLocalResidual::heatFlux(flux, problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf);
108
109 return flux;
110 }
111
113 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
114 const Element &element,
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& elemVolVars,
117 const ElementFaceVariables& elemFaceVars,
118 const SubControlVolume &scv) const
119 {
120 CellCenterPrimaryVariables result(0.0);
121
122 // get the values from the problem
123 const auto sourceValues = problem.source(element, fvGeometry, elemVolVars, elemFaceVars, scv);
124
125 // copy the respective cell center related values to the result
126 for (int i = 0; i < result.size(); ++i)
127 result[i] = sourceValues[i + cellCenterOffset];
128
129 return result;
130 }
131
132
134 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
135 const SubControlVolume& scv,
136 const VolumeVariables& volVars) const
137 {
138 CellCenterPrimaryVariables storage;
139 storage[Indices::conti0EqIdx - ModelTraits::dim()] = volVars.density();
140
141 EnergyLocalResidual::fluidPhaseStorage(storage, volVars);
142
143 return storage;
144 }
145
147 FacePrimaryVariables computeStorageForFace(const Problem& problem,
148 const SubControlVolumeFace& scvf,
149 const VolumeVariables& volVars,
150 const ElementFaceVariables& elemFaceVars) const
151 {
152 FacePrimaryVariables storage(0.0);
153 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
154 storage[0] = volVars.density() * velocity;
155 return storage;
156 }
157
159 FacePrimaryVariables computeSourceForFace(const Problem& problem,
160 const Element& element,
161 const FVElementGeometry& fvGeometry,
162 const SubControlVolumeFace& scvf,
163 const ElementVolumeVariables& elemVolVars,
164 const ElementFaceVariables& elemFaceVars) const
165 {
166 FacePrimaryVariables source(0.0);
167 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
168 source += problem.gravity()[scvf.directionIndex()] * insideVolVars.density();
169 source += problem.source(element, fvGeometry, elemVolVars, elemFaceVars, scvf)[Indices::velocity(scvf.directionIndex())];
170
171 // Axisymmetric problems in 2D feature an extra source terms arising from the transformation to cylindrical coordinates.
172 // See Ferziger/Peric: Computational methods for fluid dynamics chapter 8.
173 // https://doi.org/10.1007/978-3-540-68228-8 (page 301)
174 if constexpr (ModelTraits::dim() == 2 && Impl::isRotationalExtrusion<Extrusion>)
175 {
176 if (scvf.directionIndex() == Extrusion::radialAxis)
177 {
178 // Velocity term
179 const auto r = scvf.center()[scvf.directionIndex()];
180 source -= -2.0*insideVolVars.effectiveViscosity() * elemFaceVars[scvf].velocitySelf() / (r*r);
181
182 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
183 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
184 // is new with respect to Cartesian coordinates and handled below as a source term.
185 const Scalar pressure =
186 normalizePressure ? insideVolVars.pressure() - problem.initial(scvf)[Indices::pressureIdx]
187 : insideVolVars.pressure();
188
189 source += pressure/r;
190 }
191 }
192
193 return source;
194 }
195
197 FacePrimaryVariables computeFluxForFace(const Problem& problem,
198 const Element& element,
199 const SubControlVolumeFace& scvf,
200 const FVElementGeometry& fvGeometry,
201 const ElementVolumeVariables& elemVolVars,
202 const ElementFaceVariables& elemFaceVars,
203 const ElementFluxVariablesCache& elemFluxVarsCache) const
204 {
205 FluxVariables fluxVars;
206 return fluxVars.computeMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache.gridFluxVarsCache());
207 }
208
212 CellCenterResidual computeBoundaryFluxForCellCenter(const Problem& problem,
213 const Element& element,
214 const FVElementGeometry& fvGeometry,
215 const SubControlVolumeFace& scvf,
216 const ElementVolumeVariables& elemVolVars,
217 const ElementFaceVariables& elemFaceVars,
218 const ElementBoundaryTypes& elemBcTypes,
219 const ElementFluxVariablesCache& elemFluxVarsCache) const
220 {
221 CellCenterResidual result(0.0);
222
223 if (scvf.boundary())
224 {
225 const auto bcTypes = problem.boundaryTypes(element, scvf);
226
227 // no fluxes occur over symmetry boundaries
228 if (bcTypes.isSymmetry())
229 return result;
230
231 // treat Dirichlet and outflow BCs
232 result = computeFluxForCellCenter(problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
233
234 // treat Neumann BCs, i.e. overwrite certain fluxes by user-specified values
235 static constexpr auto numEqCellCenter = CellCenterResidual::dimension;
236 if (bcTypes.hasNeumann())
237 {
238 const auto extrusionFactor = elemVolVars[scvf.insideScvIdx()].extrusionFactor();
239 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
240
241 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
242 {
243 if (bcTypes.isNeumann(eqIdx + cellCenterOffset))
244 result[eqIdx] = neumannFluxes[eqIdx + cellCenterOffset] * extrusionFactor * Extrusion::area(fvGeometry, scvf);
245 }
246 }
247
248 // account for wall functions, if used
249 incorporateWallFunction_(result, problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
250 }
251 return result;
252 }
253
257 void evalDirichletBoundariesForFace(FaceResidual& residual,
258 const Problem& problem,
259 const Element& element,
260 const FVElementGeometry& fvGeometry,
261 const SubControlVolumeFace& scvf,
262 const ElementVolumeVariables& elemVolVars,
263 const ElementFaceVariables& elemFaceVars,
264 const ElementBoundaryTypes& elemBcTypes,
265 const ElementFluxVariablesCache& elemFluxVarsCache) const
266 {
267 if (scvf.boundary())
268 {
269 // handle the actual boundary conditions:
270 const auto bcTypes = problem.boundaryTypes(element, scvf);
271
272 if(bcTypes.isDirichlet(Indices::velocity(scvf.directionIndex())))
273 {
274 // set a fixed value for the velocity for Dirichlet boundary conditions
275 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
276 const Scalar dirichletValue = problem.dirichlet(element, scvf)[Indices::velocity(scvf.directionIndex())];
277 residual = velocity - dirichletValue;
278 }
279 else if(bcTypes.isSymmetry())
280 {
281 // for symmetry boundary conditions, there is no flow across the boundary and
282 // we therefore treat it like a Dirichlet boundary conditions with zero velocity
283 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
284 const Scalar fixedValue = 0.0;
285 residual = velocity - fixedValue;
286 }
287 }
288 }
289
293 FaceResidual computeBoundaryFluxForFace(const Problem& problem,
294 const Element& element,
295 const FVElementGeometry& fvGeometry,
296 const SubControlVolumeFace& scvf,
297 const ElementVolumeVariables& elemVolVars,
298 const ElementFaceVariables& elemFaceVars,
299 const ElementBoundaryTypes& elemBcTypes,
300 const ElementFluxVariablesCache& elemFluxVarsCache) const
301 {
302 FaceResidual result(0.0);
303
304 if (scvf.boundary())
305 {
306 FluxVariables fluxVars;
307
308 // handle the actual boundary conditions:
309 const auto bcTypes = problem.boundaryTypes(element, scvf);
310 if (bcTypes.isNeumann(Indices::velocity(scvf.directionIndex())))
311 {
312 // the source term has already been accounted for, here we
313 // add a given Neumann flux for the face on the boundary itself ...
314 const auto extrusionFactor = elemVolVars[scvf.insideScvIdx()].extrusionFactor();
315 result = problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, scvf)[Indices::velocity(scvf.directionIndex())]
316 * extrusionFactor * Extrusion::area(fvGeometry, scvf);
317
318 // ... and treat the fluxes of the remaining (frontal and lateral) faces of the staggered control volume
319 result += fluxVars.computeMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache.gridFluxVarsCache());
320 }
321 else if(bcTypes.isDirichlet(Indices::pressureIdx))
322 {
323 // we are at an "fixed pressure" boundary for which the resdiual of the momentum balance needs to be assembled
324 // as if it where inside the domain and not on the boundary (source term has already been acounted for)
325 result = fluxVars.computeMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache.gridFluxVarsCache());
326
327 // incorporate the inflow or outflow contribution
328 result += fluxVars.inflowOutflowBoundaryFlux(problem, scvf, fvGeometry, elemVolVars, elemFaceVars);
329 }
330 }
331 return result;
332 }
333
334private:
335
337 template<class ...Args, bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<!turbulenceModel, int> = 0>
338 void incorporateWallFunction_(Args&&... args) const
339 {}
340
342 template<bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<turbulenceModel, int> = 0>
343 void incorporateWallFunction_(CellCenterResidual& boundaryFlux,
344 const Problem& problem,
345 const Element& element,
346 const FVElementGeometry& fvGeometry,
347 const SubControlVolumeFace& scvf,
348 const ElementVolumeVariables& elemVolVars,
349 const ElementFaceVariables& elemFaceVars) const
350 {
351 static constexpr auto numEqCellCenter = CellCenterResidual::dimension;
352 const auto extrusionFactor = elemVolVars[scvf.insideScvIdx()].extrusionFactor();
353
354 // account for wall functions, if used
355 for(int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
356 {
357 // use a wall function
358 if(problem.useWallFunction(element, scvf, eqIdx + cellCenterOffset))
359 {
360 boundaryFlux[eqIdx] = problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf)[eqIdx]
361 * extrusionFactor * Extrusion::area(fvGeometry, scvf);
362 }
363 }
364 }
365
367 Implementation &asImp_()
368 { return *static_cast<Implementation *>(this); }
369
371 const Implementation &asImp_() const
372 { return *static_cast<const Implementation *>(this); }
373};
374
375} // end namespace Dumux
376
377#endif
Definition: freeflow/nonisothermal/localresidual.hh:24
FacePrimaryVariables computeFluxForFace(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the momentum flux for the face control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:197
void evalDirichletBoundariesForFace(FaceResidual &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate Dirichlet (fixed value) boundary conditions for a face dof.
Definition: freeflow/navierstokes/staggered/localresidual.hh:257
FacePrimaryVariables computeStorageForFace(const Problem &problem, const SubControlVolumeFace &scvf, const VolumeVariables &volVars, const ElementFaceVariables &elemFaceVars) const
Evaluate the storage term for the face control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:147
FacePrimaryVariables computeSourceForFace(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluate the source term for the face control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:159
FaceResidual computeBoundaryFluxForFace(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate boundary fluxes for a face dof.
Definition: freeflow/navierstokes/staggered/localresidual.hh:293
CellCenterResidual computeBoundaryFluxForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate boundary conditions for a cell center dof.
Definition: freeflow/navierstokes/staggered/localresidual.hh:212
CellCenterPrimaryVariables computeFluxForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate fluxes entering or leaving the cell center control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:95
CellCenterPrimaryVariables computeStorageForCellCenter(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate the storage term for the cell center control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:134
CellCenterPrimaryVariables computeSourceForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolume &scv) const
Evaluate the source term for the cell center control volume.
Definition: freeflow/navierstokes/staggered/localresidual.hh:113
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/localresidual.hh:23
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:29
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
constexpr bool isRotationalExtrusion
Convenience trait to check whether the extrusion is rotational.
Definition: extrusion.hh:172
Calculates the element-wise residual for the staggered FV scheme.