12#ifndef DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
13#define DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
28template<
class Problem,
31 class ElementVolumeVariables,
32 class ElementFluxVariablesCache,
33 class UpwindScheme = UpwindScheme<typename ProblemTraits<Problem>::GridGeometry>>
36 typename ProblemTraits<Problem>::GridGeometry::LocalView,
37 ElementVolumeVariables,
38 ElementFluxVariablesCache>
40 using Scalar =
typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
42 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
43 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
44 using FVElementGeometry =
typename GridGeometry::LocalView;
55 const SubControlVolumeFace& scvf,
60 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(
fvGeometry, scvf)*extrusionFactor_(
elemVolVars, scvf);
68 template<
typename FunctionType>
71 if constexpr (ModelTraits::enableAdvection())
73 const auto& scvf = this->
scvFace();
75 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(this->
fvGeometry(), scvf)*extrusionFactor_(this->
elemVolVars(), scvf);
76 return UpwindScheme::apply(*
this, upwindTerm, volumeFlux, 0);
87 if constexpr (ModelTraits::enableEnergyBalance())
89 using HeatConductionType =
typename FluxTypes::HeatConductionType;
90 return HeatConductionType::flux(this->
problem(),
106 if constexpr (ModelTraits::enableEnergyBalance())
108 const auto upwindTerm = [](
const auto& volVars) {
return volVars.density() * volVars.enthalpy(); };
126 template<
class NumEqVector>
129 if constexpr (ModelTraits::enableEnergyBalance())
130 flux[ModelTraits::Indices::energyEqIdx] =
heatFlux();
135 static Scalar extrusionFactor_(
const ElementVolumeVariables&
elemVolVars,
const SubControlVolumeFace& scvf)
137 const auto& insideVolVars =
elemVolVars[scvf.insideScvIdx()];
138 const auto& outsideVolVars =
elemVolVars[scvf.outsideScvIdx()];
139 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:33
const ElementVolumeVariables & elemVolVars() const
Definition: fluxvariablesbase.hh:69
const SubControlVolumeFace & scvFace() const
Definition: fluxvariablesbase.hh:63
const ElementFluxVariablesCache & elemFluxVarsCache() const
Definition: fluxvariablesbase.hh:72
const ProblemTraits< Problem >::GridGeometry::LocalView & fvGeometry() const
Definition: fluxvariablesbase.hh:66
const Element & element() const
Definition: fluxvariablesbase.hh:60
const Problem & problem() const
Definition: fluxvariablesbase.hh:57
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.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: scalarfluxvariables.hh:50
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