24#ifndef DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
40template<
class Problem,
43 class ElementVolumeVariables,
44 class ElementFluxVariablesCache,
45 class UpwindScheme = UpwindScheme<typename ProblemTraits<Problem>::GridGeometry>>
48 typename ProblemTraits<Problem>::GridGeometry::LocalView,
49 ElementVolumeVariables,
50 ElementFluxVariablesCache>
52 using Scalar =
typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
54 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
55 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
67 const SubControlVolumeFace& scvf,
72 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(
fvGeometry, scvf)*extrusionFactor_(
elemVolVars, scvf);
80 template<
typename FunctionType>
83 if constexpr (ModelTraits::enableAdvection())
85 const auto& scvf = this->
scvFace();
87 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(this->
fvGeometry(), scvf)*extrusionFactor_(this->
elemVolVars(), scvf);
88 return UpwindScheme::apply(*
this, upwindTerm, volumeFlux, 0);
99 if constexpr (ModelTraits::enableEnergyBalance())
101 using HeatConductionType =
typename FluxTypes::HeatConductionType;
102 return HeatConductionType::flux(this->
problem(),
118 if constexpr (ModelTraits::enableEnergyBalance())
120 const auto upwindTerm = [](
const auto& volVars) {
return volVars.density() * volVars.enthalpy(); };
138 template<
class NumEqVector>
141 if constexpr (ModelTraits::enableEnergyBalance())
142 flux[ModelTraits::Indices::energyEqIdx] =
heatFlux();
147 static Scalar extrusionFactor_(
const ElementVolumeVariables&
elemVolVars,
const SubControlVolumeFace& scvf)
149 const auto& insideVolVars =
elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars =
elemVolVars[scvf.outsideScvIdx()];
151 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
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:42
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:69
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:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:45
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:45
const ElementVolumeVariables & elemVolVars() const
Definition: fluxvariablesbase.hh:81
const SubControlVolumeFace & scvFace() const
Definition: fluxvariablesbase.hh:75
const ElementFluxVariablesCache & elemFluxVarsCache() const
Definition: fluxvariablesbase.hh:84
const ProblemTraits< Problem >::GridGeometry::LocalView & fvGeometry() const
Definition: fluxvariablesbase.hh:78
const Element & element() const
Definition: fluxvariablesbase.hh:72
const Problem & problem() const
Definition: fluxvariablesbase.hh:69
The flux variables base class for scalar quantities balanced in the Navier-Stokes model.
Definition: scalarfluxvariables.hh:51
Scalar heatConductionFlux() const
Returns the conductive energy flux computed by the respective law.
Definition: scalarfluxvariables.hh:97
Scalar heatFlux() const
Returns the total energy flux.
Definition: scalarfluxvariables.hh:130
void addHeatFlux(NumEqVector &flux) const
Adds the energy flux to a given flux vector.
Definition: scalarfluxvariables.hh:139
Scalar heatAdvectionFlux() const
Returns the advective energy flux.
Definition: scalarfluxvariables.hh:116
Scalar getAdvectiveFlux(const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: scalarfluxvariables.hh:81
Definition: scalarfluxvariables.hh:62
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:63
Type traits for problem classes.
Base class for the upwind scheme.