Loading [MathJax]/extensions/tex2jax.js
3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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/*****************************************************************************
3 * See the file COPYING for full copying permissions. *
4 * *
5 * This program is free software: you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation, either version 3 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
17 *****************************************************************************/
25#ifndef DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
26#define DUMUX_ENERGY_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
27
28#include <cmath>
34
35namespace Dumux {
36
42// forward declaration
43template <class TypeTag, int numEnergyEqFluid>
45
46template<class TypeTag>
47class EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
48{
52 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
53 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60
62 using Indices = typename ModelTraits::Indices;
63
64 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
65 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
66 static constexpr auto energyEq0Idx = Indices::energyEq0Idx;
67 static constexpr auto energyEqSolidIdx = Indices::energyEqSolidIdx;
68
69 static constexpr auto numPhases = ModelTraits::numFluidPhases();
70 static constexpr auto numComponents = ModelTraits::numFluidComponents();
71
72public:
74 static void fluidPhaseStorage(NumEqVector& storage,
75 const SubControlVolume& scv,
76 const VolumeVariables& volVars,
77 int phaseIdx)
78 {
79 //in case we have one energy equation for more than one fluid phase, add up parts on the one energy equation
80 storage[energyEq0Idx] += volVars.porosity()
81 * volVars.density(phaseIdx)
82 * volVars.internalEnergy(phaseIdx)
83 * volVars.saturation(phaseIdx);
84
85 }
86
87
89 static void solidPhaseStorage(NumEqVector& storage,
90 const SubControlVolume& scv,
91 const VolumeVariables& volVars)
92 {
93 // heat conduction for the fluid phases
94 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
95 {
96 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
97 * volVars.solidHeatCapacity()
98 * volVars.solidDensity()
99 * (1.0 - volVars.porosity());
100 }
101 }
102
109 static void heatDispersionFlux(NumEqVector& flux,
110 FluxVariables& fluxVars)
111 {}
112
114 static void heatConvectionFlux(NumEqVector& flux,
115 FluxVariables& fluxVars,
116 int phaseIdx)
117 {
118 auto upwindTerm = [phaseIdx](const auto& volVars)
119 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
120
121 //in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
122 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
123
124 //now add the diffusive part
125 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
126 const auto& elemVolVars = fluxVars.elemVolVars();
127 const auto& scvf = fluxVars.scvFace();
128 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
129 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
130
131 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
132 {
133 //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
134 if (phaseIdx == compIdx)
135 continue;
136 //we need the upwind enthalpy. Even better would be the componentEnthalpy
137 auto enthalpy = 0.0;
138 if (diffusiveFluxes[compIdx] > 0)
139 enthalpy += insideVolVars.enthalpy(phaseIdx);
140 else
141 enthalpy += outsideVolVars.enthalpy(phaseIdx);
142
143 //check for the reference system and adapt units of the diffusive flux accordingly.
144 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
145 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
146 else
147 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
148 }
149 }
150
152 static void heatConductionFlux(NumEqVector& flux,
153 FluxVariables& fluxVars)
154 {
155 //in case we have one energy equation for more than one fluid phase we use an effective law in the nonequilibrium fourierslaw
156 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
157 //heat conduction for the solid phases
158 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
159 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
160 }
161
171 static void computeSourceEnergy(NumEqVector& source,
172 const Element& element,
173 const FVElementGeometry& fvGeometry,
174 const ElementVolumeVariables& elemVolVars,
175 const SubControlVolume &scv)
176 {
177 // specialization for 2 fluid phases
178 const auto& volVars = elemVolVars[scv];
179 const Scalar characteristicLength = volVars.characteristicLength() ;
180
181 // interfacial area
182 // Shi & Wang, Transport in porous media (2011)
183 const Scalar as = volVars.fluidSolidInterfacialArea();
184
185 // temperature fluid is the same for both fluids
186 const Scalar TFluid = volVars.temperatureFluid(0);
187 const Scalar TSolid = volVars.temperatureSolid();
188
189 Scalar solidToFluidEnergyExchange ;
190
191 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
192
193 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
194
195 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
196
197 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
198
199 for(int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
200 {
201 switch (energyEqIdx)
202 {
203 case 0 :
204 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
205 break;
206 case 1 :
207 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
208 break;
209 default:
210 DUNE_THROW(Dune::NotImplemented,
211 "wrong index");
212 } // end switch
213 } // end energyEqIdx
214 } // end source
215};
216
217template<class TypeTag>
218class EnergyLocalResidualNonEquilibrium<TypeTag, 2/*numEnergyEqFluid*/>
219: public EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
220{
224 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
225 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
230 using Element = typename GridView::template Codim<0>::Entity;
231 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
232 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
233
235 using Indices = typename ModelTraits::Indices;
236
237 enum { numPhases = ModelTraits::numFluidPhases() };
238 enum { numEnergyEqFluid = ModelTraits::numEnergyEqFluid() };
239 enum { numEnergyEqSolid = ModelTraits::numEnergyEqSolid() };
240 enum { energyEq0Idx = Indices::energyEq0Idx };
241 enum { energyEqSolidIdx = Indices::energyEqSolidIdx};
242 enum { conti0EqIdx = Indices::conti0EqIdx };
243
244 enum { numComponents = ModelTraits::numFluidComponents() };
245 enum { phase0Idx = FluidSystem::phase0Idx};
246 enum { phase1Idx = FluidSystem::phase1Idx};
247 enum { sPhaseIdx = numPhases};
248
249 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
250
251public:
252
254 static void fluidPhaseStorage(NumEqVector& storage,
255 const SubControlVolume& scv,
256 const VolumeVariables& volVars,
257 int phaseIdx)
258 {
259 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
260 * volVars.density(phaseIdx)
261 * volVars.internalEnergy(phaseIdx)
262 * volVars.saturation(phaseIdx);
263 }
264
266 static void heatConvectionFlux(NumEqVector& flux,
267 FluxVariables& fluxVars,
268 int phaseIdx)
269 {
270 auto upwindTerm = [phaseIdx](const auto& volVars)
271 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
272
273 // in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
274 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
275
276 // add the diffusiv part
277 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
278 const auto& elemVolVars = fluxVars.elemVolVars();
279 const auto& scvf = fluxVars.scvFace();
280 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
281 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
282
283 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
284 {
285 // 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
286 if (phaseIdx == compIdx)
287 continue;
288 // we need the upwind enthalpy. Even better would be the componentEnthalpy
289 auto enthalpy = 0.0;
290 if (diffusiveFluxes[compIdx] > 0)
291 enthalpy += insideVolVars.enthalpy(phaseIdx);
292 else
293 enthalpy += outsideVolVars.enthalpy(phaseIdx);
294 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
295 }
296 }
297
299 static void heatConductionFlux(NumEqVector& flux,
300 FluxVariables& fluxVars)
301 {
302 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
303 {
304 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
305 }
306 for(int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
307 {
308 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
309 }
310 }
311
318 static void heatDispersionFlux(NumEqVector& flux,
319 FluxVariables& fluxVars)
320 {}
321
331 static void computeSourceEnergy(NumEqVector& source,
332 const Element& element,
333 const FVElementGeometry& fvGeometry,
334 const ElementVolumeVariables& elemVolVars,
335 const SubControlVolume &scv)
336 {
337 // specialization for 2 fluid phases
338 const auto &volVars = elemVolVars[scv];
339
340 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
341 const Scalar aws = volVars.interfacialArea(phase0Idx, sPhaseIdx);
342 const Scalar ans = volVars.interfacialArea(phase1Idx, sPhaseIdx);
343
344 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
345 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
346 const Scalar Ts = volVars.temperatureSolid();
347
348 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
349 const Scalar lambdaNonwetting = volVars.fluidThermalConductivity(phase1Idx);
350 const Scalar lambdaSolid = volVars.solidThermalConductivity();
351
352 const Scalar lambdaWN = harmonicMean(lambdaWetting, lambdaNonwetting);
353 const Scalar lambdaWS = harmonicMean(lambdaWetting, lambdaSolid);
354 const Scalar lambdaNS = harmonicMean(lambdaNonwetting, lambdaSolid);
355
356 const Scalar characteristicLength = volVars.characteristicLength() ;
357 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
358
359 const Scalar nusseltWN = harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
360 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
361 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
362
363 const Scalar wettingToNonwettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * awn * lambdaWN * nusseltWN ;
364 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * aws * lambdaWS * nusseltWS ;
365 const Scalar nonwettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * ans * lambdaNS * nusseltNS ;
366
367 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
368 {
369 switch (phaseIdx)
370 {
371 case phase0Idx:
372 source[energyEq0Idx + phaseIdx] += ( - wettingToNonwettingEnergyExchange - wettingToSolidEnergyExchange);
373 break;
374 case phase1Idx:
375 source[energyEq0Idx + phaseIdx] += (+ wettingToNonwettingEnergyExchange - nonwettingToSolidEnergyExchange);
376 break;
377 case sPhaseIdx:
378 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonwettingToSolidEnergyExchange);
379 break;
380 default:
381 DUNE_THROW(Dune::NotImplemented,
382 "wrong index");
383 } // end switch
384
385
386 using std::isfinite;
387 if (!isfinite(source[energyEq0Idx + phaseIdx]))
388 DUNE_THROW(NumericalProblem, "Calculated non-finite source, " << "Tw="<< Tw << " Tn="<< Tn<< " Ts="<< Ts);
389 }// end phases
390
391 // we only need to do this for when there is more than 1 fluid phase
392 if (enableChemicalNonEquilibrium)
393 {
394 // Here comes the catch: We are not doing energy conservation for the whole
395 // system, but rather for each individual phase.
396 // -> Therefore the energy fluxes over each phase boundary need be
397 // individually accounted for.
398 // -> Each particle crossing a phase boundary does carry some mass and
399 // thus energy!
400 // -> Therefore, this contribution needs to be added.
401 // -> the particle always brings the energy of the originating phase.
402 // -> Energy advectivly transported into a phase = the moles of a component that go into a phase
403 // * molMass * enthalpy of the component in the *originating* phase
404
405 const auto& fluidState = volVars.fluidState();
406
407 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
408 {
409 switch (phaseIdx)
410 {
411 case phase0Idx:
412 //sum up the transfered energy by the components into the wetting phase
413 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
414 {
415 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
416 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
417 * FluidSystem::molarMass(compIdx)
418 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
419 }
420 break;
421 case phase1Idx:
422 //sum up the transfered energy by the components into the nonwetting phase
423 for(int compIdx =0; compIdx<numComponents; ++compIdx)
424 {
425 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
426 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
427 * FluidSystem::molarMass(compIdx)
428 *FluidSystem::componentEnthalpy(fluidState, phase0Idx, compIdx));
429 }
430 break;
431 case sPhaseIdx:
432 break; // no sorption
433 default:
434 DUNE_THROW(Dune::NotImplemented,
435 "wrong index");
436 } // end switch
437 } // end phases
438 } // EnableChemicalNonEquilibrium
439 } // end source
440};
441} // end namespace Dumux
442
443#endif
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
A helper to deduce a vector with the same size as numbers of equations.
Some exceptions thrown in DuMux
Provides 3rd order polynomial splines.
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:69
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:46
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:44
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:171
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:89
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:109
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:74
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:152
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:114
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:318
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:331
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:266
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:299
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:254
Declares all properties used in Dumux.