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/concepts/variables_.hh>
21#include <dumux/common/typetraits/localdofs_.hh>
23
30
33
34namespace Dumux {
35
36namespace Detail {
37
39template<class P, class FVG, class EV, class IPD>
41 std::declval<P>().source(std::declval<FVG>(), std::declval<EV>(), std::declval<IPD>())
42);
43
44template<class P, class FVG, class EV, class IPD>
46{ return Dune::Std::is_detected<SourceWithIpDataInterface, P, FVG, EV, IPD>::value; }
47
48} // end namespace Detail
49
54template<class TypeTag>
57{
59
61
62 using GridVariablesCache = Concept::GridVariablesCache_t<GridVariables>;
63 using ElementVariables = typename GridVariablesCache::LocalView;
64 using Variables = Concept::Variables_t<GridVariables>;
65
69 using FVElementGeometry = typename GridGeometry::LocalView;
70 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
71 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
72 using GridView = typename GridGeometry::GridView;
73 using Element = typename GridView::template Codim<0>::Entity;
77 using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
78
79 using Extrusion = Extrusion_t<GridGeometry>;
80
82
83 static constexpr auto dim = GridView::dimension;
84
85 using LocalBasis = typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
86 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
90
92
93public:
95 using ElementResidualVector = typename ParentType::ElementResidualVector;
96 using ParentType::ParentType;
97
108 NumEqVector computeStorage(const Problem& problem,
109 const FVElementGeometry& fvGeometry,
110 const SubControlVolume& scv,
111 const Variables& vars,
112 const bool isPreviousStorage) const
113 {
114 return problem.density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scv), isPreviousStorage) * vars.velocity();
115 }
116
126 NumEqVector storageIntegral(const FVElementGeometry& fvGeometry,
127 const ElementVariables& elemVars,
128 const SubControlVolume& scv,
129 bool isPreviousTimeLevel) const
130 {
131 const auto& vars = elemVars[scv];
132 // We apply mass lumping
133 NumEqVector storage = this->asImp().problem().density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scv), isPreviousTimeLevel)
134 * vars.velocity();
135
136 storage *= Extrusion::volume(fvGeometry, scv) * vars.extrusionFactor();
137
138 return storage;
139 }
140
152 NumEqVector computeSource(const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVariables& elemVars,
156 const SubControlVolume& scv) const
157 {
158 NumEqVector source;
159
160 if constexpr (Detail::hasProblemSourceWithIpDataInterface<Problem, FVElementGeometry, ElementVariables, BaseIpData>())
161 {
162 source = problem.source(fvGeometry, elemVars, ipData(fvGeometry, scv.center()));
163
164 // ToDo: point source data with ipData
165 // add contribution from possible point sources
166 if (!problem.pointSourceMap().empty())
167 source += problem.scvPointSources(element, fvGeometry, elemVars, scv);
168 }
169 else
170 source = ParentType::computeSource(problem, element, fvGeometry, elemVars, scv);
171
172
173 // add rho*g (note that gravity might be zero in case it's disabled in the problem)
174 const auto& data = ipData(fvGeometry, scv);
175 source += problem.density(element, fvGeometry, data) * problem.gravity();
176
177 // Axisymmetric problems in 2D feature an extra source term arising from the transformation to cylindrical coordinates.
178 // See Ferziger/Peric: Computational methods for Fluid Dynamics (2020)
179 // https://doi.org/10.1007/978-3-319-99693-6
180 // Chapter 9.9 and Eq. (9.81) and comment on finite volume methods
181 if constexpr (dim == 2 && isRotationalExtrusion<Extrusion>)
182 {
183 // the radius with respect to the rotation axis
184 const auto r = scv.center()[Extrusion::radialAxis] - fvGeometry.gridGeometry().bBoxMin()[Extrusion::radialAxis];
185
186 // The velocity term is new with respect to Cartesian coordinates and handled below as a source term
187 // It only enters the balance of the momentum balance in radial direction
188 source[Extrusion::radialAxis] += -2.0*problem.effectiveViscosity(element, fvGeometry, data)
189 * elemVars[scv].velocity(Extrusion::radialAxis) / (r*r);
190
191 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
192 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
193 // is new with respect to Cartesian coordinates and handled below as a source term.
194 source[Extrusion::radialAxis] += problem.pressure(element, fvGeometry, data)/r;
195 }
196
197 return source;
198 }
199
208 NumEqVector sourceIntegral(const FVElementGeometry& fvGeometry,
209 const ElementVariables& elemVars,
210 const SubControlVolume& scv) const
211 {
212 static_assert(!(dim == 2 && isRotationalExtrusion<Extrusion>), "Rotational extrusion source terms are not implemented for integral interface.");
213
214 const auto& problem = this->asImp().problem();
215
216 NumEqVector source(0.0);
217 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scv))
218 {
219 source += qpData.weight() * (problem.source(fvGeometry, elemVars, qpData.ipData())
220 + problem.density(fvGeometry.element(), fvGeometry, qpData.ipData()) * problem.gravity());
221 }
222
223 // ToDo: point source data with ipData
224 // add contribution from possible point sources
225 if (!problem.pointSourceMap().empty())
226 source += Extrusion::volume(fvGeometry, scv) * problem.scvPointSources(fvGeometry.element(), fvGeometry, elemVars, scv);
227
228 source *= elemVars[scv].extrusionFactor();
229
230 return source;
231 }
232
243 template<class ElementFluxVariablesCache>
244 NumEqVector computeFlux(const Problem& problem,
245 const Element& element,
246 const FVElementGeometry& fvGeometry,
247 const ElementVariables& elemVars,
248 const SubControlVolumeFace& scvf,
249 const ElementFluxVariablesCache& elemFluxVarsCache) const
250 {
252 FluxContext context(problem, fvGeometry, elemVars, elemFluxVarsCache, scvf);
253 FluxHelper fluxHelper;
254
255 NumEqVector flux(0.0);
256 flux += fluxHelper.advectiveMomentumFlux(context);
257 flux += fluxHelper.diffusiveMomentumFlux(context);
258 flux += fluxHelper.pressureContribution(context);
259 return flux;
260 }
261
271 template<class ElementFluxVariablesCache>
272 [[deprecated("This function is deprecated and will be removed after release 3.11. "
273 "Use fluxIntegral without elemFluxVarsCache instead.")]]
274 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
275 const ElementVariables& elemVars,
276 const ElementFluxVariablesCache& elemFluxVarsCache,
277 const SubControlVolumeFace& scvf) const
278 {
279 using FluxVariablesCache = typename ElementFluxVariablesCache::FluxVariablesCache;
281
282 const auto& problem = this->asImp().problem();
283
284 NumEqVector flux(0.0);
285 GlobalPosition velIntegral(0.0);
286 FluxFunctionHelper fluxFunctionHelper;
287
288 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
289 {
290 const auto& fluxVarsCache = elemFluxVarsCache[qpData.ipData()];
291 FluxFunctionContext context(this->problem(), fvGeometry, elemVars, fluxVarsCache);
292
293 velIntegral += context.velocity() * qpData.weight();
294 flux += qpData.weight() * ( fluxFunctionHelper.diffusiveMomentumFluxIntegrand(context, qpData.ipData())
295 + fluxFunctionHelper.pressureFluxIntegrand(context, qpData.ipData()) );
296 }
297 flux += fluxFunctionHelper.advectiveMomentumFluxIntegral(problem, fvGeometry, elemVars, scvf, velIntegral);
298
299 flux *= elemVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
300
301 return flux;
302 }
303
312 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
313 const ElementVariables& elemVars,
314 const SubControlVolumeFace& scvf) const
315 {
316 const auto& problem = this->asImp().problem();
317
318 NumEqVector flux(0.0);
319 GlobalPosition velIntegral(0.0);
320 FluxFunctionHelper fluxFunctionHelper;
322
323 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
324 {
325 const auto& ipCache = cache(elemVars, qpData.ipData());
326 FluxFunctionContext context(this->problem(), fvGeometry, elemVars, ipCache);
327
328 velIntegral += context.velocity() * qpData.weight();
329 flux += qpData.weight() * ( fluxFunctionHelper.diffusiveMomentumFluxIntegrand(context, qpData.ipData())
330 + fluxFunctionHelper.pressureFluxIntegrand(context, qpData.ipData()) );
331 }
332 flux += fluxFunctionHelper.advectiveMomentumFluxIntegral(problem, fvGeometry, elemVars, scvf, velIntegral);
333
334 flux *= elemVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
335
336 return flux;
337 }
338
340 const Problem& problem,
341 const Element& element,
342 const FVElementGeometry& fvGeometry,
343 const ElementVariables& prevElemVolVars,
344 const ElementVariables& curElemVolVars) const
345 {
347 residual, problem, fvGeometry, prevElemVolVars, curElemVolVars, this->timeLoop().timeStepSize()
348 );
349 }
350
351 template<class ElementFluxVariablesCache>
352 [[deprecated("This function is deprecated and will be removed after release 3.11. "
353 "Use addToElementFluxAndSourceResidual without elemFluxVarsCache instead.")]]
355 const Problem& problem,
356 const Element& element,
357 const FVElementGeometry& fvGeometry,
358 const ElementVariables& elemVars,
359 const ElementFluxVariablesCache& elemFluxVarsCache,
360 const ElementBoundaryTypes &elemBcTypes) const
361 {
363 residual, problem, fvGeometry, elemVars, elemFluxVarsCache, elemBcTypes
364 );
365 }
366
368 const Problem& problem,
369 const Element& element,
370 const FVElementGeometry& fvGeometry,
371 const ElementVariables& elemVars,
372 const ElementBoundaryTypes &elemBcTypes) const
373 {
375 residual, problem, fvGeometry, elemVars, elemBcTypes
376 );
377 }
378
379};
380
381} // end namespace Dumux
382
383#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:31
Element-wise calculation of the Navier-Stokes residual for models using CVFE discretizations.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:57
NumEqVector sourceIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source integral.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:208
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:152
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:244
void addToElementFluxAndSourceResidual(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementBoundaryTypes &elemBcTypes) const
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:367
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:354
NumEqVector fluxIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Calculates the flux integral over a sub control volume face.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:312
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 &vars, const bool isPreviousStorage) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:108
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:339
typename ParentType::ElementResidualVector ElementResidualVector
Use the parent type's constructor.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:95
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:126
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:111
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:373
NumEqVector diffusiveMomentumFluxIntegrand(const Context &context, const IpData &ipData) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:346
NumEqVector advectiveMomentumFluxIntegral(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, 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 Problem & problem() const
Definition: flux.hh:129
Defines all properties used in Dumux.
Classes representing interpolation point data for control-volume finite element schemes.
The default local operator than can be specialized for each discretization scheme.
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:45
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:42
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition: defaultlocaloperator.hh:27
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.