version 3.11-dev
freeflow/navierstokes/momentum/cvfe/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_NAVIERSTOKES_MOMENTUM_CVFE_LOCAL_RESIDUAL_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
16#include <dune/geometry/quadraturerules.hh>
17
20#include <dumux/common/typetraits/localdofs_.hh>
22
29
32
33namespace Dumux {
34
35namespace Detail {
36
38template<class P, class FVG, class EV, class IPD>
40 std::declval<P>().source(std::declval<FVG>(), std::declval<EV>(), std::declval<IPD>())
41);
42
43template<class P, class FVG, class EV, class IPD>
45{ return Dune::Std::is_detected<SourceWithIpDataInterface, P, FVG, EV, IPD>::value; }
46
47} // end namespace Detail
48
53template<class TypeTag>
55: public CVFELocalResidual<TypeTag>
56{
58
60
61 using GridVariablesCache = typename GridVariables::GridVolumeVariables;
62 using ElementVariables = typename GridVariablesCache::LocalView;
63 using Variables = typename GridVariablesCache::VolumeVariables;
64
65 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
66 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
67
71 using FVElementGeometry = typename GridGeometry::LocalView;
72 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
73 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
74 using GridView = typename GridGeometry::GridView;
75 using Element = typename GridView::template Codim<0>::Entity;
79 using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
80
81 using Extrusion = Extrusion_t<GridGeometry>;
82
84
85 static constexpr auto dim = GridView::dimension;
86
87 using LocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
88 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
93
94 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
97
98public:
101 using ParentType::ParentType;
102
113 NumEqVector computeStorage(const Problem& problem,
114 const FVElementGeometry& fvGeometry,
115 const SubControlVolume& scv,
116 const Variables& volVars,
117 const bool isPreviousStorage) const
118 {
119 return problem.density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scv), isPreviousStorage) * volVars.velocity();
120 }
121
131 NumEqVector storageIntegral(const FVElementGeometry& fvGeometry,
132 const ElementVariables& elemVars,
133 const SubControlVolume& scv,
134 bool isPreviousTimeLevel) const
135 {
136 const auto& volVars = elemVars[scv];
137 // We apply mass lumping
138 NumEqVector storage = this->asImp().problem().density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scv), isPreviousTimeLevel)
139 * volVars.velocity();
140
141 storage *= Extrusion::volume(fvGeometry, scv) * volVars.extrusionFactor();
142
143 return storage;
144 }
145
157 NumEqVector computeSource(const Problem& problem,
158 const Element& element,
159 const FVElementGeometry& fvGeometry,
160 const ElementVariables& elemVars,
161 const SubControlVolume& scv) const
162 {
163 NumEqVector source;
164
165 if constexpr (Detail::hasProblemSourceWithIpDataInterface<Problem, FVElementGeometry, ElementVariables, BaseIpData>())
166 {
167 source = problem.source(fvGeometry, elemVars, ipData(fvGeometry, scv.center()));
168
169 // ToDo: point source data with ipData
170 // add contribution from possible point sources
171 if (!problem.pointSourceMap().empty())
172 source += problem.scvPointSources(element, fvGeometry, elemVars, scv);
173 }
174 else
175 source = ParentType::computeSource(problem, element, fvGeometry, elemVars, scv);
176
177
178 // add rho*g (note that gravity might be zero in case it's disabled in the problem)
179 const auto& data = ipData(fvGeometry, scv);
180 source += problem.density(element, fvGeometry, data) * problem.gravity();
181
182 // Axisymmetric problems in 2D feature an extra source term arising from the transformation to cylindrical coordinates.
183 // See Ferziger/Peric: Computational methods for Fluid Dynamics (2020)
184 // https://doi.org/10.1007/978-3-319-99693-6
185 // Chapter 9.9 and Eq. (9.81) and comment on finite volume methods
186 if constexpr (dim == 2 && isRotationalExtrusion<Extrusion>)
187 {
188 // the radius with respect to the rotation axis
189 const auto r = scv.center()[Extrusion::radialAxis] - fvGeometry.gridGeometry().bBoxMin()[Extrusion::radialAxis];
190
191 // The velocity term is new with respect to Cartesian coordinates and handled below as a source term
192 // It only enters the balance of the momentum balance in radial direction
193 source[Extrusion::radialAxis] += -2.0*problem.effectiveViscosity(element, fvGeometry, data)
194 * elemVars[scv].velocity(Extrusion::radialAxis) / (r*r);
195
196 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
197 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
198 // is new with respect to Cartesian coordinates and handled below as a source term.
199 source[Extrusion::radialAxis] += problem.pressure(element, fvGeometry, data)/r;
200 }
201
202 return source;
203 }
204
213 NumEqVector sourceIntegral(const FVElementGeometry& fvGeometry,
214 const ElementVariables& elemVars,
215 const SubControlVolume& scv) const
216 {
217 static_assert(!(dim == 2 && isRotationalExtrusion<Extrusion>), "Rotational extrusion source terms are not implemented for integral interface.");
218
219 const auto& problem = this->asImp().problem();
220
221 NumEqVector source(0.0);
222 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scv))
223 {
224 source += qpData.weight() * (problem.source(fvGeometry, elemVars, qpData.ipData())
225 + problem.density(fvGeometry.element(), fvGeometry, qpData.ipData()) * problem.gravity());
226 }
227
228 // ToDo: point source data with ipData
229 // add contribution from possible point sources
230 if (!problem.pointSourceMap().empty())
231 source += Extrusion::volume(fvGeometry, scv) * problem.scvPointSources(fvGeometry.element(), fvGeometry, elemVars, scv);
232
233 source *= elemVars[scv].extrusionFactor();
234
235 return source;
236 }
237
248 NumEqVector computeFlux(const Problem& problem,
249 const Element& element,
250 const FVElementGeometry& fvGeometry,
251 const ElementVariables& elemVars,
252 const SubControlVolumeFace& scvf,
253 const ElementFluxVariablesCache& elemFluxVarsCache) const
254 {
255 FluxContext context(problem, fvGeometry, elemVars, elemFluxVarsCache, scvf);
256 FluxHelper fluxHelper;
257
258 NumEqVector flux(0.0);
259 flux += fluxHelper.advectiveMomentumFlux(context);
260 flux += fluxHelper.diffusiveMomentumFlux(context);
261 flux += fluxHelper.pressureContribution(context);
262 return flux;
263 }
264
274 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
275 const ElementVariables& elemVars,
276 const ElementFluxVariablesCache& elemFluxVarsCache,
277 const SubControlVolumeFace& scvf) const
278 {
279 const auto& problem = this->asImp().problem();
280
281 NumEqVector flux(0.0);
282 GlobalPosition velIntegral(0.0);
283 FluxFunctionHelper fluxFunctionHelper;
284
285 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
286 {
287 const auto& fluxVarsCache = elemFluxVarsCache[qpData.ipData()];
288 FluxFunctionContext context(this->problem(), fvGeometry, elemVars, fluxVarsCache);
289
290 velIntegral += context.velocity() * qpData.weight();
291 flux += qpData.weight() * ( fluxFunctionHelper.diffusiveMomentumFluxIntegrand(context, qpData.ipData())
292 + fluxFunctionHelper.pressureFluxIntegrand(context, qpData.ipData()) );
293 }
294 flux += fluxFunctionHelper.advectiveMomentumFluxIntegral(problem, fvGeometry, elemVars, elemFluxVarsCache, scvf, velIntegral);
295
296 flux *= elemVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
297
298 return flux;
299 }
300
302 const Problem& problem,
303 const Element& element,
304 const FVElementGeometry& fvGeometry,
305 const ElementVariables& prevElemVolVars,
306 const ElementVariables& curElemVolVars) const
307 {
309 residual, problem, fvGeometry, prevElemVolVars, curElemVolVars, this->timeLoop().timeStepSize()
310 );
311 }
312
314 const Problem& problem,
315 const Element& element,
316 const FVElementGeometry& fvGeometry,
317 const ElementVariables& elemVars,
318 const ElementFluxVariablesCache& elemFluxVarsCache,
319 const ElementBoundaryTypes &elemBcTypes) const
320 {
322 residual, problem, fvGeometry, elemVars, elemFluxVarsCache, elemBcTypes
323 );
324 }
325};
326
327} // end namespace Dumux
328
329#endif
Boundary flag to store e.g. in sub control volume faces.
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:77
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cvfelocalresidual.hh:98
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:36
Implementation & asImp()
Definition: fvlocalresidual.hh:489
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: fvlocalresidual.hh:209
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:474
const TimeLoop & timeLoop() const
Definition: fvlocalresidual.hh:479
Element-wise calculation of the Navier-Stokes residual for models using CVFE discretizations.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:56
NumEqVector sourceIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source integral.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:213
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source term of the equation.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:157
void addToElementFluxAndSourceResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes) const
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:313
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the mass flux over a face of a sub control volume.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:248
void addToElementStorageResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &prevElemVolVars, const ElementVariables &curElemVolVars) const
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:301
NumEqVector storageIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv, bool isPreviousTimeLevel) const
Calculate the storage integral.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:131
NumEqVector fluxIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Calculates the flux integral over a sub control volume face.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:274
NumEqVector computeStorage(const Problem &problem, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const Variables &volVars, const bool isPreviousStorage) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:113
Helper class for evaluating FE-based local residuals.
Definition: felocalresidual.hh:31
static void addFluxAndSourceTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes)
Add flux and source residual contribution for non-CV local dofs.
Definition: felocalresidual.hh:109
static void addStorageTerms(ResidualVector &residual, const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &prevElemVars, const ElementVariables &curElemVars, const Scalar timeStepSize)
Add storage residual contribution for non-CV local dofs.
Definition: felocalresidual.hh:46
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
Definition: flux.hh:179
NumEqVector advectiveMomentumFlux(const Context &context) const
Returns the advective momentum flux.
Definition: flux.hh:197
NumEqVector diffusiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:231
NumEqVector pressureContribution(const Context &context) const
Definition: flux.hh:270
Context for computing fluxes.
Definition: flux.hh:39
The flux function class for the Navier-Stokes model using control-volume finite element schemes.
Definition: flux.hh:300
NumEqVector pressureFluxIntegrand(const Context &context, const IpData &ipData) const
Definition: flux.hh:374
NumEqVector diffusiveMomentumFluxIntegrand(const Context &context, const IpData &ipData) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:347
NumEqVector advectiveMomentumFluxIntegral(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf, const VelocityVector &integratedVelocity) const
Returns the advective momentum flux contribution for a given integrated velocity at the face.
Definition: flux.hh:319
Context for interpolating data on interpolation points.
Definition: flux.hh:99
const GlobalPosition & velocity() const
Definition: flux.hh:141
Defines all properties used in Dumux.
Classes representing interpolation point data for control-volume finite element schemes.
Calculates the element-wise residual for control-volume finite element schemes.
Helper classes to compute the integration elements.
Helper functions for assembling FE-based local residuals.
Shape functions and gradients at an interpolation point.
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:34
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
constexpr bool hasProblemSourceWithIpDataInterface()
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:44
decltype(std::declval< P >().source(std::declval< FVG >(), std::declval< EV >(), std::declval< IPD >())) SourceWithIpDataInterface
helper struct detecting if a problem has new source interface
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:41
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.
Quadrature rules over sub-control volumes and sub-control volume faces.