version 3.8
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-FileCopyrightInfo: 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
19
23
25
26namespace Dumux {
27
32template<class TypeTag>
34: public CVFELocalResidual<TypeTag>
35{
37
39
40 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
41 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
42 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
43
44 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
45 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
46
51 using FVElementGeometry = typename GridGeometry::LocalView;
52 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using GridView = typename GridGeometry::GridView;
55 using Element = typename GridView::template Codim<0>::Entity;
59 using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
60
61 using Extrusion = Extrusion_t<GridGeometry>;
62
64
65 static constexpr auto dim = GridView::dimension;
66
69
70public:
72 using ParentType::ParentType;
73
84 NumEqVector computeStorage(const Problem& problem,
85 const FVElementGeometry& fvGeometry,
86 const SubControlVolume& scv,
87 const VolumeVariables& volVars,
88 const bool isPreviousStorage) const
89 {
90 return problem.density(fvGeometry.element(), fvGeometry, scv, isPreviousStorage) * volVars.velocity();
91 }
92
104 NumEqVector computeSource(const Problem& problem,
105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const SubControlVolume& scv) const
109 {
110 NumEqVector source = ParentType::computeSource(problem, element, fvGeometry, elemVolVars, scv);
111
112 // add rho*g (note that gravity might be zero in case it's disabled in the problem)
113 source += problem.density(element, fvGeometry, scv) * problem.gravity();
114
115 // Axisymmetric problems in 2D feature an extra source terms arising from the transformation to cylindrical coordinates.
116 // See Ferziger/Peric: Computational methods for Fluid Dynamics (2020)
117 // https://doi.org/10.1007/978-3-319-99693-6
118 // Chapter 9.9 and Eq. (9.81) and comment on finite volume methods
119 if constexpr (dim == 2 && isRotationalExtrusion<Extrusion>)
120 {
121 // the radius with respect to the rotation axis
122 const auto r = scv.center()[Extrusion::radialAxis] - fvGeometry.gridGeometry().bBoxMin()[Extrusion::radialAxis];
123
124 // The velocity term is new with respect to Cartesian coordinates and handled below as a source term
125 // it only enters the balance of the momentum balance in radial direction
126 source[Extrusion::radialAxis] += -2.0*problem.effectiveViscosity(element, fvGeometry, scv)
127 * elemVolVars[scv].velocity(Extrusion::radialAxis) / (r*r);
128
129 // Pressure term (needed because we incorporate pressure in terms of a surface integral).
130 // grad(p) becomes div(pI) + (p/r)*n_r in cylindrical coordinates. The second term
131 // is new with respect to Cartesian coordinates and handled below as a source term.
132 source[Extrusion::radialAxis] += problem.pressure(element, fvGeometry, scv)/r;
133 }
134
135 return source;
136 }
137
148 NumEqVector computeFlux(const Problem& problem,
149 const Element& element,
150 const FVElementGeometry& fvGeometry,
151 const ElementVolumeVariables& elemVolVars,
152 const SubControlVolumeFace& scvf,
153 const ElementFluxVariablesCache& elemFluxVarsCache) const
154 {
155 FluxContext context(problem, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
156 FluxHelper fluxHelper;
157
158 NumEqVector flux(0.0);
159 flux += fluxHelper.advectiveMomentumFlux(context);
160 flux += fluxHelper.diffusiveMomentumFlux(context);
161 flux += fluxHelper.pressureContribution(context);
162 return flux;
163 }
164};
165
166} // end namespace Dumux
167
168#endif
The element-wise residual for control-volume finite element schemes.
Definition: cvfelocalresidual.hh:60
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:35
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:208
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:472
Element-wise calculation of the Navier-Stokes residual for models using the CVFE discretizations.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:35
NumEqVector computeStorage(const Problem &problem, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const VolumeVariables &volVars, const bool isPreviousStorage) const
Calculate the source term of the equation.
Definition: freeflow/navierstokes/momentum/cvfe/localresidual.hh:84
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, 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:148
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/cvfe/localresidual.hh:104
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
Definition: flux.hh:91
NumEqVector advectiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:109
NumEqVector diffusiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:143
NumEqVector pressureContribution(const Context &context) const
Definition: flux.hh:182
Context for computing fluxes.
Definition: flux.hh:39
Defines all properties used in Dumux.
Calculates the element-wise residual for control-volume finite element schemes.
Helper classes to compute the integration elements.
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
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: 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.