25#ifndef DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
26#define DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
42template <
class TypeTag,
int numEnergyEqFlu
id>
45template<
class TypeTag>
52 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
56 using Element =
typename GridView::template Codim<0>::Entity;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
61 using Indices =
typename ModelTraits::Indices;
63 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
64 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
65 static constexpr auto energyEq0Idx = Indices::energyEq0Idx;
66 static constexpr auto energyEqSolidIdx = Indices::energyEqSolidIdx;
68 static constexpr auto numPhases = ModelTraits::numFluidPhases();
69 static constexpr auto numComponents = ModelTraits::numFluidComponents();
74 const SubControlVolume& scv,
75 const VolumeVariables& volVars,
79 storage[energyEq0Idx] += volVars.porosity()
80 * volVars.density(phaseIdx)
81 * volVars.internalEnergy(phaseIdx)
82 * volVars.saturation(phaseIdx);
89 const SubControlVolume& scv,
90 const VolumeVariables& volVars)
93 for(
int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
95 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
96 * volVars.solidHeatCapacity()
97 * volVars.solidDensity()
98 * (1.0 - volVars.porosity());
105 FluxVariables& fluxVars,
108 auto upwindTerm = [phaseIdx](
const auto& volVars)
109 {
return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
112 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
115 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
116 const auto& elemVolVars = fluxVars.elemVolVars();
117 const auto& scvf = fluxVars.scvFace();
118 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
119 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
121 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
124 if (phaseIdx == compIdx)
128 if (diffusiveFluxes[compIdx] > 0)
129 enthalpy += insideVolVars.enthalpy(phaseIdx);
131 enthalpy += outsideVolVars.enthalpy(phaseIdx);
135 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
137 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
143 FluxVariables& fluxVars)
146 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
148 for(
int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
149 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
162 const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const SubControlVolume &scv)
168 const auto& volVars = elemVolVars[scv];
169 const Scalar characteristicLength = volVars.characteristicLength() ;
173 const Scalar as = volVars.fluidSolidInterfacialArea();
176 const Scalar TFluid = volVars.temperatureFluid(0);
177 const Scalar TSolid = volVars.temperatureSolid();
179 Scalar solidToFluidEnergyExchange ;
181 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
183 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
185 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
187 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
189 for(
int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
194 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
197 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
200 DUNE_THROW(Dune::NotImplemented,
207template<
class TypeTag>
215 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
220 using Element =
typename GridView::template Codim<0>::Entity;
222 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
225 using Indices =
typename ModelTraits::Indices;
227 enum { numPhases = ModelTraits::numFluidPhases() };
228 enum { numEnergyEqFluid = ModelTraits::numEnergyEqFluid() };
229 enum { numEnergyEqSolid = ModelTraits::numEnergyEqSolid() };
230 enum { energyEq0Idx = Indices::energyEq0Idx };
231 enum { energyEqSolidIdx = Indices::energyEqSolidIdx};
232 enum { conti0EqIdx = Indices::conti0EqIdx };
234 enum { numComponents = ModelTraits::numFluidComponents() };
235 enum { phase0Idx = FluidSystem::phase0Idx};
236 enum { phase1Idx = FluidSystem::phase1Idx};
237 enum { sPhaseIdx = numPhases};
239 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
245 const SubControlVolume& scv,
246 const VolumeVariables& volVars,
249 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
250 * volVars.density(phaseIdx)
251 * volVars.internalEnergy(phaseIdx)
252 * volVars.saturation(phaseIdx);
257 FluxVariables& fluxVars,
260 auto upwindTerm = [phaseIdx](
const auto& volVars)
261 {
return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
264 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
267 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
268 const auto& elemVolVars = fluxVars.elemVolVars();
269 const auto& scvf = fluxVars.scvFace();
270 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
271 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
273 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
276 if (phaseIdx == compIdx)
280 if (diffusiveFluxes[compIdx] > 0)
281 enthalpy += insideVolVars.enthalpy(phaseIdx);
283 enthalpy += outsideVolVars.enthalpy(phaseIdx);
284 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
290 FluxVariables& fluxVars)
292 for(
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
294 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
296 for(
int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
298 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
311 const Element& element,
312 const FVElementGeometry& fvGeometry,
313 const ElementVolumeVariables& elemVolVars,
314 const SubControlVolume &scv)
317 const auto &volVars = elemVolVars[scv];
319 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
320 const Scalar aws = volVars.interfacialArea(phase0Idx, sPhaseIdx);
321 const Scalar ans = volVars.interfacialArea(phase1Idx, sPhaseIdx);
323 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
324 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
325 const Scalar Ts = volVars.temperatureSolid();
327 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
328 const Scalar lambdaNonWetting = volVars.fluidThermalConductivity(phase1Idx);
329 const Scalar lambdaSolid = volVars.solidThermalConductivity();
331 const Scalar lambdaWN =
harmonicMean(lambdaWetting, lambdaNonWetting);
332 const Scalar lambdaWS =
harmonicMean(lambdaWetting, lambdaSolid);
333 const Scalar lambdaNS =
harmonicMean(lambdaNonWetting, lambdaSolid);
335 const Scalar characteristicLength = volVars.characteristicLength() ;
336 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
338 const Scalar nusseltWN =
harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
339 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
340 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
342 const Scalar wettingToNonWettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * awn * lambdaWN * nusseltWN ;
343 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * aws * lambdaWS * nusseltWS ;
344 const Scalar nonWettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * ans * lambdaNS * nusseltNS ;
346 for(
int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
351 source[energyEq0Idx + phaseIdx] += ( - wettingToNonWettingEnergyExchange - wettingToSolidEnergyExchange);
354 source[energyEq0Idx + phaseIdx] += (+ wettingToNonWettingEnergyExchange - nonWettingToSolidEnergyExchange);
357 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonWettingToSolidEnergyExchange);
360 DUNE_THROW(Dune::NotImplemented,
366 if (!isfinite(source[energyEq0Idx + phaseIdx]))
367 DUNE_THROW(
NumericalProblem,
"Calculated non-finite source, " <<
"Tw="<< Tw <<
" Tn="<< Tn<<
" Ts="<< Ts);
371 if (enableChemicalNonEquilibrium)
384 const auto& fluidState = volVars.fluidState();
386 for(
int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
392 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
394 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
395 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
396 * FluidSystem::molarMass(compIdx)
397 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
402 for(
int compIdx =0; compIdx<numComponents; ++compIdx)
404 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
405 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
406 * FluidSystem::molarMass(compIdx)
407 *FluidSystem::componentEnthalpy(fluidState, phase0Idx, compIdx));
413 DUNE_THROW(Dune::NotImplemented,
Some exceptions thrown in DuMux
Provides 3rd order polynomial splines.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:68
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:43
static void computeSourceEnergy(NumEqVector &source, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
heat transfer between the phases for nonequilibrium models
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:161
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:88
static void fluidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
The energy storage in the fluid phase with index phaseIdx.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:73
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:142
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:104
static void computeSourceEnergy(NumEqVector &source, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:310
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:256
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:289
static void fluidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
The energy storage in the fluid phase with index phaseIdx.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:244
Declares all properties used in Dumux.