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 FreeflowNIModel | ||
10 | * \copydoc Dumux::FreeFlowEnergyLocalResidual | ||
11 | */ | ||
12 | #ifndef DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH | ||
13 | #define DUMUX_FREE_FLOW_ENERGY_LOCAL_RESIDUAL_HH | ||
14 | |||
15 | #include <cmath> | ||
16 | |||
17 | #include <dumux/discretization/method.hh> | ||
18 | #include <dumux/flux/referencesystemformulation.hh> | ||
19 | |||
20 | namespace Dumux { | ||
21 | |||
22 | // forward declaration | ||
23 | template<class GridGeometry, class FluxVariables, class DiscretizationMethod, bool enableEneryBalance, bool isCompositional> | ||
24 | class FreeFlowEnergyLocalResidualImplementation; | ||
25 | |||
26 | /*! | ||
27 | * \ingroup FreeflowNIModel | ||
28 | * \brief Element-wise calculation of the local residual for non-isothermal | ||
29 | * free-flow models | ||
30 | */ | ||
31 | template<class GridGeometry, class FluxVariables, bool enableEneryBalance, bool isCompositional> | ||
32 | using FreeFlowEnergyLocalResidual = | ||
33 | FreeFlowEnergyLocalResidualImplementation<GridGeometry, | ||
34 | FluxVariables, | ||
35 | typename GridGeometry::DiscretizationMethod, | ||
36 | enableEneryBalance, isCompositional>; | ||
37 | |||
38 | /*! | ||
39 | * \ingroup FreeflowNIModel | ||
40 | * \brief Specialization for isothermal models, does nothing | ||
41 | */ | ||
42 | template<class GridGeometry, class FluxVariables, class DiscretizationMethod, bool isCompositional> | ||
43 | class FreeFlowEnergyLocalResidualImplementation<GridGeometry, FluxVariables, DiscretizationMethod, false, isCompositional> | ||
44 | { | ||
45 | public: | ||
46 | |||
47 | //! do nothing for the isothermal case | ||
48 | template <typename... Args> | ||
49 | ✗ | static void fluidPhaseStorage(Args&&... args) | |
50 | ✗ | {} | |
51 | |||
52 | //! do nothing for the isothermal case | ||
53 | template <typename... Args> | ||
54 | ✗ | static void heatFlux(Args&&... args) | |
55 | ✗ | {} | |
56 | }; | ||
57 | |||
58 | /*! | ||
59 | * \ingroup FreeflowNIModel | ||
60 | * \brief Specialization for staggered one-phase, non-isothermal models | ||
61 | */ | ||
62 | template<class GridGeometry, class FluxVariables> | ||
63 | class FreeFlowEnergyLocalResidualImplementation<GridGeometry, | ||
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 | |||
72 | public: | ||
73 | |||
74 | //! The energy storage in the fluid phase | ||
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 | 72074560 | storage[localEnergyBalanceIdx] += volVars.density() * volVars.internalEnergy(); | |
81 | } | ||
82 | |||
83 | //! The convective and conductive heat fluxes in the fluid phase | ||
84 | template<class NumEqVector, class Problem, class ElementVolumeVariables, class ElementFaceVariables> | ||
85 | 61148748 | 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 | 366892488 | auto upwindTerm = [](const auto& volVars) { return volVars.density() * volVars.enthalpy(); }; | |
96 | 61148748 | flux[localEnergyBalanceIdx] += FluxVariables::advectiveFluxForCellCenter(problem, | |
97 | fvGeometry, | ||
98 | elemVolVars, | ||
99 | elemFaceVars, | ||
100 | scvf, | ||
101 | upwindTerm); | ||
102 | |||
103 | 61148748 | flux[localEnergyBalanceIdx] += FluxVariables::HeatConductionType::flux(problem, | |
104 | element, | ||
105 | fvGeometry, | ||
106 | elemVolVars, | ||
107 | scvf); | ||
108 | 61148748 | } | |
109 | }; | ||
110 | |||
111 | /*! | ||
112 | * \ingroup FreeflowNIModel | ||
113 | * \brief Specialization for staggered compositional, non-isothermal models | ||
114 | */ | ||
115 | template<class GridGeometry, class FluxVariables> | ||
116 | class FreeFlowEnergyLocalResidualImplementation<GridGeometry, | ||
117 | FluxVariables, | ||
118 | DiscretizationMethods::Staggered, | ||
119 | true, true> | ||
120 | : public FreeFlowEnergyLocalResidualImplementation<GridGeometry, | ||
121 | FluxVariables, | ||
122 | DiscretizationMethods::Staggered, | ||
123 | true, false> | ||
124 | { | ||
125 | using ParentType = FreeFlowEnergyLocalResidualImplementation<GridGeometry, | ||
126 | FluxVariables, | ||
127 | DiscretizationMethods::Staggered, | ||
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 | |||
133 | public: | ||
134 | //! The convective and conductive heat fluxes in the fluid phase | ||
135 | template<class NumEqVector, class Problem, class ElementVolumeVariables, class ElementFaceVariables> | ||
136 | 25111468 | 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 | 25111468 | ParentType::heatFlux(flux, problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf); | |
145 | |||
146 | static constexpr auto localEnergyBalanceIdx = NumEqVector::dimension - 1; | ||
147 | 25111468 | auto diffusiveFlux = FluxVariables::MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf); | |
148 |
2/2✓ Branch 0 taken 50222936 times.
✓ Branch 1 taken 25111468 times.
|
75334404 | 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 |
4/4✓ Branch 0 taken 25666269 times.
✓ Branch 1 taken 24556667 times.
✓ Branch 2 taken 25666269 times.
✓ Branch 3 taken 24556667 times.
|
100445872 | const bool insideIsUpstream = !signbit(diffusiveFlux[compIdx]); |
153 |
2/2✓ Branch 0 taken 25666269 times.
✓ Branch 1 taken 24556667 times.
|
50222936 | const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()]; |
154 | |||
155 | if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) | ||
156 | 100445872 | flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx); | |
157 | else | ||
158 | flux[localEnergyBalanceIdx] += diffusiveFlux[compIdx] * upstreamVolVars.componentEnthalpy(compIdx)* elemVolVars[scvf.insideScvIdx()].molarMass(compIdx); | ||
159 | } | ||
160 | 25111468 | } | |
161 | }; | ||
162 | |||
163 | } // end namespace Dumux | ||
164 | |||
165 | #endif | ||
166 |