version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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-FileCopyrightText: 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>
24
25namespace Dumux {
26
32// forward declaration
33template <class TypeTag, int numEnergyEqFluid>
35
36template<class TypeTag>
37class EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
38{
43 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
44 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
48 using Element = typename GridView::template Codim<0>::Entity;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
50 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
51
53 using Indices = typename ModelTraits::Indices;
54
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;
59
60 static constexpr auto numPhases = ModelTraits::numFluidPhases();
61 static constexpr auto numComponents = ModelTraits::numFluidComponents();
62
63public:
64 template <typename T = void>
65 static void fluidPhaseStorage(NumEqVector& storage,
66 const SubControlVolume& scv,
67 const VolumeVariables& volVars,
68 int phaseIdx)
69 {
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.");
71 }
72
74 static void fluidPhaseStorage(NumEqVector& storage,
75 const Problem&,
76 const SubControlVolume& scv,
77 const VolumeVariables& volVars,
78 int phaseIdx)
79 {
80 // in case we have one energy equation for more than one fluid phase,
81 // add up parts on the one energy equation
82 storage[energyEq0Idx] += volVars.porosity()
83 * volVars.density(phaseIdx)
84 * volVars.internalEnergy(phaseIdx)
85 * volVars.saturation(phaseIdx);
86
87 }
88
89
91 static void solidPhaseStorage(NumEqVector& storage,
92 const SubControlVolume& scv,
93 const VolumeVariables& volVars)
94 {
95 // heat conduction for the fluid phases
96 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
97 {
98 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
99 * volVars.solidHeatCapacity()
100 * volVars.solidDensity()
101 * (1.0 - volVars.porosity());
102 }
103 }
104
111 static void heatDispersionFlux(NumEqVector& flux,
112 FluxVariables& fluxVars)
113 {}
114
116 static void heatConvectionFlux(NumEqVector& flux,
117 FluxVariables& fluxVars,
118 int phaseIdx)
119 {
120 auto upwindTerm = [phaseIdx](const auto& volVars)
121 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
122
123 //in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
124 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
125
126 //now add the diffusive part
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()];
132
133 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
134 {
135 //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
136 if (phaseIdx == compIdx)
137 continue;
138 //we need the upwind enthalpy. Even better would be the componentEnthalpy
139 auto enthalpy = 0.0;
140 if (diffusiveFluxes[compIdx] > 0)
141 enthalpy += insideVolVars.enthalpy(phaseIdx);
142 else
143 enthalpy += outsideVolVars.enthalpy(phaseIdx);
144
145 //check for the reference system and adapt units of the diffusive flux accordingly.
146 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
147 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
148 else
149 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
150 }
151 }
152
154 static void heatConductionFlux(NumEqVector& flux,
155 FluxVariables& fluxVars)
156 {
157 //in case we have one energy equation for more than one fluid phase we use an effective law in the nonequilibrium fourierslaw
158 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
159 //heat conduction for the solid phases
160 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
161 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
162 }
163
173 static void computeSourceEnergy(NumEqVector& source,
174 const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const SubControlVolume &scv)
178 {
179 // specialization for 2 fluid phases
180 const auto& volVars = elemVolVars[scv];
181 const Scalar characteristicLength = volVars.characteristicLength() ;
182
183 // interfacial area
184 // Shi & Wang, Transport in porous media (2011)
185 const Scalar as = volVars.fluidSolidInterfacialArea();
186
187 // temperature fluid is the same for both fluids
188 const Scalar TFluid = volVars.temperatureFluid(0);
189 const Scalar TSolid = volVars.temperatureSolid();
190
191 Scalar solidToFluidEnergyExchange ;
192
193 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
194
195 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
196
197 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
198
199 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
200
201 for(int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
202 {
203 switch (energyEqIdx)
204 {
205 case 0 :
206 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
207 break;
208 case 1 :
209 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
210 break;
211 default:
212 DUNE_THROW(Dune::NotImplemented,
213 "wrong index");
214 } // end switch
215 } // end energyEqIdx
216 } // end source
217};
218
223template<class TypeTag>
224class EnergyLocalResidualNonEquilibrium<TypeTag, 2/*numEnergyEqFluid*/>
225: public EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
226{
231 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
232 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
237 using Element = typename GridView::template Codim<0>::Entity;
238 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
239 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
240
242 using Indices = typename ModelTraits::Indices;
243
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;
250
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;
255
256 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
257
258public:
259 template <typename T = void>
260 static void fluidPhaseStorage(NumEqVector& storage,
261 const SubControlVolume& scv,
262 const VolumeVariables& volVars,
263 int phaseIdx)
264 {
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.");
266 }
267
269 static void fluidPhaseStorage(NumEqVector& storage,
270 const Problem&,
271 const SubControlVolume& scv,
272 const VolumeVariables& volVars,
273 int phaseIdx)
274 {
275 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
276 * volVars.density(phaseIdx)
277 * volVars.internalEnergy(phaseIdx)
278 * volVars.saturation(phaseIdx);
279
280 }
281
283 static void heatConvectionFlux(NumEqVector& flux,
284 FluxVariables& fluxVars,
285 int phaseIdx)
286 {
287 auto upwindTerm = [phaseIdx](const auto& volVars)
288 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
289
290 // in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
291 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
292
293 // add the diffusiv part
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()];
299
300 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
301 {
302 // 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
303 if (phaseIdx == compIdx)
304 continue;
305 // we need the upwind enthalpy. Even better would be the componentEnthalpy
306 auto enthalpy = 0.0;
307 if (diffusiveFluxes[compIdx] > 0)
308 enthalpy += insideVolVars.enthalpy(phaseIdx);
309 else
310 enthalpy += outsideVolVars.enthalpy(phaseIdx);
311 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
312 }
313 }
314
316 static void heatConductionFlux(NumEqVector& flux,
317 FluxVariables& fluxVars)
318 {
319 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
320 {
321 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
322 }
323 for(int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
324 {
325 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
326 }
327 }
328
335 static void heatDispersionFlux(NumEqVector& flux,
336 FluxVariables& fluxVars)
337 {}
338
348 static void computeSourceEnergy(NumEqVector& source,
349 const Element& element,
350 const FVElementGeometry& fvGeometry,
351 const ElementVolumeVariables& elemVolVars,
352 const SubControlVolume &scv)
353 {
354 // specialization for 2 fluid phases
355 const auto &volVars = elemVolVars[scv];
356
357 const Scalar areaWN = volVars.interfacialArea(phase0Idx, phase1Idx);
358 const Scalar areaWS = volVars.interfacialArea(phase0Idx, sPhaseIdx);
359 const Scalar areaNS = volVars.interfacialArea(phase1Idx, sPhaseIdx);
360
361 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
362 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
363 const Scalar Ts = volVars.temperatureSolid();
364
365 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
366 const Scalar lambdaNonwetting = volVars.fluidThermalConductivity(phase1Idx);
367 const Scalar lambdaSolid = volVars.solidThermalConductivity();
368
369 const Scalar lambdaWN = harmonicMean(lambdaWetting, lambdaNonwetting);
370 const Scalar lambdaWS = harmonicMean(lambdaWetting, lambdaSolid);
371 const Scalar lambdaNS = harmonicMean(lambdaNonwetting, lambdaSolid);
372
373 const Scalar characteristicLength = volVars.characteristicLength() ;
374 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
375
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);
379
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 ;
383
384 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
385 {
386 switch (phaseIdx)
387 {
388 case phase0Idx:
389 source[energyEq0Idx + phaseIdx] += ( - wettingToNonwettingEnergyExchange - wettingToSolidEnergyExchange);
390 break;
391 case phase1Idx:
392 source[energyEq0Idx + phaseIdx] += (+ wettingToNonwettingEnergyExchange - nonwettingToSolidEnergyExchange);
393 break;
394 case sPhaseIdx:
395 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonwettingToSolidEnergyExchange);
396 break;
397 default:
398 DUNE_THROW(Dune::NotImplemented,
399 "wrong index");
400 } // end switch
401
402
403 using std::isfinite;
404 if (!isfinite(source[energyEq0Idx + phaseIdx]))
405 DUNE_THROW(NumericalProblem, "Calculated non-finite source, " << "Tw="<< Tw << " Tn="<< Tn<< " Ts="<< Ts);
406 }// end phases
407
408 // we only need to do this for when there is more than 1 fluid phase
409 if (enableChemicalNonEquilibrium)
410 {
411 // Here comes the catch: We are not doing energy conservation for the whole
412 // system, but rather for each individual phase.
413 // -> Therefore the energy fluxes over each phase boundary need be
414 // individually accounted for.
415 // -> Each particle crossing a phase boundary does carry some mass and
416 // thus energy!
417 // -> Therefore, this contribution needs to be added.
418 // -> the particle always brings the energy of the originating phase.
419 // -> Energy advectivly transported into a phase = the moles of a component that go into a phase
420 // * molMass * enthalpy of the component in the *originating* phase
421
422 const auto& fluidState = volVars.fluidState();
423
424 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
425 {
426 switch (phaseIdx)
427 {
428 case phase0Idx:
429 // sum up the transferred energy by the components into the wetting phase
430 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
431 {
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) );
436 }
437 break;
438 case phase1Idx:
439 // sum up the transferred energy by the components into the nonwetting phase
440 for(int compIdx =0; compIdx<numComponents; ++compIdx)
441 {
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));
446 }
447 break;
448 case sPhaseIdx:
449 break; // no sorption
450 default:
451 DUNE_THROW(Dune::NotImplemented,
452 "wrong index");
453 } // end switch
454 } // end phases
455 } // EnableChemicalNonEquilibrium
456 } // end source
457};
458} // end namespace Dumux
459
460#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: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.
Type traits.
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.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24