version 3.11-dev
Loading...
Searching...
No Matches
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
124 template<class ElementFluxVariablesCache>
125 NumEqVector computeFlux(const Problem& problem,
126 const Element& element,
127 const FVElementGeometry& fvGeometry,
128 const ElementVariables& elemVars,
129 const SubControlVolumeFace& scvf,
130 const ElementFluxVariablesCache& elemFluxVarsCache) const
131 {
132 FluxVariables fluxVars;
133 fluxVars.init(problem, element, fvGeometry, elemVars, scvf, elemFluxVarsCache);
134 auto flux = fluxVars.flux(0);
135
136 // the auxiliary flux is enabled if the trait is specialized for the problem
137 // this can be used, for example, to implement flux stabilization terms
139 flux += problem.auxiliaryFlux(element, fvGeometry, elemVars, elemFluxVarsCache, scvf);
140
141 return flux;
142 }
143
152 NumEqVector fluxIntegral(const FVElementGeometry& fvGeometry,
153 const ElementVariables& elemVars,
154 const SubControlVolumeFace& scvf) const
155 {
156 const auto& problem = this->asImp().problem();
157 NumEqVector flux(0.0);
158 for (const auto& qpData : CVFE::quadratureRule(fvGeometry, scvf))
159 {
160 const auto& faceIpData = qpData.ipData();
161 flux += qpData.weight() * (problem.velocity(fvGeometry, faceIpData) * faceIpData.unitOuterNormal());
162 }
163
164 static const auto upwindWeight
165 = getParamFromGroup<Scalar>(this->problem().paramGroup(), "Flux.UpwindWeight", 1.0);
166
167 const auto& insideVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
168 const auto& outsideVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
169
170 flux *= (upwindWeight * insideVars.density() + (1.0 - upwindWeight) * outsideVars.density());
171
172 // the auxiliary flux is enabled if the trait is specialized for the problem
173 // this can be used, for example, to implement flux stabilization terms
175 flux += problem.auxiliaryFlux(fvGeometry.element(), fvGeometry, elemVars, scvf);
176
177 return flux * insideVars.extrusionFactor();
178 }
179};
180
181} // end namespace Dumux
182
183#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 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:125
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:152
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
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
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:159
Definition adapt.hh:17
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition defaultlocaloperator.hh:26
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
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