12#ifndef DUMUX_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
13#define DUMUX_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
20#include <dumux/common/concepts/variables_.hh>
32template<
class Problem>
34:
public std::false_type
41template<
class TypeTag>
47 using GridVariablesCache = Concept::GridVariablesCache_t<GridVariables>;
48 using ElementVariables =
typename GridVariablesCache::LocalView;
49 using Variables = Concept::Variables_t<GridVariables>;
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;
67 using ParentType::ParentType;
73 const SubControlVolume& scv,
74 const Variables& vars)
const
76 NumEqVector storage(0.0);
77 storage[ModelTraits::Indices::conti0EqIdx] = vars.density();
80 if constexpr (ModelTraits::enableEnergyBalance())
81 storage[ModelTraits::Indices::energyEqIdx] = vars.density() * vars.internalEnergy();
96 const ElementVariables& elemVars,
97 const SubControlVolume& scv,
98 bool isPreviousTimeLevel)
const
100 const auto& vars = elemVars[scv];
102 NumEqVector storage(0.0);
103 storage[ModelTraits::Indices::conti0EqIdx] = vars.density();
106 if constexpr (ModelTraits::enableEnergyBalance())
107 storage[ModelTraits::Indices::energyEqIdx] = vars.density() * vars.internalEnergy();
123 const ElementVariables& elemVars,
124 const SubControlVolume& scv)
const
126 const auto& problem = this->asImp().problem();
128 NumEqVector source(0.0);
130 source += qpData.weight() * problem.source(fvGeometry, elemVars, qpData.ipData());
134 if (!problem.pointSourceMap().empty())
135 source +=
Extrusion::volume(fvGeometry, scv) * problem.scvPointSources(fvGeometry.element(), fvGeometry, elemVars, scv);
137 return source * elemVars[scv].extrusionFactor();
150 template<
class ElementFluxVariablesCache>
152 const Element& element,
153 const FVElementGeometry& fvGeometry,
154 const ElementVariables& elemVars,
155 const SubControlVolumeFace& scvf,
156 const ElementFluxVariablesCache& elemFluxVarsCache)
const
158 FluxVariables fluxVars;
159 fluxVars.init(problem, element, fvGeometry, elemVars, scvf, elemFluxVarsCache);
160 auto flux = fluxVars.flux(0);
165 flux += problem.auxiliaryFlux(element, fvGeometry, elemVars, elemFluxVarsCache, scvf);
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.")]]
183 const ElementVariables& elemVars,
184 const ElementFluxVariablesCache& elemFluxVarsCache,
185 const SubControlVolumeFace& scvf)
const
187 const auto& problem = this->asImp().problem();
188 NumEqVector flux(0.0);
191 const auto& faceIpData = qpData.ipData();
192 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
195 static const auto upwindWeight
196 = getParamFromGroup<Scalar>(this->problem().paramGroup(),
"Flux.UpwindWeight", 1.0);
198 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
199 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
201 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
206 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, elemFluxVarsCache, scvf);
208 return flux * insideVars.extrusionFactor();
220 const ElementVariables& elemVars,
221 const SubControlVolumeFace& scvf)
const
223 const auto& problem = this->asImp().problem();
224 NumEqVector flux(0.0);
227 const auto& faceIpData = qpData.ipData();
228 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
231 static const auto upwindWeight
232 = getParamFromGroup<Scalar>(this->problem().paramGroup(),
"Flux.UpwindWeight", 1.0);
234 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
235 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
237 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
242 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, scvf);
244 return flux * insideVars.extrusionFactor();
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
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