14#ifndef DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
15#define DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
33template <
class TypeTag,
int numEnergyEqFlu
id>
36template<
class TypeTag>
44 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
48 using Element =
typename GridView::template Codim<0>::Entity;
50 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
53 using Indices =
typename ModelTraits::Indices;
55 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
56 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
57 static constexpr auto energyEq0Idx = Indices::energyEq0Idx;
58 static constexpr auto energyEqSolidIdx = Indices::energyEqSolidIdx;
60 static constexpr auto numPhases = ModelTraits::numFluidPhases();
61 static constexpr auto numComponents = ModelTraits::numFluidComponents();
64 template <
typename T =
void>
66 const SubControlVolume& scv,
67 const VolumeVariables& volVars,
70 static_assert(
AlwaysFalse<T>::value,
"Deprecated interface that has been removed! Use new interface with additional argument problem instead. Will be entirely removed after release 3.10.");
76 const SubControlVolume& scv,
77 const VolumeVariables& volVars,
82 storage[energyEq0Idx] += volVars.porosity()
83 * volVars.density(phaseIdx)
84 * volVars.internalEnergy(phaseIdx)
85 * volVars.saturation(phaseIdx);
92 const SubControlVolume& scv,
93 const VolumeVariables& volVars)
96 for(
int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
98 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
99 * volVars.solidHeatCapacity()
100 * volVars.solidDensity()
101 * (1.0 - volVars.porosity());
112 FluxVariables& fluxVars)
117 FluxVariables& fluxVars,
120 auto upwindTerm = [phaseIdx](
const auto& volVars)
121 {
return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
124 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
127 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
128 const auto& elemVolVars = fluxVars.elemVolVars();
129 const auto& scvf = fluxVars.scvFace();
130 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
131 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
133 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
136 if (phaseIdx == compIdx)
140 if (diffusiveFluxes[compIdx] > 0)
141 enthalpy += insideVolVars.enthalpy(phaseIdx);
143 enthalpy += outsideVolVars.enthalpy(phaseIdx);
147 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
149 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
155 FluxVariables& fluxVars)
158 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
160 for(
int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
161 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
174 const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const SubControlVolume &scv)
180 const auto& volVars = elemVolVars[scv];
181 const Scalar characteristicLength = volVars.characteristicLength() ;
185 const Scalar as = volVars.fluidSolidInterfacialArea();
188 const Scalar TFluid = volVars.temperatureFluid(0);
189 const Scalar TSolid = volVars.temperatureSolid();
191 Scalar solidToFluidEnergyExchange ;
193 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
195 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
197 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
199 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
201 for(
int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
206 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
209 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
212 DUNE_THROW(Dune::NotImplemented,
223template<
class TypeTag>
232 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
237 using Element =
typename GridView::template Codim<0>::Entity;
239 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
242 using Indices =
typename ModelTraits::Indices;
244 static constexpr auto numPhases = ModelTraits::numFluidPhases();
245 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
246 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
247 static constexpr int energyEq0Idx = Indices::energyEq0Idx;
248 static constexpr int energyEqSolidIdx = Indices::energyEqSolidIdx;
249 static constexpr int conti0EqIdx = Indices::conti0EqIdx;
251 static constexpr auto numComponents = ModelTraits::numFluidComponents();
252 static constexpr int phase0Idx = FluidSystem::phase0Idx;
253 static constexpr int phase1Idx = FluidSystem::phase1Idx;
254 static constexpr int sPhaseIdx = numPhases;
256 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
259 template <
typename T =
void>
261 const SubControlVolume& scv,
262 const VolumeVariables& volVars,
265 static_assert(
AlwaysFalse<T>::value,
"Deprecated interface that has been removed! Use new interface with additional argument problem instead. Will be entirely removed after release 3.10.");
271 const SubControlVolume& scv,
272 const VolumeVariables& volVars,
275 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
276 * volVars.density(phaseIdx)
277 * volVars.internalEnergy(phaseIdx)
278 * volVars.saturation(phaseIdx);
284 FluxVariables& fluxVars,
287 auto upwindTerm = [phaseIdx](
const auto& volVars)
288 {
return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
291 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
294 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
295 const auto& elemVolVars = fluxVars.elemVolVars();
296 const auto& scvf = fluxVars.scvFace();
297 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
298 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
300 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
303 if (phaseIdx == compIdx)
307 if (diffusiveFluxes[compIdx] > 0)
308 enthalpy += insideVolVars.enthalpy(phaseIdx);
310 enthalpy += outsideVolVars.enthalpy(phaseIdx);
311 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
317 FluxVariables& fluxVars)
319 for(
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
321 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
323 for(
int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
325 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
336 FluxVariables& fluxVars)
349 const Element& element,
350 const FVElementGeometry& fvGeometry,
351 const ElementVolumeVariables& elemVolVars,
352 const SubControlVolume &scv)
355 const auto &volVars = elemVolVars[scv];
357 const Scalar areaWN = volVars.interfacialArea(phase0Idx, phase1Idx);
358 const Scalar areaWS = volVars.interfacialArea(phase0Idx, sPhaseIdx);
359 const Scalar areaNS = volVars.interfacialArea(phase1Idx, sPhaseIdx);
361 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
362 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
363 const Scalar Ts = volVars.temperatureSolid();
365 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
366 const Scalar lambdaNonwetting = volVars.fluidThermalConductivity(phase1Idx);
367 const Scalar lambdaSolid = volVars.solidThermalConductivity();
369 const Scalar lambdaWN =
harmonicMean(lambdaWetting, lambdaNonwetting);
370 const Scalar lambdaWS =
harmonicMean(lambdaWetting, lambdaSolid);
371 const Scalar lambdaNS =
harmonicMean(lambdaNonwetting, lambdaSolid);
373 const Scalar characteristicLength = volVars.characteristicLength() ;
374 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
376 const Scalar nusseltWN =
harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
377 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
378 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
380 const Scalar wettingToNonwettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * areaWN * lambdaWN * nusseltWN ;
381 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * areaWS * lambdaWS * nusseltWS ;
382 const Scalar nonwettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * areaNS * lambdaNS * nusseltNS ;
384 for(
int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
389 source[energyEq0Idx + phaseIdx] += ( - wettingToNonwettingEnergyExchange - wettingToSolidEnergyExchange);
392 source[energyEq0Idx + phaseIdx] += (+ wettingToNonwettingEnergyExchange - nonwettingToSolidEnergyExchange);
395 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonwettingToSolidEnergyExchange);
398 DUNE_THROW(Dune::NotImplemented,
404 if (!isfinite(source[energyEq0Idx + phaseIdx]))
405 DUNE_THROW(
NumericalProblem,
"Calculated non-finite source, " <<
"Tw="<< Tw <<
" Tn="<< Tn<<
" Ts="<< Ts);
409 if (enableChemicalNonEquilibrium)
422 const auto& fluidState = volVars.fluidState();
424 for(
int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
430 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
432 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
433 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
434 * FluidSystem::molarMass(compIdx)
435 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
440 for(
int compIdx =0; compIdx<numComponents; ++compIdx)
442 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
443 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
444 * FluidSystem::molarMass(compIdx)
445 *FluidSystem::componentEnthalpy(fluidState, phase0Idx, compIdx));
451 DUNE_THROW(Dune::NotImplemented,
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:173
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:91
static void fluidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:65
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:111
static void fluidPhaseStorage(NumEqVector &storage, const Problem &, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
The energy storage in the fluid phase with index phaseIdx.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:74
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:154
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:116
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:335
static void fluidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:260
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:348
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:283
static void fluidPhaseStorage(NumEqVector &storage, const Problem &, const SubControlVolume &scv, const VolumeVariables &volVars, int phaseIdx)
The energy storage in the fluid phase with index phaseIdx.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:269
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:316
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:34
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
Defines all properties used in Dumux.
Some exceptions thrown in DuMux
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
A helper to deduce a vector with the same size as numbers of equations.
Provides 3rd order polynomial splines.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24