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 StaggeredDiscretization | ||
10 | * \copydoc Dumux::StaggeredFaceVariables | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH | ||
13 | #define DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH | ||
14 | |||
15 | #include <array> | ||
16 | #include <vector> | ||
17 | #include <type_traits> | ||
18 | |||
19 | namespace Dumux { | ||
20 | |||
21 | namespace Detail { | ||
22 | |||
23 | template<class Scalar, int upwindSchemeOrder> | ||
24 | 36000 | struct InAxisVelocities | |
25 | { | ||
26 | Scalar self = 0.0; | ||
27 | Scalar opposite = 0.0; | ||
28 | std::array<Scalar, upwindSchemeOrder-1> forward{}; | ||
29 | std::array<Scalar, upwindSchemeOrder-1> backward{}; | ||
30 | }; | ||
31 | |||
32 | template<class Scalar> | ||
33 | 679136 | struct InAxisVelocities<Scalar, 1> | |
34 | { | ||
35 | Scalar self = 0.0; | ||
36 | Scalar opposite = 0.0; | ||
37 | }; | ||
38 | |||
39 | } // end namespace Detail | ||
40 | |||
41 | /*! | ||
42 | * \ingroup StaggeredDiscretization | ||
43 | * \brief The face variables class for free flow staggered grid models. | ||
44 | * Contains all relevant velocities for the assembly of the momentum balance. | ||
45 | * When the upwindSchemeOrder is set to 2, additional velocities located at Dofs | ||
46 | * further from the central stencil will be added and used when calculating the | ||
47 | * advective term. When the order remains at 1, these velocities will not be provided. | ||
48 | */ | ||
49 | template<class FacePrimaryVariables, int dim, int upwindSchemeOrder> | ||
50 | 715136 | class StaggeredFaceVariables | |
51 | { | ||
52 | static constexpr int numPairs = (dim == 2) ? 2 : 4; | ||
53 | static constexpr bool useHigherOrder = upwindSchemeOrder > 1; | ||
54 | |||
55 | using Scalar = typename FacePrimaryVariables::block_type; | ||
56 | using InAxisVelocities = Detail::InAxisVelocities<Scalar, upwindSchemeOrder>; | ||
57 | |||
58 | public: | ||
59 | |||
60 | /*! | ||
61 | * \brief Partial update of the face variables. Only the face itself is considered. | ||
62 | * | ||
63 | * \param priVars The face-specific primary variales | ||
64 | */ | ||
65 | ✗ | void updateOwnFaceOnly(const FacePrimaryVariables& priVars) | |
66 | { | ||
67 |
2/2✓ Branch 0 taken 1651648 times.
✓ Branch 1 taken 7294256 times.
|
8945904 | inAxisVelocities_.self = priVars[0]; |
68 | ✗ | } | |
69 | |||
70 | /*! | ||
71 | * \brief Complete update of the face variables (i.e. velocities for free flow) | ||
72 | * for a given face | ||
73 | * | ||
74 | * \param faceSol The face-specific solution vector | ||
75 | * \param problem The problem | ||
76 | * \param element The element | ||
77 | * \param fvGeometry The finite-volume geometry | ||
78 | * \param scvf The sub-control volume face of interest | ||
79 | */ | ||
80 | template<class FaceSolution, class Problem, class Element, | ||
81 | class FVElementGeometry, class SubControlVolumeFace> | ||
82 | 76964816 | void update(const FaceSolution& faceSol, | |
83 | const Problem& problem, | ||
84 | const Element& element, | ||
85 | const FVElementGeometry& fvGeometry, | ||
86 | const SubControlVolumeFace& scvf) | ||
87 | { | ||
88 | static_assert(std::decay_t<decltype(faceSol[0])>::dimension == 1, | ||
89 | "\n\n\nVelocity primary variable must be a scalar value. \n\n Make sure to use\n\n ffSol = partial(sol, ffFaceIdx, ffCellCenterIdx);\n\n"); | ||
90 | |||
91 | 153929632 | inAxisVelocities_.self = faceSol[scvf.dofIndex()]; | |
92 | 153929632 | inAxisVelocities_.opposite = faceSol[scvf.dofIndexOpposingFace()]; | |
93 | |||
94 | 76964816 | addHigherOrderInAxisVelocities_(faceSol, scvf, std::integral_constant<bool, useHigherOrder>{}); | |
95 | |||
96 | // handle all sub faces | ||
97 |
3/3✓ Branch 0 taken 75318816 times.
✓ Branch 1 taken 152283632 times.
✓ Branch 2 taken 3292000 times.
|
230894448 | for (int i = 0; i < velocityParallel_.size(); ++i) |
98 | 461788896 | velocityParallel_[i].fill(0.0); | |
99 | |||
100 |
4/4✓ Branch 0 taken 153929632 times.
✓ Branch 1 taken 76964816 times.
✓ Branch 2 taken 153929632 times.
✓ Branch 3 taken 76964816 times.
|
307859264 | for (int i = 0; i < scvf.pairData().size(); ++i) |
101 | { | ||
102 | 153929632 | const auto& subFaceData = scvf.pairData(i); | |
103 | |||
104 | // treat the velocities normal to the face | ||
105 |
2/2✓ Branch 1 taken 146783280 times.
✓ Branch 2 taken 7146352 times.
|
153929632 | velocityLateralInside_[i] = faceSol[subFaceData.lateralPair.first]; |
106 | |||
107 |
4/4✓ Branch 0 taken 146783280 times.
✓ Branch 1 taken 7146352 times.
✓ Branch 2 taken 146783280 times.
✓ Branch 3 taken 7146352 times.
|
307859264 | if (scvf.hasOuterLateral(i)) |
108 | 146783280 | velocityLateralOutside_[i] = faceSol[subFaceData.lateralPair.second]; | |
109 | |||
110 | // treat the velocities parallel to the self face | ||
111 |
2/2✓ Branch 0 taken 157221632 times.
✓ Branch 1 taken 153929632 times.
|
311151264 | for (int j = 0; j < upwindSchemeOrder; j++) |
112 | { | ||
113 |
4/4✓ Branch 0 taken 14563892 times.
✓ Branch 1 taken 142657740 times.
✓ Branch 2 taken 14563892 times.
✓ Branch 3 taken 142657740 times.
|
314443264 | if (scvf.hasParallelNeighbor(i,j)) |
114 | 297720824 | velocityParallel_[i][j] = faceSol[subFaceData.parallelDofs[j]]; | |
115 | } | ||
116 | } | ||
117 | 76964816 | } | |
118 | |||
119 | /*! | ||
120 | * \brief Returns the velocity at the face itself | ||
121 | */ | ||
122 | ✗ | Scalar velocitySelf() const | |
123 | { | ||
124 | ✗ | return inAxisVelocities_.self; | |
125 | } | ||
126 | |||
127 | /*! | ||
128 | * \brief Returns the velocity at the opposing face | ||
129 | */ | ||
130 | ✗ | Scalar velocityOpposite() const | |
131 | { | ||
132 | ✗ | return inAxisVelocities_.opposite; | |
133 | } | ||
134 | |||
135 | /*! | ||
136 | * \brief Returns the velocity at a backward face | ||
137 | * | ||
138 | * \param backwardIdx The index describing how many faces backward this dof is from the opposite face | ||
139 | */ | ||
140 | template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0> | ||
141 | Scalar velocityBackward(const int backwardIdx) const | ||
142 | { | ||
143 | 1669136 | return inAxisVelocities_.backward[backwardIdx]; | |
144 | } | ||
145 | |||
146 | /*! | ||
147 | * \brief Returns the velocity at a forward face | ||
148 | * | ||
149 | * \param forwardIdx The index describing how many faces forward this dof is of the self face | ||
150 | */ | ||
151 | template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0> | ||
152 | Scalar velocityForward(const int forwardIdx) const | ||
153 | { | ||
154 | 2181152 | return inAxisVelocities_.forward[forwardIdx]; | |
155 | } | ||
156 | |||
157 | /*! | ||
158 | * \brief Returns the velocity at a parallel face | ||
159 | * | ||
160 | * \param localSubFaceIdx The local index of the subface | ||
161 | * \param parallelDegreeIdx The index describing how many faces parallel this dof is of the parallel face | ||
162 | */ | ||
163 | Scalar velocityParallel(const int localSubFaceIdx, const int parallelDegreeIdx = 0) const | ||
164 | { | ||
165 |
18/18✓ Branch 0 taken 123825619 times.
✓ Branch 1 taken 128256115 times.
✓ Branch 2 taken 123825619 times.
✓ Branch 3 taken 128256115 times.
✓ Branch 4 taken 123825619 times.
✓ Branch 5 taken 128256115 times.
✓ Branch 6 taken 1186089 times.
✓ Branch 7 taken 467792 times.
✓ Branch 8 taken 1186089 times.
✓ Branch 9 taken 467792 times.
✓ Branch 10 taken 1186089 times.
✓ Branch 11 taken 467792 times.
✓ Branch 12 taken 2115753 times.
✓ Branch 13 taken 53878 times.
✓ Branch 14 taken 2115753 times.
✓ Branch 15 taken 53878 times.
✓ Branch 16 taken 2115753 times.
✓ Branch 17 taken 53878 times.
|
1598415540 | return velocityParallel_[localSubFaceIdx][parallelDegreeIdx]; |
166 | } | ||
167 | |||
168 | /*! | ||
169 | * \brief Returns the velocity at the inner normal face | ||
170 | * | ||
171 | * \param localSubFaceIdx The local index of the subface | ||
172 | */ | ||
173 | Scalar velocityLateralInside(const int localSubFaceIdx) const | ||
174 | { | ||
175 | 1103976190 | return velocityLateralInside_[localSubFaceIdx]; | |
176 | } | ||
177 | |||
178 | /*! | ||
179 | * \brief Returns the velocity at the outer normal face | ||
180 | * | ||
181 | * \param localSubFaceIdx The local index of the subface | ||
182 | */ | ||
183 | Scalar velocityLateralOutside(const int localSubFaceIdx) const | ||
184 | { | ||
185 | 565763780 | return velocityLateralOutside_[localSubFaceIdx]; | |
186 | } | ||
187 | |||
188 | private: | ||
189 | |||
190 | template<class SolVector, class SubControlVolumeFace> | ||
191 | ✗ | void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::false_type) {} | |
192 | |||
193 | template<class SolVector, class SubControlVolumeFace> | ||
194 | 1646000 | void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::true_type) | |
195 | { | ||
196 | |||
197 | // treat the velocity forward of the self face i.e. the face that is | ||
198 | // forward wrt the self face by degree i | ||
199 |
2/2✓ Branch 0 taken 1646000 times.
✓ Branch 1 taken 1646000 times.
|
4938000 | for (int i = 0; i < scvf.axisData().inAxisForwardDofs.size(); i++) |
200 | { | ||
201 |
4/4✓ Branch 0 taken 27666 times.
✓ Branch 1 taken 1618334 times.
✓ Branch 2 taken 27666 times.
✓ Branch 3 taken 1618334 times.
|
3292000 | if (scvf.hasForwardNeighbor(i)) |
202 | 3236668 | inAxisVelocities_.forward[i]= faceSol[scvf.axisData().inAxisForwardDofs[i]]; | |
203 | } | ||
204 | |||
205 | // treat the velocity at the first backward face i.e. the face that is | ||
206 | // behind the opposite face by degree i | ||
207 |
2/2✓ Branch 0 taken 1646000 times.
✓ Branch 1 taken 1646000 times.
|
3292000 | for (int i = 0; i < scvf.axisData().inAxisBackwardDofs.size(); i++) |
208 | { | ||
209 |
4/4✓ Branch 0 taken 33266 times.
✓ Branch 1 taken 1612734 times.
✓ Branch 2 taken 33266 times.
✓ Branch 3 taken 1612734 times.
|
3292000 | if (scvf.hasBackwardNeighbor(i)) |
210 | 3225468 | inAxisVelocities_.backward[i] = faceSol[scvf.axisData().inAxisBackwardDofs[i]]; | |
211 | } | ||
212 | 1646000 | } | |
213 | |||
214 | InAxisVelocities inAxisVelocities_; | ||
215 | std::array<std::array<Scalar, upwindSchemeOrder>, numPairs> velocityParallel_; | ||
216 | std::array<Scalar, numPairs> velocityLateralInside_; | ||
217 | std::array<Scalar, numPairs> velocityLateralOutside_; | ||
218 | |||
219 | }; | ||
220 | |||
221 | } // end namespace Dumux | ||
222 | |||
223 | #endif | ||
224 |