version 3.8
freeflow/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//
12#ifndef DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH
13#define DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH
14
15#include <cmath>
16
19
20namespace Dumux {
21
22// forward declaration
23template<class GridGeometry, class FluxVariables, class DiscretizationMethod, bool enableEneryBalance, bool isCompositional>
25
31template<class GridGeometry, class FluxVariables, bool enableEneryBalance, bool isCompositional>
34 FluxVariables,
35 typename GridGeometry::DiscretizationMethod,
36 enableEneryBalance, isCompositional>;
37
42template<class GridGeometry, class FluxVariables, class DiscretizationMethod, bool isCompositional>
43class FreeFlowEnergyLocalResidualImplementation<GridGeometry, FluxVariables, DiscretizationMethod, false, isCompositional>
44{
45public:
46
48 template <typename... Args>
49 static void fluidPhaseStorage(Args&&... args)
50 {}
51
53 template <typename... Args>
54 static void heatFlux(Args&&... args)
55 {}
56};
57
62template<class GridGeometry, class FluxVariables>
64 FluxVariables,
65 DiscretizationMethods::Staggered,
66 true, false>
67{
68 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
69 using FVElementGeometry = typename GridGeometry::LocalView;
70 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
71
72public:
73
75 template<class NumEqVector, class VolumeVariables>
76 static void fluidPhaseStorage(NumEqVector& storage,
77 const VolumeVariables& volVars)
78 {
79 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
80 storage[localEnergyBalanceIdx] += volVars.density() * volVars.internalEnergy();
81 }
82
84 template<class NumEqVector, class Problem, class ElementVolumeVariables, class ElementFaceVariables>
85 static void heatFlux(NumEqVector& flux,
86 const Problem& problem,
87 const Element &element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const ElementFaceVariables& elemFaceVars,
91 const SubControlVolumeFace& scvf)
92 {
93 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
94
95 auto upwindTerm = [](const auto& volVars) { return volVars.density() * volVars.enthalpy(); };
96 flux[localEnergyBalanceIdx] += FluxVariables::advectiveFluxForCellCenter(problem,
97 fvGeometry,
98 elemVolVars,
99 elemFaceVars,
100 scvf,
101 upwindTerm);
102
103 flux[localEnergyBalanceIdx] += FluxVariables::HeatConductionType::flux(problem,
104 element,
105 fvGeometry,
106 elemVolVars,
107 scvf);
108 }
109};
110
115template<class GridGeometry, class FluxVariables>
117 FluxVariables,
118 DiscretizationMethods::Staggered,
119 true, true>
120 : public FreeFlowEnergyLocalResidualImplementation<GridGeometry,
121 FluxVariables,
122 DiscretizationMethods::Staggered,
123 true, false>
124{
126 FluxVariables,
128 true, false>;
129 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
130 using FVElementGeometry = typename GridGeometry::LocalView;
131 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
132
133public:
135 template<class NumEqVector, class Problem, class ElementVolumeVariables, class ElementFaceVariables>
136 static void heatFlux(NumEqVector& flux,
137 const Problem& problem,
138 const Element &element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const ElementFaceVariables& elemFaceVars,
142 const SubControlVolumeFace& scvf)
143 {
144 ParentType::heatFlux(flux, problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf);
145
146 static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1;
147 auto diffusiveFlux = FluxVariables::MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf);
148 for (int compIdx = 0; compIdx < FluxVariables::numComponents; ++compIdx)
149 {
150 // define the upstream direction according to the sign of the diffusive flux
151 using std::signbit;
152 const bool insideIsUpstream = !signbit(diffusiveFlux[compIdx]);
153 const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()];
154
155 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
156 flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx);
157 else
158 flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx)* elemVolVars[scvf.insideScvIdx()].molarMass(compIdx);
159 }
160 }
161};
162
163} // end namespace Dumux
164
165#endif
static void heatFlux(NumEqVector &flux, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf)
The convective and conductive heat fluxes in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:136
Specialization for staggered one-phase, non-isothermal models.
Definition: freeflow/nonisothermal/localresidual.hh:67
static void heatFlux(NumEqVector &flux, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf)
The convective and conductive heat fluxes in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:85
static void fluidPhaseStorage(NumEqVector &storage, const VolumeVariables &volVars)
The energy storage in the fluid phase.
Definition: freeflow/nonisothermal/localresidual.hh:76
static void heatFlux(Args &&... args)
do nothing for the isothermal case
Definition: freeflow/nonisothermal/localresidual.hh:54
static void fluidPhaseStorage(Args &&... args)
do nothing for the isothermal case
Definition: freeflow/nonisothermal/localresidual.hh:49
Definition: freeflow/nonisothermal/localresidual.hh:24
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
The available discretization methods in Dumux.
Definition: adapt.hh:17
The reference frameworks and formulations available for splitting total fluxes into a advective and d...