version 3.11-dev
freeflow/navierstokes/mass/1p/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_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
13#define DUMUX_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
14
15#include <type_traits>
16
20#include <dumux/common/concepts/variables_.hh>
24
25namespace Dumux {
26
32template<class Problem>
34: public std::false_type
35{};
36
41template<class TypeTag>
44{
47 using GridVariablesCache = Concept::GridVariablesCache_t<GridVariables>;
48 using ElementVariables = typename GridVariablesCache::LocalView;
49 using Variables = Concept::Variables_t<GridVariables>;
50
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
57 using GridView = typename GridGeometry::GridView;
58 using Element = typename GridView::template Codim<0>::Entity;
62
63 using Extrusion = Extrusion_t<GridGeometry>;
64
65public:
67 using ParentType::ParentType;
68
72 NumEqVector computeStorage(const Problem& problem,
73 const SubControlVolume& scv,
74 const Variables& vars) const
75 {
76 NumEqVector storage(0.0);
77 storage[ModelTraits::Indices::conti0EqIdx] = vars.density();
78
79 // consider energy storage for non-isothermal models
80 if constexpr (ModelTraits::enableEnergyBalance())
81 storage[ModelTraits::Indices::energyEqIdx] = vars.density() * vars.internalEnergy();
82
83 return storage;
84 }
85
95 NumEqVector storageIntegral(const FVElementGeometry& fvGeometry,
96 const ElementVariables& elemVars,
97 const SubControlVolume& scv,
98 bool isPreviousTimeLevel) const
99 {
100 const auto& vars = elemVars[scv];
101 // We apply mass lumping here
102 NumEqVector storage(0.0);
103 storage[ModelTraits::Indices::conti0EqIdx] = vars.density();
104
105 // consider energy storage for non-isothermal models
106 if constexpr (ModelTraits::enableEnergyBalance())
107 storage[ModelTraits::Indices::energyEqIdx] = vars.density() * vars.internalEnergy();
108
109 storage *= Extrusion::volume(fvGeometry, scv) * vars.extrusionFactor();
110
111 return storage;
112 }
113
122 NumEqVector sourceIntegral(const FVElementGeometry& fvGeometry,
123 const ElementVariables& elemVars,
124 const SubControlVolume& scv) const
125 {
126 const auto& problem = this->asImp().problem();
127
128 NumEqVector source(0.0);
129 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scv))
130 source += qpData.weight() * problem.source(fvGeometry, elemVars, qpData.ipData());
131
132 // ToDo: point source data with ipData
133 // add contribution from possible point sources
134 if (!problem.pointSourceMap().empty())
135 source += Extrusion::volume(fvGeometry, scv) * problem.scvPointSources(fvGeometry.element(), fvGeometry, elemVars, scv);
136
137 return source * elemVars[scv].extrusionFactor();
138 }
139
150 template<class ElementFluxVariablesCache>
151 NumEqVector computeFlux(const Problem& problem,
152 const Element& element,
153 const FVElementGeometry& fvGeometry,
154 const ElementVariables& elemVars,
155 const SubControlVolumeFace& scvf,
156 const ElementFluxVariablesCache& elemFluxVarsCache) const
157 {
158 FluxVariables fluxVars;
159 fluxVars.init(problem, element, fvGeometry, elemVars, scvf, elemFluxVarsCache);
160 auto flux = fluxVars.flux(0);
161
162 // the auxiliary flux is enabled if the trait is specialized for the problem
163 // this can be used, for example, to implement flux stabilization terms
165 flux += problem.auxiliaryFlux(element, fvGeometry, elemVars, elemFluxVarsCache, scvf);
166
167 return flux;
168 }
169
179 template<class ElementFluxVariablesCache>
180 [[deprecated("This function is deprecated and will be removed after release 3.11. "
181 "Use fluxIntegral(fvGeometry, elemVars, scvf) instead.")]]
182 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
183 const ElementVariables& elemVars,
184 const ElementFluxVariablesCache& elemFluxVarsCache,
185 const SubControlVolumeFace& scvf) const
186 {
187 const auto& problem = this->asImp().problem();
188 NumEqVector flux(0.0);
189 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
190 {
191 const auto& faceIpData = qpData.ipData();
192 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
193 }
194
195 static const auto upwindWeight
196 = getParamFromGroup<Scalar>(this->problem().paramGroup(), "Flux.UpwindWeight", 1.0);
197
198 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
199 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
200
201 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
202
203 // the auxiliary flux is enabled if the trait is specialized for the problem
204 // this can be used, for example, to implement flux stabilization terms
206 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, elemFluxVarsCache, scvf);
207
208 return flux * insideVars.extrusionFactor();
209 }
210
219 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
220 const ElementVariables& elemVars,
221 const SubControlVolumeFace& scvf) const
222 {
223 const auto& problem = this->asImp().problem();
224 NumEqVector flux(0.0);
225 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
226 {
227 const auto& faceIpData = qpData.ipData();
228 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
229 }
230
231 static const auto upwindWeight
232 = getParamFromGroup<Scalar>(this->problem().paramGroup(), "Flux.UpwindWeight", 1.0);
233
234 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
235 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
236
237 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
238
239 // the auxiliary flux is enabled if the trait is specialized for the problem
240 // this can be used, for example, to implement flux stabilization terms
242 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, scvf);
243
244 return flux * insideVars.extrusionFactor();
245 }
246};
247
248} // end namespace Dumux
249
250#endif
Element-wise calculation of the Navier-Stokes residual for single-phase flow.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:44
NumEqVector storageIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv, bool isPreviousTimeLevel) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:95
NumEqVector fluxIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Calculates the flux integral over a face of a sub control volume.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:182
NumEqVector sourceIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolume &scv) const
Calculate the source integral.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:122
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const Variables &vars) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:72
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the mass flux over a face of a sub control volume.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:151
NumEqVector fluxIntegral(const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolumeFace &scvf) const
Calculates the flux integral over a face of a sub control volume.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:219
Defines all properties used in Dumux.
The default local operator than can be specialized for each discretization scheme.
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
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
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
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.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Quadrature rules over sub-control volumes and sub-control volume faces.
Traits class to be specialized for problems to add auxiliary fluxes.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:35