version 3.7
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
206template<class TypeTag>
207class EnergyLocalResidualNonEquilibrium<TypeTag, 2/*numEnergyEqFluid*/>
208: public EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
209{
213 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
214 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
219 using Element = typename GridView::template Codim<0>::Entity;
220 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
221 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
222
224 using Indices = typename ModelTraits::Indices;
225
226 enum { numPhases = ModelTraits::numFluidPhases() };
227 enum { numEnergyEqFluid = ModelTraits::numEnergyEqFluid() };
228 enum { numEnergyEqSolid = ModelTraits::numEnergyEqSolid() };
229 enum { energyEq0Idx = Indices::energyEq0Idx };
230 enum { energyEqSolidIdx = Indices::energyEqSolidIdx};
231 enum { conti0EqIdx = Indices::conti0EqIdx };
232
233 enum { numComponents = ModelTraits::numFluidComponents() };
234 enum { phase0Idx = FluidSystem::phase0Idx};
235 enum { phase1Idx = FluidSystem::phase1Idx};
236 enum { sPhaseIdx = numPhases};
237
238 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
239
240public:
241
243 static void fluidPhaseStorage(NumEqVector& storage,
244 const SubControlVolume& scv,
245 const VolumeVariables& volVars,
246 int phaseIdx)
247 {
248 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
249 * volVars.density(phaseIdx)
250 * volVars.internalEnergy(phaseIdx)
251 * volVars.saturation(phaseIdx);
252 }
253
255 static void heatConvectionFlux(NumEqVector& flux,
256 FluxVariables& fluxVars,
257 int phaseIdx)
258 {
259 auto upwindTerm = [phaseIdx](const auto& volVars)
260 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
261
262 // in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
263 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
264
265 // add the diffusiv part
266 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
267 const auto& elemVolVars = fluxVars.elemVolVars();
268 const auto& scvf = fluxVars.scvFace();
269 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
270 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
271
272 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
273 {
274 // 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
275 if (phaseIdx == compIdx)
276 continue;
277 // we need the upwind enthalpy. Even better would be the componentEnthalpy
278 auto enthalpy = 0.0;
279 if (diffusiveFluxes[compIdx] > 0)
280 enthalpy += insideVolVars.enthalpy(phaseIdx);
281 else
282 enthalpy += outsideVolVars.enthalpy(phaseIdx);
283 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
284 }
285 }
286
288 static void heatConductionFlux(NumEqVector& flux,
289 FluxVariables& fluxVars)
290 {
291 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
292 {
293 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
294 }
295 for(int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
296 {
297 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
298 }
299 }
300
307 static void heatDispersionFlux(NumEqVector& flux,
308 FluxVariables& fluxVars)
309 {}
310
320 static void computeSourceEnergy(NumEqVector& source,
321 const Element& element,
322 const FVElementGeometry& fvGeometry,
323 const ElementVolumeVariables& elemVolVars,
324 const SubControlVolume &scv)
325 {
326 // specialization for 2 fluid phases
327 const auto &volVars = elemVolVars[scv];
328
329 const Scalar areaWN = volVars.interfacialArea(phase0Idx, phase1Idx);
330 const Scalar areaWS = volVars.interfacialArea(phase0Idx, sPhaseIdx);
331 const Scalar areaNS = volVars.interfacialArea(phase1Idx, sPhaseIdx);
332
333 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
334 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
335 const Scalar Ts = volVars.temperatureSolid();
336
337 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
338 const Scalar lambdaNonwetting = volVars.fluidThermalConductivity(phase1Idx);
339 const Scalar lambdaSolid = volVars.solidThermalConductivity();
340
341 const Scalar lambdaWN = harmonicMean(lambdaWetting, lambdaNonwetting);
342 const Scalar lambdaWS = harmonicMean(lambdaWetting, lambdaSolid);
343 const Scalar lambdaNS = harmonicMean(lambdaNonwetting, lambdaSolid);
344
345 const Scalar characteristicLength = volVars.characteristicLength() ;
346 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
347
348 const Scalar nusseltWN = harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
349 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
350 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
351
352 const Scalar wettingToNonwettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * areaWN * lambdaWN * nusseltWN ;
353 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * areaWS * lambdaWS * nusseltWS ;
354 const Scalar nonwettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * areaNS * lambdaNS * nusseltNS ;
355
356 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
357 {
358 switch (phaseIdx)
359 {
360 case phase0Idx:
361 source[energyEq0Idx + phaseIdx] += ( - wettingToNonwettingEnergyExchange - wettingToSolidEnergyExchange);
362 break;
363 case phase1Idx:
364 source[energyEq0Idx + phaseIdx] += (+ wettingToNonwettingEnergyExchange - nonwettingToSolidEnergyExchange);
365 break;
366 case sPhaseIdx:
367 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonwettingToSolidEnergyExchange);
368 break;
369 default:
370 DUNE_THROW(Dune::NotImplemented,
371 "wrong index");
372 } // end switch
373
374
375 using std::isfinite;
376 if (!isfinite(source[energyEq0Idx + phaseIdx]))
377 DUNE_THROW(NumericalProblem, "Calculated non-finite source, " << "Tw="<< Tw << " Tn="<< Tn<< " Ts="<< Ts);
378 }// end phases
379
380 // we only need to do this for when there is more than 1 fluid phase
381 if (enableChemicalNonEquilibrium)
382 {
383 // Here comes the catch: We are not doing energy conservation for the whole
384 // system, but rather for each individual phase.
385 // -> Therefore the energy fluxes over each phase boundary need be
386 // individually accounted for.
387 // -> Each particle crossing a phase boundary does carry some mass and
388 // thus energy!
389 // -> Therefore, this contribution needs to be added.
390 // -> the particle always brings the energy of the originating phase.
391 // -> Energy advectivly transported into a phase = the moles of a component that go into a phase
392 // * molMass * enthalpy of the component in the *originating* phase
393
394 const auto& fluidState = volVars.fluidState();
395
396 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
397 {
398 switch (phaseIdx)
399 {
400 case phase0Idx:
401 // sum up the transferred energy by the components into the wetting phase
402 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
403 {
404 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
405 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
406 * FluidSystem::molarMass(compIdx)
407 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
408 }
409 break;
410 case phase1Idx:
411 // sum up the transferred energy by the components into the nonwetting phase
412 for(int compIdx =0; compIdx<numComponents; ++compIdx)
413 {
414 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
415 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
416 * FluidSystem::molarMass(compIdx)
417 *FluidSystem::componentEnthalpy(fluidState, phase0Idx, compIdx));
418 }
419 break;
420 case sPhaseIdx:
421 break; // no sorption
422 default:
423 DUNE_THROW(Dune::NotImplemented,
424 "wrong index");
425 } // end switch
426 } // end phases
427 } // EnableChemicalNonEquilibrium
428 } // end source
429};
430} // end namespace Dumux
431
432#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:307
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:320
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:255
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:288
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:243
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:267
Definition: adapt.hh:17
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
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.