version 3.8
porousmediumflow/nonequilibrium/thermal/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
15#define DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
16
17#include <cmath>
23
24namespace Dumux {
25
31// forward declaration
32template <class TypeTag, int numEnergyEqFluid>
34
35template<class TypeTag>
36class EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
37{
41 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
42 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
46 using Element = typename GridView::template Codim<0>::Entity;
47 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
48 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
49
51 using Indices = typename ModelTraits::Indices;
52
53 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
54 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
55 static constexpr auto energyEq0Idx = Indices::energyEq0Idx;
56 static constexpr auto energyEqSolidIdx = Indices::energyEqSolidIdx;
57
58 static constexpr auto numPhases = ModelTraits::numFluidPhases();
59 static constexpr auto numComponents = ModelTraits::numFluidComponents();
60
61public:
63 static void fluidPhaseStorage(NumEqVector& storage,
64 const SubControlVolume& scv,
65 const VolumeVariables& volVars,
66 int phaseIdx)
67 {
68 //in case we have one energy equation for more than one fluid phase, add up parts on the one energy equation
69 storage[energyEq0Idx] += volVars.porosity()
70 * volVars.density(phaseIdx)
71 * volVars.internalEnergy(phaseIdx)
72 * volVars.saturation(phaseIdx);
73
74 }
75
76
78 static void solidPhaseStorage(NumEqVector& storage,
79 const SubControlVolume& scv,
80 const VolumeVariables& volVars)
81 {
82 // heat conduction for the fluid phases
83 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
84 {
85 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
86 * volVars.solidHeatCapacity()
87 * volVars.solidDensity()
88 * (1.0 - volVars.porosity());
89 }
90 }
91
98 static void heatDispersionFlux(NumEqVector& flux,
99 FluxVariables& fluxVars)
100 {}
101
103 static void heatConvectionFlux(NumEqVector& flux,
104 FluxVariables& fluxVars,
105 int phaseIdx)
106 {
107 auto upwindTerm = [phaseIdx](const auto& volVars)
108 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
109
110 //in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
111 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
112
113 //now add the diffusive part
114 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
115 const auto& elemVolVars = fluxVars.elemVolVars();
116 const auto& scvf = fluxVars.scvFace();
117 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
118 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
119
120 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
121 {
122 //no diffusion of the main component, this is a hack to use normal fick's law which computes both diffusions (main and component). We only add the part from the component here
123 if (phaseIdx == compIdx)
124 continue;
125 //we need the upwind enthalpy. Even better would be the componentEnthalpy
126 auto enthalpy = 0.0;
127 if (diffusiveFluxes[compIdx] > 0)
128 enthalpy += insideVolVars.enthalpy(phaseIdx);
129 else
130 enthalpy += outsideVolVars.enthalpy(phaseIdx);
131
132 //check for the reference system and adapt units of the diffusive flux accordingly.
133 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
134 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
135 else
136 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
137 }
138 }
139
141 static void heatConductionFlux(NumEqVector& flux,
142 FluxVariables& fluxVars)
143 {
144 //in case we have one energy equation for more than one fluid phase we use an effective law in the nonequilibrium fourierslaw
145 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
146 //heat conduction for the solid phases
147 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
148 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
149 }
150
160 static void computeSourceEnergy(NumEqVector& source,
161 const Element& element,
162 const FVElementGeometry& fvGeometry,
163 const ElementVolumeVariables& elemVolVars,
164 const SubControlVolume &scv)
165 {
166 // specialization for 2 fluid phases
167 const auto& volVars = elemVolVars[scv];
168 const Scalar characteristicLength = volVars.characteristicLength() ;
169
170 // interfacial area
171 // Shi & Wang, Transport in porous media (2011)
172 const Scalar as = volVars.fluidSolidInterfacialArea();
173
174 // temperature fluid is the same for both fluids
175 const Scalar TFluid = volVars.temperatureFluid(0);
176 const Scalar TSolid = volVars.temperatureSolid();
177
178 Scalar solidToFluidEnergyExchange ;
179
180 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
181
182 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
183
184 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
185
186 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
187
188 for(int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
189 {
190 switch (energyEqIdx)
191 {
192 case 0 :
193 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
194 break;
195 case 1 :
196 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
197 break;
198 default:
199 DUNE_THROW(Dune::NotImplemented,
200 "wrong index");
201 } // end switch
202 } // end energyEqIdx
203 } // end source
204};
205
210template<class TypeTag>
211class EnergyLocalResidualNonEquilibrium<TypeTag, 2/*numEnergyEqFluid*/>
212: public EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
213{
217 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
218 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
223 using Element = typename GridView::template Codim<0>::Entity;
224 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
225 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
226
228 using Indices = typename ModelTraits::Indices;
229
230 enum { numPhases = ModelTraits::numFluidPhases() };
231 enum { numEnergyEqFluid = ModelTraits::numEnergyEqFluid() };
232 enum { numEnergyEqSolid = ModelTraits::numEnergyEqSolid() };
233 enum { energyEq0Idx = Indices::energyEq0Idx };
234 enum { energyEqSolidIdx = Indices::energyEqSolidIdx};
235 enum { conti0EqIdx = Indices::conti0EqIdx };
236
237 enum { numComponents = ModelTraits::numFluidComponents() };
238 enum { phase0Idx = FluidSystem::phase0Idx};
239 enum { phase1Idx = FluidSystem::phase1Idx};
240 enum { sPhaseIdx = numPhases};
241
242 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
243
244public:
245
247 static void fluidPhaseStorage(NumEqVector& storage,
248 const SubControlVolume& scv,
249 const VolumeVariables& volVars,
250 int phaseIdx)
251 {
252 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
253 * volVars.density(phaseIdx)
254 * volVars.internalEnergy(phaseIdx)
255 * volVars.saturation(phaseIdx);
256 }
257
259 static void heatConvectionFlux(NumEqVector& flux,
260 FluxVariables& fluxVars,
261 int phaseIdx)
262 {
263 auto upwindTerm = [phaseIdx](const auto& volVars)
264 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
265
266 // in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
267 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
268
269 // add the diffusiv part
270 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
271 const auto& elemVolVars = fluxVars.elemVolVars();
272 const auto& scvf = fluxVars.scvFace();
273 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
274 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
275
276 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
277 {
278 // no diffusion of the main component, this is a hack to use normal fick's law which computes both diffusions (main and component). We only add the part from the component here
279 if (phaseIdx == compIdx)
280 continue;
281 // we need the upwind enthalpy. Even better would be the componentEnthalpy
282 auto enthalpy = 0.0;
283 if (diffusiveFluxes[compIdx] > 0)
284 enthalpy += insideVolVars.enthalpy(phaseIdx);
285 else
286 enthalpy += outsideVolVars.enthalpy(phaseIdx);
287 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
288 }
289 }
290
292 static void heatConductionFlux(NumEqVector& flux,
293 FluxVariables& fluxVars)
294 {
295 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
296 {
297 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
298 }
299 for(int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
300 {
301 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
302 }
303 }
304
311 static void heatDispersionFlux(NumEqVector& flux,
312 FluxVariables& fluxVars)
313 {}
314
324 static void computeSourceEnergy(NumEqVector& source,
325 const Element& element,
326 const FVElementGeometry& fvGeometry,
327 const ElementVolumeVariables& elemVolVars,
328 const SubControlVolume &scv)
329 {
330 // specialization for 2 fluid phases
331 const auto &volVars = elemVolVars[scv];
332
333 const Scalar areaWN = volVars.interfacialArea(phase0Idx, phase1Idx);
334 const Scalar areaWS = volVars.interfacialArea(phase0Idx, sPhaseIdx);
335 const Scalar areaNS = volVars.interfacialArea(phase1Idx, sPhaseIdx);
336
337 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
338 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
339 const Scalar Ts = volVars.temperatureSolid();
340
341 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
342 const Scalar lambdaNonwetting = volVars.fluidThermalConductivity(phase1Idx);
343 const Scalar lambdaSolid = volVars.solidThermalConductivity();
344
345 const Scalar lambdaWN = harmonicMean(lambdaWetting, lambdaNonwetting);
346 const Scalar lambdaWS = harmonicMean(lambdaWetting, lambdaSolid);
347 const Scalar lambdaNS = harmonicMean(lambdaNonwetting, lambdaSolid);
348
349 const Scalar characteristicLength = volVars.characteristicLength() ;
350 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
351
352 const Scalar nusseltWN = harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
353 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
354 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
355
356 const Scalar wettingToNonwettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * areaWN * lambdaWN * nusseltWN ;
357 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * areaWS * lambdaWS * nusseltWS ;
358 const Scalar nonwettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * areaNS * lambdaNS * nusseltNS ;
359
360 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
361 {
362 switch (phaseIdx)
363 {
364 case phase0Idx:
365 source[energyEq0Idx + phaseIdx] += ( - wettingToNonwettingEnergyExchange - wettingToSolidEnergyExchange);
366 break;
367 case phase1Idx:
368 source[energyEq0Idx + phaseIdx] += (+ wettingToNonwettingEnergyExchange - nonwettingToSolidEnergyExchange);
369 break;
370 case sPhaseIdx:
371 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonwettingToSolidEnergyExchange);
372 break;
373 default:
374 DUNE_THROW(Dune::NotImplemented,
375 "wrong index");
376 } // end switch
377
378
379 using std::isfinite;
380 if (!isfinite(source[energyEq0Idx + phaseIdx]))
381 DUNE_THROW(NumericalProblem, "Calculated non-finite source, " << "Tw="<< Tw << " Tn="<< Tn<< " Ts="<< Ts);
382 }// end phases
383
384 // we only need to do this for when there is more than 1 fluid phase
385 if (enableChemicalNonEquilibrium)
386 {
387 // Here comes the catch: We are not doing energy conservation for the whole
388 // system, but rather for each individual phase.
389 // -> Therefore the energy fluxes over each phase boundary need be
390 // individually accounted for.
391 // -> Each particle crossing a phase boundary does carry some mass and
392 // thus energy!
393 // -> Therefore, this contribution needs to be added.
394 // -> the particle always brings the energy of the originating phase.
395 // -> Energy advectivly transported into a phase = the moles of a component that go into a phase
396 // * molMass * enthalpy of the component in the *originating* phase
397
398 const auto& fluidState = volVars.fluidState();
399
400 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
401 {
402 switch (phaseIdx)
403 {
404 case phase0Idx:
405 // sum up the transferred energy by the components into the wetting phase
406 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
407 {
408 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
409 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
410 * FluidSystem::molarMass(compIdx)
411 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
412 }
413 break;
414 case phase1Idx:
415 // sum up the transferred energy by the components into the nonwetting phase
416 for(int compIdx =0; compIdx<numComponents; ++compIdx)
417 {
418 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
419 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
420 * FluidSystem::molarMass(compIdx)
421 *FluidSystem::componentEnthalpy(fluidState, phase0Idx, compIdx));
422 }
423 break;
424 case sPhaseIdx:
425 break; // no sorption
426 default:
427 DUNE_THROW(Dune::NotImplemented,
428 "wrong index");
429 } // end switch
430 } // end phases
431 } // EnableChemicalNonEquilibrium
432 } // end source
433};
434} // end namespace Dumux
435
436#endif
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:160
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:78
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:98
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:63
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:141
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:103
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:311
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:324
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:259
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:292
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:247
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:33
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
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Provides 3rd order polynomial splines.