3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/navierstokes/momentum/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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH
25#define DUMUX_NAVIERSTOKES_MOMENTUM_LOCAL_RESIDUAL_HH
26
32#include <dune/common/hybridutilities.hh>
33
34namespace Dumux {
35
36namespace Impl {
37template<class T>
38static constexpr bool isRotationalExtrusion = false;
39
40template<int radialAxis>
41static constexpr bool isRotationalExtrusion<RotationalExtrusion<radialAxis>> = true;
42} // end namespace Impl
43
48template<class TypeTag>
50: public FaceCenteredLocalResidual<TypeTag>
51{
53 friend class FaceCenteredLocalResidual<TypeTag>;
54
56
57 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
58 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
59 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
60
61 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
62 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
63
68 using FVElementGeometry = typename GridGeometry::LocalView;
69 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
70 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
71 using GridView = typename GridGeometry::GridView;
72 using Element = typename GridView::template Codim<0>::Entity;
77
78 using Extrusion = Extrusion_t<GridGeometry>;
79
81
82 static constexpr auto dim = GridView::dimension;
83
84public:
86 using ParentType::ParentType;
87
98 NumEqVector computeStorage(const Problem& problem,
99 const SubControlVolume& scv,
100 const VolumeVariables& volVars,
101 const bool isPreviousStorage = false) const
102 {
103 const auto& element = problem.gridGeometry().element(scv.elementIndex());
104 return problem.density(element, scv, isPreviousStorage) * volVars.velocity();
105 }
106
120 NumEqVector computeSource(const Problem& problem,
121 const Element& element,
122 const FVElementGeometry& fvGeometry,
123 const ElementVolumeVariables& elemVolVars,
124 const SubControlVolume& scv) const
125 {
126 NumEqVector source = ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
127 source += problem.gravity()[scv.dofAxis()] * problem.density(element, scv);
128
129 // Axisymmetric problems in 2D feature an extra source terms arising from the transformation to cylindrical coordinates.
130 // See Ferziger/Peric: Computational methods for fluid dynamics chapter 8.
131 // https://doi.org/10.1007/978-3-540-68228-8 (page 301)
132 if constexpr (dim == 2 && Impl::isRotationalExtrusion<Extrusion>)
133 {
134 if (scv.dofAxis() == Extrusion::radialAxis)
135 {
136 const auto r = scv.center()[scv.dofAxis()] - fvGeometry.gridGeometry().bBoxMin()[scv.dofAxis()];
137 const auto& scvf = (*scvfs(fvGeometry, scv).begin()); // the frontal scvf belonging to the scv
138
139 // Velocity term
140 source -= -2.0*problem.effectiveViscosity(element, fvGeometry, scvf) * elemVolVars[scv].velocity() / (r*r);
141
142 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
143 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
144 // is new with respect to Cartesian coordinates and handled below as a source term.
145 source += problem.pressure(element, fvGeometry, scvf)/r;
146 }
147 }
148
149 return source;
150 }
151
163 NumEqVector computeFlux(const Problem& problem,
164 const Element& element,
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const SubControlVolumeFace& scvf,
168 const ElementFluxVariablesCache& elemFluxVarsCache,
169 const ElementBoundaryTypes& elemBcTypes) const
170 {
171 FluxVariables fluxVars(problem, element, fvGeometry, scvf, elemVolVars, elemFluxVarsCache, elemBcTypes);
172
173 NumEqVector flux(0.0);
174 flux += fluxVars.advectiveMomentumFlux();
175 flux += fluxVars.diffusiveMomentumFlux();
176 flux += fluxVars.pressureContribution();
177 return flux;
178 }
179
182 NumEqVector maybeHandleDirichletBoundary(const Problem& problem,
183 const Element& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars,
186 const ElementBoundaryTypes& elemBcTypes,
187 const ElementFluxVariablesCache& elemFluxVarsCache,
188 const SubControlVolumeFace& scvf) const
189 {
190 if (const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); scv.boundary())
191 {
192 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
193 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
194 if (bcTypes.isDirichlet(scv.dofAxis()))
195 return NumEqVector(0.0); // skip calculation as Dirichlet BC will be incorporated explicitly by assembler.
196 }
197
198 // Default: The scv does not lie on a boundary with a Dirichlet condition for the velocity in normal direction.
199 // Check for a Neumann condition or calculate the flux as normal.
200 if (elemBcTypes.hasNeumann())
201 return this->asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
202 else
203 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
204 }
205
208 NumEqVector maybeHandleNeumannBoundary(const Problem& problem,
209 const Element& element,
210 const FVElementGeometry& fvGeometry,
211 const ElementVolumeVariables& elemVolVars,
212 const ElementBoundaryTypes& elemBcTypes,
213 const ElementFluxVariablesCache& elemFluxVarsCache,
214 const SubControlVolumeFace& scvf) const
215 {
216 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered); // TODO overload this method for different discretizations
217 static_assert(
218 std::decay_t<decltype(
219 problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf)
220 )>::size() == ModelTraits::dim(),
221 "The momentum model expects problem.neumann to return a vector of size dim. "
222 "When in doubt you should be able to use 'using NumEqVector = typename ParentType::NumEqVector;'."
223 );
224 assert(elemBcTypes.hasNeumann());
225
226 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
227
228 if (scv.boundary())
229 {
230 // This is the most simple case for a frontal scvf lying directly on the boundary.
231 // Just evaluate the Neumann flux for the normal velocity component.
232 //
233 // ---------------- * || frontal face of staggered half-control-volume
234 // | || # *
235 // | || # * D dof position (coincides with integration point where
236 // | || scv D * the Neumann flux is evaluated, here for x component)
237 // | || # *
238 // | || # * -- element boundaries
239 // ---------------- *
240 // y * domain boundary
241 // ^
242 // | # current scvf over which Neumann flux is evaluated
243 // ----> x
244 //
245 if (scvf.boundary() && scvf.isFrontal())
246 {
247 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
248 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
249 {
250 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
251 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
252 }
253 }
254 else if (scvf.isLateral())
255 {
256 // If a sub control volume lies on a boundary, Neumann fluxes also need to be considered for its lateral faces,
257 // even though those do not necessarily lie on a boundary (only their edges / integration points are touching the boundary).
258 // We cannot simply calculate momentum fluxes across the lateral boundary in this situation because we might not
259 // have enough information at the boundary to do so.
260 // We first need to find the scv's frontal face on the boundary to check if there actually is a Neumann BC.
261 // Afterwards, we use the actual (lateral) scvf to retrieve the Neumann flux at its integration point intersecting with the boundary.
262 //
263 //
264 // ---------######O * || frontal face of staggered half-control-volume
265 // | || : *
266 // | || : * D dof position
267 // | || scv D *
268 // | || : * -- element boundaries
269 // | || : *
270 // ---------------- * * domain boundary
271 //
272 // y # current scvf over which Neumann flux is evaluated
273 // ^
274 // | : frontal scvf on boundary
275 // ----> x
276 // O integration point at which Neumann flux is evaluated (here for y component)
277 //
278 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
279 assert(frontalScvfOnBoundary.isFrontal() && frontalScvfOnBoundary.boundary());
280 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
281 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scvf.normalAxis()))
282 {
283 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
284 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
285 }
286 }
287 }
288 else if (scvf.boundary())
289 {
290 // Here, the lateral face does lie on a boundary. Retrive the Neumann flux component for
291 // the corresponding scvs' DOF orientation.
292 //
293 //
294 // *****************
295 // ---------######O || frontal face of staggered half-control-volume
296 // | || :
297 // | || : D dof position
298 // | || scv D
299 // | || : -- element boundaries
300 // | || :
301 // ---------------- * domain boundary
302 //
303 // y # current scvf over which Neumann flux is evaluated
304 // ^
305 // | : frontal scvf on boundary
306 // ----> x
307 // O integration point at which Neumann flux is evaluated (here for x component)
308 //
309 assert(scvf.isLateral());
310 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
311 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
312 {
313 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
314 return neumannFluxes[scv.dofAxis()] * Extrusion::area(scvf) * elemVolVars[scv].extrusionFactor();
315 }
316 }
317
318 // Default: Neither the scvf itself nor a lateral one considers a Neumann flux. We just calculate the flux as normal.
319 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
320 }
321
323 Implementation &asImp_()
324 { return *static_cast<Implementation *>(this); }
325
327 const Implementation &asImp_() const
328 { return *static_cast<const Implementation *>(this); }
329};
330
331} // end namespace Dumux
332
333#endif
Calculates the element-wise residual for the box scheme.
A helper to deduce a vector with the same size as numbers of equations.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
constexpr FCStaggered fcstaggered
Definition: method.hh:142
The element-wise residual for the box scheme.
Definition: fclocalresidual.hh:46
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: fclocalresidual.hh:112
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:501
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:486
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/momentum/localresidual.hh:51
NumEqVector maybeHandleNeumannBoundary(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition: freeflow/navierstokes/momentum/localresidual.hh:208
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: freeflow/navierstokes/momentum/localresidual.hh:120
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes) const
Evaluates the mass flux over a face of a sub control volume.
Definition: freeflow/navierstokes/momentum/localresidual.hh:163
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:323
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:327
NumEqVector maybeHandleDirichletBoundary(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition: freeflow/navierstokes/momentum/localresidual.hh:182
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars, const bool isPreviousStorage=false) const
Calculate the source term of the equation.
Definition: freeflow/navierstokes/momentum/localresidual.hh:98
Declares all properties used in Dumux.