version 3.10-dev
scalarfluxvariables.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_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
13#define DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
14
15#include <dumux/common/math.hh>
21
22namespace Dumux {
23
28template<class Problem,
29 class ModelTraits,
30 class FluxTypes,
31 class ElementVolumeVariables,
32 class ElementFluxVariablesCache,
33 class UpwindScheme = UpwindScheme<typename ProblemTraits<Problem>::GridGeometry>>
35: public FluxVariablesBase<Problem,
36 typename ProblemTraits<Problem>::GridGeometry::LocalView,
37 ElementVolumeVariables,
38 ElementFluxVariablesCache>
39{
40 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
41 using GridGeometry = typename ProblemTraits<Problem>::GridGeometry;
42 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
43 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
44 using FVElementGeometry = typename GridGeometry::LocalView;
46
47public:
48
50 {
51 static Scalar flux(const Problem& problem,
52 const Element& element,
53 const FVElementGeometry& fvGeometry,
54 const ElementVolumeVariables& elemVolVars,
55 const SubControlVolumeFace& scvf,
56 const int phaseIdx,
57 const ElementFluxVariablesCache& elemFluxVarsCache)
58 {
59 const auto velocity = problem.faceVelocity(element, fvGeometry, scvf);
60 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(fvGeometry, scvf)*extrusionFactor_(elemVolVars, scvf);
61 return volumeFlux;
62 }
63 };
64
68 template<typename FunctionType>
69 Scalar getAdvectiveFlux(const FunctionType& upwindTerm) const
70 {
71 if constexpr (ModelTraits::enableAdvection())
72 {
73 const auto& scvf = this->scvFace();
74 const auto velocity = this->problem().faceVelocity(this->element(), this->fvGeometry(), scvf);
75 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(this->fvGeometry(), scvf)*extrusionFactor_(this->elemVolVars(), scvf);
76 return UpwindScheme::apply(*this, upwindTerm, volumeFlux, 0/*phaseIdx*/);
77 }
78 else
79 return 0.0;
80 }
81
85 Scalar heatConductionFlux() const
86 {
87 if constexpr (ModelTraits::enableEnergyBalance())
88 {
89 using HeatConductionType = typename FluxTypes::HeatConductionType;
90 return HeatConductionType::flux(this->problem(),
91 this->element(),
92 this->fvGeometry(),
93 this->elemVolVars(),
94 this->scvFace(),
95 this->elemFluxVarsCache());
96 }
97 else
98 return 0.0;
99 }
100
104 Scalar heatAdvectionFlux() const
105 {
106 if constexpr (ModelTraits::enableEnergyBalance())
107 {
108 const auto upwindTerm = [](const auto& volVars) { return volVars.density() * volVars.enthalpy(); };
109 return getAdvectiveFlux(upwindTerm);
110 }
111 else
112 return 0.0;
113 }
114
118 Scalar heatFlux() const
119 {
121 }
122
126 template<class NumEqVector>
127 void addHeatFlux(NumEqVector& flux) const
128 {
129 if constexpr (ModelTraits::enableEnergyBalance())
130 flux[ModelTraits::Indices::energyEqIdx] = heatFlux();
131 }
132
133private:
134
135 static Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf)
136 {
137 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
138 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
139 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
140 }
141};
142
143} // end namespace Dumux
144
145#endif
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:33
const ProblemTraits< Problem >::GridGeometry::LocalView & fvGeometry() const
Definition: fluxvariablesbase.hh:66
The flux variables base class for scalar quantities balanced in the Navier-Stokes model.
Definition: scalarfluxvariables.hh:39
Scalar heatConductionFlux() const
Returns the conductive energy flux computed by the respective law.
Definition: scalarfluxvariables.hh:85
Scalar heatFlux() const
Returns the total energy flux.
Definition: scalarfluxvariables.hh:118
void addHeatFlux(NumEqVector &flux) const
Adds the energy flux to a given flux vector.
Definition: scalarfluxvariables.hh:127
Scalar heatAdvectionFlux() const
Returns the advective energy flux.
Definition: scalarfluxvariables.hh:104
Scalar getAdvectiveFlux(const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: scalarfluxvariables.hh:69
Type traits for problem classes.
Helper classes to compute the integration elements.
Base class for the upwind scheme.
Base class for the flux variables living on a sub control volume face.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:57
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
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:30
Define some often used mathematical functions.
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
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: scalarfluxvariables.hh:51
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:33