3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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>
33
34namespace Dumux {
35
41// forward declaration
42template <class TypeTag, int numEnergyEqFluid>
44
45template<class TypeTag>
46class EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
47{
51 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
52 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using Element = typename GridView::template Codim<0>::Entity;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
59
61 using Indices = typename ModelTraits::Indices;
62
63 static constexpr auto numEnergyEqFluid = ModelTraits::numEnergyEqFluid();
64 static constexpr auto numEnergyEqSolid = ModelTraits::numEnergyEqSolid();
65 static constexpr auto energyEq0Idx = Indices::energyEq0Idx;
66 static constexpr auto energyEqSolidIdx = Indices::energyEqSolidIdx;
67
68 static constexpr auto numPhases = ModelTraits::numFluidPhases();
69 static constexpr auto numComponents = ModelTraits::numFluidComponents();
70
71public:
73 static void fluidPhaseStorage(NumEqVector& storage,
74 const SubControlVolume& scv,
75 const VolumeVariables& volVars,
76 int phaseIdx)
77 {
78 //in case we have one energy equation for more than one fluid phase, add up parts on the one energy equation
79 storage[energyEq0Idx] += volVars.porosity()
80 * volVars.density(phaseIdx)
81 * volVars.internalEnergy(phaseIdx)
82 * volVars.saturation(phaseIdx);
83
84 }
85
86
88 static void solidPhaseStorage(NumEqVector& storage,
89 const SubControlVolume& scv,
90 const VolumeVariables& volVars)
91 {
92 // heat conduction for the fluid phases
93 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
94 {
95 storage[energyEqSolidIdx+sPhaseIdx] += volVars.temperatureSolid()
96 * volVars.solidHeatCapacity()
97 * volVars.solidDensity()
98 * (1.0 - volVars.porosity());
99 }
100 }
101
102
104 static void heatConvectionFlux(NumEqVector& flux,
105 FluxVariables& fluxVars,
106 int phaseIdx)
107 {
108 auto upwindTerm = [phaseIdx](const auto& volVars)
109 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
110
111 //in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
112 flux[energyEq0Idx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
113
114 //now add the diffusive part
115 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
116 const auto& elemVolVars = fluxVars.elemVolVars();
117 const auto& scvf = fluxVars.scvFace();
118 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
119 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
120
121 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
122 {
123 //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
124 if (phaseIdx == compIdx)
125 continue;
126 //we need the upwind enthalpy. Even better would be the componentEnthalpy
127 auto enthalpy = 0.0;
128 if (diffusiveFluxes[compIdx] > 0)
129 enthalpy += insideVolVars.enthalpy(phaseIdx);
130 else
131 enthalpy += outsideVolVars.enthalpy(phaseIdx);
132
133 //check for the reference system and adapt units of the diffusive flux accordingly.
134 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
135 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*enthalpy;
136 else
137 flux[energyEq0Idx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
138 }
139 }
140
142 static void heatConductionFlux(NumEqVector& flux,
143 FluxVariables& fluxVars)
144 {
145 //in case we have one energy equation for more than one fluid phase we use an effective law in the nonequilibrium fourierslaw
146 flux[energyEq0Idx] += fluxVars.heatConductionFlux(0);
147 //heat conduction for the solid phases
148 for(int sPhaseIdx = 0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
149 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
150 }
151
161 static void computeSourceEnergy(NumEqVector& source,
162 const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const SubControlVolume &scv)
166 {
167 // specialization for 2 fluid phases
168 const auto& volVars = elemVolVars[scv];
169 const Scalar characteristicLength = volVars.characteristicLength() ;
170
171 // interfacial area
172 // Shi & Wang, Transport in porous media (2011)
173 const Scalar as = volVars.fluidSolidInterfacialArea();
174
175 // temperature fluid is the same for both fluids
176 const Scalar TFluid = volVars.temperatureFluid(0);
177 const Scalar TSolid = volVars.temperatureSolid();
178
179 Scalar solidToFluidEnergyExchange ;
180
181 const Scalar fluidConductivity = volVars.fluidThermalConductivity(0) ;
182
183 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
184
185 solidToFluidEnergyExchange = factorEnergyTransfer * (TSolid - TFluid) / characteristicLength * as * fluidConductivity;
186
187 solidToFluidEnergyExchange *= volVars.nusseltNumber(0);
188
189 for(int energyEqIdx = 0; energyEqIdx < numEnergyEqFluid+numEnergyEqSolid; ++energyEqIdx)
190 {
191 switch (energyEqIdx)
192 {
193 case 0 :
194 source[energyEq0Idx + energyEqIdx] += solidToFluidEnergyExchange;
195 break;
196 case 1 :
197 source[energyEq0Idx + energyEqIdx] -= solidToFluidEnergyExchange;
198 break;
199 default:
200 DUNE_THROW(Dune::NotImplemented,
201 "wrong index");
202 } // end switch
203 } // end energyEqIdx
204 } // end source
205};
206
207template<class TypeTag>
208class EnergyLocalResidualNonEquilibrium<TypeTag, 2/*numEnergyEqFluid*/>
209: public EnergyLocalResidualNonEquilibrium<TypeTag, 1/*numEnergyEqFluid*/>
210{
214 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
215 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
220 using Element = typename GridView::template Codim<0>::Entity;
221 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
222 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
223
225 using Indices = typename ModelTraits::Indices;
226
227 enum { numPhases = ModelTraits::numFluidPhases() };
228 enum { numEnergyEqFluid = ModelTraits::numEnergyEqFluid() };
229 enum { numEnergyEqSolid = ModelTraits::numEnergyEqSolid() };
230 enum { energyEq0Idx = Indices::energyEq0Idx };
231 enum { energyEqSolidIdx = Indices::energyEqSolidIdx};
232 enum { conti0EqIdx = Indices::conti0EqIdx };
233
234 enum { numComponents = ModelTraits::numFluidComponents() };
235 enum { phase0Idx = FluidSystem::phase0Idx};
236 enum { phase1Idx = FluidSystem::phase1Idx};
237 enum { sPhaseIdx = numPhases};
238
239 static constexpr bool enableChemicalNonEquilibrium = ModelTraits::enableChemicalNonEquilibrium();
240
241public:
242
244 static void fluidPhaseStorage(NumEqVector& storage,
245 const SubControlVolume& scv,
246 const VolumeVariables& volVars,
247 int phaseIdx)
248 {
249 storage[energyEq0Idx+phaseIdx] += volVars.porosity()
250 * volVars.density(phaseIdx)
251 * volVars.internalEnergy(phaseIdx)
252 * volVars.saturation(phaseIdx);
253 }
254
256 static void heatConvectionFlux(NumEqVector& flux,
257 FluxVariables& fluxVars,
258 int phaseIdx)
259 {
260 auto upwindTerm = [phaseIdx](const auto& volVars)
261 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
262
263 // in case we have one energy equation for more than one fluid phase, add up advective parts on the one energy equation
264 flux[energyEq0Idx+phaseIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
265
266 // add the diffusiv part
267 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
268 const auto& elemVolVars = fluxVars.elemVolVars();
269 const auto& scvf = fluxVars.scvFace();
270 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
271 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
272
273 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
274 {
275 // 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
276 if (phaseIdx == compIdx)
277 continue;
278 // we need the upwind enthalpy. Even better would be the componentEnthalpy
279 auto enthalpy = 0.0;
280 if (diffusiveFluxes[compIdx] > 0)
281 enthalpy += insideVolVars.enthalpy(phaseIdx);
282 else
283 enthalpy += outsideVolVars.enthalpy(phaseIdx);
284 flux[energyEq0Idx+phaseIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx)*enthalpy;
285 }
286 }
287
289 static void heatConductionFlux(NumEqVector& flux,
290 FluxVariables& fluxVars)
291 {
292 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
293 {
294 flux[energyEq0Idx+phaseIdx] += fluxVars.heatConductionFlux(phaseIdx);
295 }
296 for(int sPhaseIdx=0; sPhaseIdx<numEnergyEqSolid; ++sPhaseIdx)
297 {
298 flux[energyEqSolidIdx+sPhaseIdx] += fluxVars.heatConductionFlux(numPhases + sPhaseIdx);
299 }
300 }
310 static void computeSourceEnergy(NumEqVector& source,
311 const Element& element,
312 const FVElementGeometry& fvGeometry,
313 const ElementVolumeVariables& elemVolVars,
314 const SubControlVolume &scv)
315 {
316 // specialization for 2 fluid phases
317 const auto &volVars = elemVolVars[scv];
318
319 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
320 const Scalar aws = volVars.interfacialArea(phase0Idx, sPhaseIdx);
321 const Scalar ans = volVars.interfacialArea(phase1Idx, sPhaseIdx);
322
323 const Scalar Tw = volVars.temperatureFluid(phase0Idx);
324 const Scalar Tn = volVars.temperatureFluid(phase1Idx);
325 const Scalar Ts = volVars.temperatureSolid();
326
327 const Scalar lambdaWetting = volVars.fluidThermalConductivity(phase0Idx);
328 const Scalar lambdaNonWetting = volVars.fluidThermalConductivity(phase1Idx);
329 const Scalar lambdaSolid = volVars.solidThermalConductivity();
330
331 const Scalar lambdaWN = harmonicMean(lambdaWetting, lambdaNonWetting);
332 const Scalar lambdaWS = harmonicMean(lambdaWetting, lambdaSolid);
333 const Scalar lambdaNS = harmonicMean(lambdaNonWetting, lambdaSolid);
334
335 const Scalar characteristicLength = volVars.characteristicLength() ;
336 const Scalar factorEnergyTransfer = volVars.factorEnergyTransfer() ;
337
338 const Scalar nusseltWN = harmonicMean(volVars.nusseltNumber(phase0Idx), volVars.nusseltNumber(phase1Idx));
339 const Scalar nusseltWS = volVars.nusseltNumber(phase0Idx);
340 const Scalar nusseltNS = volVars.nusseltNumber(phase1Idx);
341
342 const Scalar wettingToNonWettingEnergyExchange = factorEnergyTransfer * (Tw - Tn) / characteristicLength * awn * lambdaWN * nusseltWN ;
343 const Scalar wettingToSolidEnergyExchange = factorEnergyTransfer * (Tw - Ts) / characteristicLength * aws * lambdaWS * nusseltWS ;
344 const Scalar nonWettingToSolidEnergyExchange = factorEnergyTransfer * (Tn - Ts) / characteristicLength * ans * lambdaNS * nusseltNS ;
345
346 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
347 {
348 switch (phaseIdx)
349 {
350 case phase0Idx:
351 source[energyEq0Idx + phaseIdx] += ( - wettingToNonWettingEnergyExchange - wettingToSolidEnergyExchange);
352 break;
353 case phase1Idx:
354 source[energyEq0Idx + phaseIdx] += (+ wettingToNonWettingEnergyExchange - nonWettingToSolidEnergyExchange);
355 break;
356 case sPhaseIdx:
357 source[energyEq0Idx + phaseIdx] += (+ wettingToSolidEnergyExchange + nonWettingToSolidEnergyExchange);
358 break;
359 default:
360 DUNE_THROW(Dune::NotImplemented,
361 "wrong index");
362 } // end switch
363
364
365 using std::isfinite;
366 if (!isfinite(source[energyEq0Idx + phaseIdx]))
367 DUNE_THROW(NumericalProblem, "Calculated non-finite source, " << "Tw="<< Tw << " Tn="<< Tn<< " Ts="<< Ts);
368 }// end phases
369
370 // we only need to do this for when there is more than 1 fluid phase
371 if (enableChemicalNonEquilibrium)
372 {
373 // Here comes the catch: We are not doing energy conservation for the whole
374 // system, but rather for each individual phase.
375 // -> Therefore the energy fluxes over each phase boundary need be
376 // individually accounted for.
377 // -> Each particle crossing a phase boundary does carry some mass and
378 // thus energy!
379 // -> Therefore, this contribution needs to be added.
380 // -> the particle always brings the energy of the originating phase.
381 // -> Energy advectivly transported into a phase = the moles of a component that go into a phase
382 // * molMass * enthalpy of the component in the *originating* phase
383
384 const auto& fluidState = volVars.fluidState();
385
386 for(int phaseIdx = 0; phaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++phaseIdx)
387 {
388 switch (phaseIdx)
389 {
390 case phase0Idx:
391 //sum up the transfered energy by the components into the wetting phase
392 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
393 {
394 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
395 source[energyEq0Idx + phaseIdx] += (source[eqIdx]
396 * FluidSystem::molarMass(compIdx)
397 * FluidSystem::componentEnthalpy(fluidState, phase1Idx, compIdx) );
398 }
399 break;
400 case phase1Idx:
401 //sum up the transfered energy by the components into the nonwetting 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, phase0Idx, compIdx));
408 }
409 break;
410 case sPhaseIdx:
411 break; // no sorption
412 default:
413 DUNE_THROW(Dune::NotImplemented,
414 "wrong index");
415 } // end switch
416 } // end phases
417 } // EnableChemicalNonEquilibrium
418 } // end source
419};
420} // end namespace Dumux
421
422#endif //
Some exceptions thrown in DuMux
Provides 3rd order polynomial splines.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:68
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:43
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:161
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:88
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:73
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:142
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:104
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:310
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:256
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonequilibrium/thermal/localresidual.hh:289
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:244
Declares all properties used in Dumux.