version 3.8
porousmediumflow/nonisothermal/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_LOCAL_RESIDUAL_HH
15#define DUMUX_ENERGY_LOCAL_RESIDUAL_HH
16
19
20namespace Dumux {
21
22// forward declaration
23template<class TypeTag, bool enableEneryBalance>
25
26template<class TypeTag>
28
33template<class TypeTag>
35{
39 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
40 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
42
43public:
52 static void fluidPhaseStorage(NumEqVector& storage,
53 const SubControlVolume& scv,
54 const VolumeVariables& volVars,
55 int phaseIdx)
56 {}
57
65 static void solidPhaseStorage(NumEqVector& storage,
66 const SubControlVolume& scv,
67 const VolumeVariables& volVars)
68 {}
69
77 static void heatConvectionFlux(NumEqVector& flux,
78 FluxVariables& fluxVars,
79 int phaseIdx)
80 {}
81
88 static void heatConductionFlux(NumEqVector& flux,
89 FluxVariables& fluxVars)
90 {}
91
98 static void heatDispersionFlux(NumEqVector& flux,
99 FluxVariables& fluxVars)
100 {}
101};
102
107template<class TypeTag>
109{
113 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
114 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
117 using Element = typename GridView::template Codim<0>::Entity;
118 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
120 using Indices = typename ModelTraits::Indices;
121
122 static constexpr int numPhases = ModelTraits::numFluidPhases();
123 enum { energyEqIdx = Indices::energyEqIdx };
124
125public:
126
135 static void fluidPhaseStorage(NumEqVector& storage,
136 const SubControlVolume& scv,
137 const VolumeVariables& volVars,
138 int phaseIdx)
139 {
140 storage[energyEqIdx] += volVars.porosity()
141 * volVars.density(phaseIdx)
142 * volVars.internalEnergy(phaseIdx)
143 * volVars.saturation(phaseIdx);
144 }
145
153 static void solidPhaseStorage(NumEqVector& storage,
154 const SubControlVolume& scv,
155 const VolumeVariables& volVars)
156 {
157 storage[energyEqIdx] += volVars.temperature()
158 * volVars.solidHeatCapacity()
159 * volVars.solidDensity()
160 * (1.0 - volVars.porosity());
161 }
162
170 static void heatConvectionFlux(NumEqVector& flux,
171 FluxVariables& fluxVars,
172 int phaseIdx)
173 {
174 auto upwindTerm = [phaseIdx](const auto& volVars)
175 { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.enthalpy(phaseIdx); };
176
177 flux[energyEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
178 }
179
186 static void heatConductionFlux(NumEqVector& flux,
187 FluxVariables& fluxVars)
188 {
189 flux[energyEqIdx] += fluxVars.heatConductionFlux();
190 }
191
198 static void heatDispersionFlux(NumEqVector& flux,
199 FluxVariables& fluxVars)
200 {
201
202 if constexpr (ModelTraits::enableThermalDispersion())
203 {
204 flux[energyEqIdx] += fluxVars.thermalDispersionFlux();
205 }
206 }
207
208
218 static void computeSourceEnergy(NumEqVector& source,
219 const Element& element,
220 const FVElementGeometry& fvGeometry,
221 const ElementVolumeVariables& elemVolVars,
222 const SubControlVolume &scv)
223 {}
224};
225
226} // end namespace Dumux
227
228#endif
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonisothermal/localresidual.hh:98
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonisothermal/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/nonisothermal/localresidual.hh:52
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonisothermal/localresidual.hh:77
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonisothermal/localresidual.hh:65
static void solidPhaseStorage(NumEqVector &storage, const SubControlVolume &scv, const VolumeVariables &volVars)
The energy storage in the solid matrix.
Definition: porousmediumflow/nonisothermal/localresidual.hh:153
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/nonisothermal/localresidual.hh:135
static void heatConductionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The diffusive energy fluxes.
Definition: porousmediumflow/nonisothermal/localresidual.hh:186
static void heatDispersionFlux(NumEqVector &flux, FluxVariables &fluxVars)
The dispersive energy fluxes.
Definition: porousmediumflow/nonisothermal/localresidual.hh:198
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/nonisothermal/localresidual.hh:218
static void heatConvectionFlux(NumEqVector &flux, FluxVariables &fluxVars, int phaseIdx)
The advective phase energy fluxes.
Definition: porousmediumflow/nonisothermal/localresidual.hh:170
Definition: porousmediumflow/nonisothermal/localresidual.hh:24
Defines all properties used in Dumux.
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.