12#ifndef DUMUX_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
13#define DUMUX_FREEFLOW_NAVIERSTOKES_MASS_1P_LOCAL_RESIDUAL_HH
31template<
class Problem>
33:
public std::false_type
40template<
class TypeTag>
46 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
47 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
48 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
50 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
51 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using GridView =
typename GridGeometry::GridView;
60 using Element =
typename GridView::template Codim<0>::Entity;
69 using ParentType::ParentType;
75 const SubControlVolume& scv,
76 const VolumeVariables& volVars)
const
78 NumEqVector storage(0.0);
79 storage[ModelTraits::Indices::conti0EqIdx] = volVars.density();
82 if constexpr (ModelTraits::enableEnergyBalance())
83 storage[ModelTraits::Indices::energyEqIdx] = volVars.density() * volVars.internalEnergy();
98 const ElementVolumeVariables& elemVars,
99 const SubControlVolume& scv,
100 bool isPreviousTimeLevel)
const
102 const auto& volVars = elemVars[scv];
104 NumEqVector storage(0.0);
105 storage[ModelTraits::Indices::conti0EqIdx] = volVars.density();
108 if constexpr (ModelTraits::enableEnergyBalance())
109 storage[ModelTraits::Indices::energyEqIdx] = volVars.density() * volVars.internalEnergy();
125 const ElementVolumeVariables& elemVars,
126 const SubControlVolume& scv)
const
128 const auto& problem = this->asImp().problem();
130 NumEqVector source(0.0);
132 source += qpData.weight() * problem.source(fvGeometry, elemVars, qpData.ipData());
136 if (!problem.pointSourceMap().empty())
137 source +=
Extrusion::volume(fvGeometry, scv) * problem.scvPointSources(fvGeometry.element(), fvGeometry, elemVars, scv);
139 return source * elemVars[scv].extrusionFactor();
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& elemVolVars,
156 const SubControlVolumeFace& scvf,
157 const ElementFluxVariablesCache& elemFluxVarsCache)
const
159 FluxVariables fluxVars;
160 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
161 auto flux = fluxVars.flux(0);
166 flux += problem.auxiliaryFlux(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
181 const ElementVolumeVariables& elemVars,
182 const ElementFluxVariablesCache& elemFluxVarsCache,
183 const SubControlVolumeFace& scvf)
const
185 const auto& problem = this->asImp().problem();
186 NumEqVector flux(0.0);
189 const auto& faceIpData = qpData.ipData();
190 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
193 static const auto upwindWeight
194 = getParamFromGroup<Scalar>(this->problem().paramGroup(),
"Flux.UpwindWeight", 1.0);
196 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
197 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
199 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
204 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, elemFluxVarsCache, scvf);
206 return flux * insideVars.extrusionFactor();
Element-wise calculation of the Navier-Stokes residual for single-phase flow.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:43
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, 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:152
NumEqVector storageIntegral(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVars, const SubControlVolume &scv, bool isPreviousTimeLevel) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:97
NumEqVector sourceIntegral(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVars, const SubControlVolume &scv) const
Calculate the source integral.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:124
NumEqVector fluxIntegral(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &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:180
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculate the storage term of the equation.
Definition: freeflow/navierstokes/mass/1p/localresidual.hh:74
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:34