3.6-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
40template<class TypeTag>
42: public FaceCenteredLocalResidual<TypeTag>
43{
45 friend class FaceCenteredLocalResidual<TypeTag>;
46
48
49 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
50 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
51 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
52
53 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
54 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
55
60 using FVElementGeometry = typename GridGeometry::LocalView;
61 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
62 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
63 using GridView = typename GridGeometry::GridView;
64 using Element = typename GridView::template Codim<0>::Entity;
69
70 using Extrusion = Extrusion_t<GridGeometry>;
71
73
74 static constexpr auto dim = GridView::dimension;
75
76public:
78 using ParentType::ParentType;
79
90 NumEqVector computeStorage(const Problem& problem,
91 const SubControlVolume& scv,
92 const VolumeVariables& volVars,
93 const bool isPreviousStorage = false) const
94 {
95 const auto& element = problem.gridGeometry().element(scv.elementIndex());
96 return problem.density(element, scv, isPreviousStorage) * volVars.velocity();
97 }
98
112 NumEqVector computeSource(const Problem& problem,
113 const Element& element,
114 const FVElementGeometry& fvGeometry,
115 const ElementVolumeVariables& elemVolVars,
116 const SubControlVolume& scv) const
117 {
118 NumEqVector source = ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
119 source += problem.gravity()[scv.dofAxis()] * problem.density(element, scv);
120
121 // Axisymmetric problems in 2D feature an extra source terms arising from the transformation to cylindrical coordinates.
122 // See Ferziger/Peric: Computational methods for fluid dynamics chapter 8.
123 // https://doi.org/10.1007/978-3-540-68228-8 (page 301)
124 if constexpr (dim == 2 && isRotationalExtrusion<Extrusion>)
125 {
126 if (scv.dofAxis() == Extrusion::radialAxis)
127 {
128 const auto r = scv.center()[scv.dofAxis()] - fvGeometry.gridGeometry().bBoxMin()[scv.dofAxis()];
129 const auto& scvf = (*scvfs(fvGeometry, scv).begin()); // the frontal scvf belonging to the scv
130
131 // Velocity term
132 source -= -2.0*problem.effectiveViscosity(element, fvGeometry, scvf) * elemVolVars[scv].velocity() / (r*r);
133
134 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
135 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
136 // is new with respect to Cartesian coordinates and handled below as a source term.
137 source += problem.pressure(element, fvGeometry, scvf)/r;
138 }
139 }
140
141 return source;
142 }
143
155 NumEqVector computeFlux(const Problem& problem,
156 const Element& element,
157 const FVElementGeometry& fvGeometry,
158 const ElementVolumeVariables& elemVolVars,
159 const SubControlVolumeFace& scvf,
160 const ElementFluxVariablesCache& elemFluxVarsCache,
161 const ElementBoundaryTypes& elemBcTypes) const
162 {
163 FluxVariables fluxVars(problem, element, fvGeometry, scvf, elemVolVars, elemFluxVarsCache, elemBcTypes);
164
165 NumEqVector flux(0.0);
166 flux += fluxVars.advectiveMomentumFlux();
167 flux += fluxVars.diffusiveMomentumFlux();
168 flux += fluxVars.pressureContribution();
169 return flux;
170 }
171
174 NumEqVector maybeHandleDirichletBoundary(const Problem& problem,
175 const Element& element,
176 const FVElementGeometry& fvGeometry,
177 const ElementVolumeVariables& elemVolVars,
178 const ElementBoundaryTypes& elemBcTypes,
179 const ElementFluxVariablesCache& elemFluxVarsCache,
180 const SubControlVolumeFace& scvf) const
181 {
182 if (const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); scv.boundary())
183 {
184 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
185 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
186 if (bcTypes.isDirichlet(scv.dofAxis()))
187 return NumEqVector(0.0); // skip calculation as Dirichlet BC will be incorporated explicitly by assembler.
188 }
189
190 // Default: The scv does not lie on a boundary with a Dirichlet condition for the velocity in normal direction.
191 // Check for a Neumann condition or calculate the flux as normal.
192 if (elemBcTypes.hasNeumann())
193 return this->asImp().maybeHandleNeumannBoundary(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
194 else
195 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
196 }
197
200 NumEqVector maybeHandleNeumannBoundary(const Problem& problem,
201 const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const ElementVolumeVariables& elemVolVars,
204 const ElementBoundaryTypes& elemBcTypes,
205 const ElementFluxVariablesCache& elemFluxVarsCache,
206 const SubControlVolumeFace& scvf) const
207 {
208 static_assert(FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::fcstaggered); // TODO overload this method for different discretizations
209 static_assert(
210 std::decay_t<decltype(
211 problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf)
212 )>::size() == ModelTraits::dim(),
213 "The momentum model expects problem.neumann to return a vector of size dim. "
214 "When in doubt you should be able to use 'using NumEqVector = typename ParentType::NumEqVector;'."
215 );
216 assert(elemBcTypes.hasNeumann());
217
218 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
219
220 if (scv.boundary())
221 {
222 // This is the most simple case for a frontal scvf lying directly on the boundary.
223 // Just evaluate the Neumann flux for the normal velocity component.
224 //
225 // ---------------- * || frontal face of staggered half-control-volume
226 // | || # *
227 // | || # * D dof position (coincides with integration point where
228 // | || scv D * the Neumann flux is evaluated, here for x component)
229 // | || # *
230 // | || # * -- element boundaries
231 // ---------------- *
232 // y * domain boundary
233 // ^
234 // | # current scvf over which Neumann flux is evaluated
235 // ----> x
236 //
237 if (scvf.boundary() && scvf.isFrontal())
238 {
239 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
240 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
241 {
242 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
243 return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor();
244 }
245 }
246 else if (scvf.isLateral())
247 {
248 // If a sub control volume lies on a boundary, Neumann fluxes also need to be considered for its lateral faces,
249 // even though those do not necessarily lie on a boundary (only their edges / integration points are touching the boundary).
250 // We cannot simply calculate momentum fluxes across the lateral boundary in this situation because we might not
251 // have enough information at the boundary to do so.
252 // We first need to find the scv's frontal face on the boundary to check if there actually is a Neumann BC.
253 // Afterwards, we use the actual (lateral) scvf to retrieve the Neumann flux at its integration point intersecting with the boundary.
254 //
255 //
256 // ---------######O * || frontal face of staggered half-control-volume
257 // | || : *
258 // | || : * D dof position
259 // | || scv D *
260 // | || : * -- element boundaries
261 // | || : *
262 // ---------------- * * domain boundary
263 //
264 // y # current scvf over which Neumann flux is evaluated
265 // ^
266 // | : frontal scvf on boundary
267 // ----> x
268 // O integration point at which Neumann flux is evaluated (here for y component)
269 //
270 const auto& frontalScvfOnBoundary = fvGeometry.frontalScvfOnBoundary(scv);
271 assert(frontalScvfOnBoundary.isFrontal() && frontalScvfOnBoundary.boundary());
272 const auto& bcTypes = elemBcTypes[frontalScvfOnBoundary.localIndex()];
273 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scvf.normalAxis()))
274 {
275 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
276 return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor();
277 }
278 }
279 }
280 else if (scvf.boundary())
281 {
282 // Here, the lateral face does lie on a boundary. Retrieve the Neumann flux component for
283 // the corresponding scvs' DOF orientation.
284 //
285 //
286 // *****************
287 // ---------######O || frontal face of staggered half-control-volume
288 // | || :
289 // | || : D dof position
290 // | || scv D
291 // | || : -- element boundaries
292 // | || :
293 // ---------------- * domain boundary
294 //
295 // y # current scvf over which Neumann flux is evaluated
296 // ^
297 // | : frontal scvf on boundary
298 // ----> x
299 // O integration point at which Neumann flux is evaluated (here for x component)
300 //
301 assert(scvf.isLateral());
302 const auto& bcTypes = elemBcTypes[scvf.localIndex()];
303 if (bcTypes.hasNeumann() && bcTypes.isNeumann(scv.dofAxis()))
304 {
305 const auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
306 return neumannFluxes[scv.dofAxis()] * Extrusion::area(fvGeometry, scvf) * elemVolVars[scv].extrusionFactor();
307 }
308 }
309
310 // Default: Neither the scvf itself nor a lateral one considers a Neumann flux. We just calculate the flux as normal.
311 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, elemBcTypes);
312 }
313
315 Implementation &asImp_()
316 { return *static_cast<Implementation *>(this); }
317
319 const Implementation &asImp_() const
320 { return *static_cast<const Implementation *>(this); }
321};
322
323} // end namespace Dumux
324
325#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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr FCStaggered fcstaggered
Definition: method.hh:140
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:499
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:484
Element-wise calculation of the Navier-Stokes residual for models using the staggered discretization.
Definition: freeflow/navierstokes/momentum/localresidual.hh:43
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:200
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:112
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:155
Implementation & asImp_()
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:315
const Implementation & asImp_() const
Returns the implementation of the problem (i.e. static polymorphism)
Definition: freeflow/navierstokes/momentum/localresidual.hh:319
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:174
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:90
Declares all properties used in Dumux.