24#ifndef DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH
25#define DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH
33template<
class Gr
idGeometry,
class FluxVariables,
class DiscretizationMethod,
bool enableEneryBalance,
bool isCompositional>
41template<
class Gr
idGeometry,
class FluxVariables,
bool enableEneryBalance,
bool isCompositional>
45 typename GridGeometry::DiscretizationMethod,
46 enableEneryBalance, isCompositional>;
52template<
class Gr
idGeometry,
class FluxVariables,
class DiscretizationMethod,
bool isCompositional>
58 template <
typename... Args>
63 template <
typename... Args>
72template<
class Gr
idGeometry,
class FluxVariables>
75 DiscretizationMethods::Staggered,
78 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
79 using FVElementGeometry =
typename GridGeometry::LocalView;
80 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
85 template<
class NumEqVector,
class VolumeVariables>
87 const VolumeVariables& volVars)
89 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
90 storage[localEnergyBalanceIdx] += volVars.density() * volVars.internalEnergy();
94 template<
class NumEqVector,
class Problem,
class ElementVolumeVariables,
class ElementFaceVariables>
96 const Problem& problem,
97 const Element &element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const ElementFaceVariables& elemFaceVars,
101 const SubControlVolumeFace& scvf)
103 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
105 auto upwindTerm = [](
const auto& volVars) {
return volVars.density() * volVars.enthalpy(); };
106 flux[localEnergyBalanceIdx] += FluxVariables::advectiveFluxForCellCenter(problem,
113 flux[localEnergyBalanceIdx] += FluxVariables::HeatConductionType::flux(problem,
125template<
class Gr
idGeometry,
class FluxVariables>
128 DiscretizationMethods::Staggered,
132 DiscretizationMethods::Staggered,
139 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
140 using FVElementGeometry =
typename GridGeometry::LocalView;
141 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
145 template<
class NumEqVector,
class Problem,
class ElementVolumeVariables,
class ElementFaceVariables>
147 const Problem& problem,
148 const Element &element,
149 const FVElementGeometry& fvGeometry,
150 const ElementVolumeVariables& elemVolVars,
151 const ElementFaceVariables& elemFaceVars,
152 const SubControlVolumeFace& scvf)
154 ParentType::heatFlux(flux, problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf);
156 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
157 auto diffusiveFlux = FluxVariables::MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf);
158 for (
int compIdx = 0; compIdx < FluxVariables::numComponents; ++compIdx)
160 const bool insideIsUpstream = scvf.directionSign() ==
sign(diffusiveFlux[compIdx]);
161 const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()];
164 flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx);
166 flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx)* elemVolVars[scvf.insideScvIdx()].molarMass(compIdx);
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
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
Definition: method.hh:103
Definition: freeflow/nonisothermal/localresidual.hh:34
static void heatFlux(Args &&... args)
do nothing for the isothermal case
Definition: freeflow/nonisothermal/localresidual.hh:64
static void fluidPhaseStorage(Args &&... args)
do nothing for the isothermal case
Definition: freeflow/nonisothermal/localresidual.hh:59
Specialization for staggered one-phase, non-isothermal models.
Definition: freeflow/nonisothermal/localresidual.hh:77
static void heatFlux(NumEqVector &flux, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf)
The convective and conductive heat fluxes in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:95
static void fluidPhaseStorage(NumEqVector &storage, const VolumeVariables &volVars)
The energy storage in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:86
static void heatFlux(NumEqVector &flux, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf)
The convective and conductive heat fluxes in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:146