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