Line |
Branch |
Exec |
Source |
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 |
|
|
// |
7 |
|
|
/*! |
8 |
|
|
* \file |
9 |
|
|
* \ingroup NIModel |
10 |
|
|
* \brief Element-wise calculation of the local residual for non-isothermal |
11 |
|
|
* fully implicit models. |
12 |
|
|
*/ |
13 |
|
|
|
14 |
|
|
#ifndef DUMUX_ENERGY_LOCAL_RESIDUAL_INCOMPRESSIBLE_HH |
15 |
|
|
#define DUMUX_ENERGY_LOCAL_RESIDUAL_INCOMPRESSIBLE_HH |
16 |
|
|
|
17 |
|
|
#include <dumux/common/properties.hh> |
18 |
|
|
#include <dumux/common/numeqvector.hh> |
19 |
|
|
#include "localresidual.hh" |
20 |
|
|
|
21 |
|
|
namespace Dumux { |
22 |
|
|
|
23 |
|
|
// forward declaration |
24 |
|
|
template<class TypeTag, bool enableEneryBalance> |
25 |
|
|
class EnergyLocalResidualIncompressibleImplementation; |
26 |
|
|
|
27 |
|
|
template<class TypeTag> |
28 |
|
|
using EnergyLocalResidualIncompressible = EnergyLocalResidualIncompressibleImplementation<TypeTag, GetPropType<TypeTag, Properties::ModelTraits>::enableEnergyBalance()>; |
29 |
|
|
|
30 |
|
|
/*! |
31 |
|
|
* \ingroup NIModel |
32 |
|
|
* \brief Element-wise calculation of the energy residual for non-isothermal problems. |
33 |
|
|
*/ |
34 |
|
|
template<class TypeTag> |
35 |
|
|
class EnergyLocalResidualIncompressibleImplementation<TypeTag, false>: public EnergyLocalResidualImplementation<TypeTag, false> |
36 |
|
|
{}; //energy balance not enabled |
37 |
|
|
|
38 |
|
|
/*! |
39 |
|
|
* \ingroup NIModel |
40 |
|
|
* \brief TODO docme! |
41 |
|
|
*/ |
42 |
|
|
template<class TypeTag> |
43 |
|
|
class EnergyLocalResidualIncompressibleImplementation<TypeTag, true>: public EnergyLocalResidualImplementation<TypeTag, true> |
44 |
|
|
{ |
45 |
|
|
using Scalar = GetPropType<TypeTag, Properties::Scalar>; |
46 |
|
|
using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; |
47 |
|
|
using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; |
48 |
|
|
using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; |
49 |
|
|
using SubControlVolume = typename FVElementGeometry::SubControlVolume; |
50 |
|
|
using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; |
51 |
|
|
using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; |
52 |
|
|
using Element = typename GridView::template Codim<0>::Entity; |
53 |
|
|
using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; |
54 |
|
|
using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; |
55 |
|
|
using Indices = typename ModelTraits::Indices; |
56 |
|
|
|
57 |
|
|
static constexpr int numPhases = ModelTraits::numFluidPhases(); |
58 |
|
|
enum { energyEqIdx = Indices::energyEqIdx }; |
59 |
|
|
|
60 |
|
|
public: |
61 |
|
|
/*! |
62 |
|
|
* \brief The advective phase energy fluxes for incompressible flow. |
63 |
|
|
* |
64 |
|
|
* Using specific internal energy $u$ instead of specific enthalpy \f$h\f$ for incompressible flow in convective flux |
65 |
|
|
* to account for otherwise neglected pressure work term (\f$\nabla p \cdot v\f$). |
66 |
|
|
* |
67 |
|
|
* Compressible formulation in EnergyLocalResidual (neglecting pressure work term (\f$\nabla p \cdot v\f$)) |
68 |
|
|
* is |
69 |
|
|
\f{align*}{ |
70 |
|
|
\frac{\partial}{\partial t} (\rho u) &= -\nabla \cdot (\rho v h) + \nabla \cdot (\lambda \nabla T) |
71 |
|
|
\f} |
72 |
|
|
* |
73 |
|
|
* Incompressible energy formulation is |
74 |
|
|
\f{align*}{ |
75 |
|
|
\frac{\partial}{\partial t} (\rho u) = -\nabla \cdot (\rho v u) + \nabla \cdot (\lambda \nabla T) |
76 |
|
|
\f} |
77 |
|
|
* |
78 |
|
|
* \param flux The flux |
79 |
|
|
* \param fluxVars The flux variables. |
80 |
|
|
* \param phaseIdx The phase index |
81 |
|
|
*/ |
82 |
|
|
static void heatConvectionFlux(NumEqVector& flux, |
83 |
|
|
FluxVariables& fluxVars, |
84 |
|
|
int phaseIdx) |
85 |
|
|
{ |
86 |
|
|
// internal energy used instead of enthalpy for incompressible flow |
87 |
|
69483 |
auto upwindTerm = [phaseIdx](const auto& volVars) |
88 |
|
185288 |
{ return volVars.density(phaseIdx)*volVars.mobility(phaseIdx)*volVars.internalEnergy(phaseIdx); }; |
89 |
|
|
|
90 |
|
23161 |
flux[energyEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm); |
91 |
|
|
} |
92 |
|
|
}; |
93 |
|
|
|
94 |
|
|
} // end namespace Dumux |
95 |
|
|
|
96 |
|
|
#endif |
97 |
|
|
|