GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/test/freeflow/navierstokes/analyticalsolutionvectors.hh
Date: 2025-04-19 19:19:10
Exec Total Coverage
Lines: 58 58 100.0%
Functions: 18 18 100.0%
Branches: 54 104 51.9%

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-FileCopyrightText: 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 NavierStokesTests
10 * \copydoc Dumux::NavierStokesAnalyticalSolutionVectors
11 */
12 #ifndef DUMUX_TEST_FREEFLOW_NAVIERSTOKES_ANALYTICALSOLVECTORS_HH
13 #define DUMUX_TEST_FREEFLOW_NAVIERSTOKES_ANALYTICALSOLVECTORS_HH
14
15 #include <dune/common/fvector.hh>
16
17 #include <dumux/discretization/method.hh>
18
19 namespace Dumux {
20 /*!
21 * \ingroup NavierStokesTests
22 * \brief Creates and provides the analytical solution in a form that can be written into the vtk output files
23 */
24 template<class Problem, class Scalar = double>
25 class NavierStokesAnalyticalSolutionVectors
26 {
27 using GridGeometry = std::decay_t<decltype(std::declval<Problem>().gridGeometry())>;
28 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
29 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
30 using Indices = typename Problem::Indices;
31
32 public:
33 9 NavierStokesAnalyticalSolutionVectors(std::shared_ptr<const Problem> problem, Scalar tInitial = 0.0)
34
1/2
✓ Branch 2 taken 9 times.
✗ Branch 3 not taken.
9 : problem_(problem)
35
1/10
✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
9 { update(tInitial); }
36
37 /*!
38 * \brief Creates the analytical solution in a form that can be written into the vtk output files
39 */
40 9 void update(Scalar time = 0.0)
41 {
42 9 analyticalPressure_.resize(problem_->gridGeometry().numCellCenterDofs());
43 9 analyticalVelocity_.resize(problem_->gridGeometry().numCellCenterDofs());
44 9 analyticalVelocityAtDofs_.resize(problem_->gridGeometry().numFaceDofs());
45
46 9 auto fvGeometry = localView(problem_->gridGeometry());
47
48
1/2
✓ Branch 2 taken 74209 times.
✗ Branch 3 not taken.
222609 for (const auto& element : elements(problem_->gridGeometry().gridView()))
49 {
50 74200 fvGeometry.bindElement(element);
51
2/2
✓ Branch 0 taken 74200 times.
✓ Branch 1 taken 74200 times.
222600 for (const auto& scv : scvs(fvGeometry))
52 {
53 // velocities on faces
54
2/2
✓ Branch 0 taken 296800 times.
✓ Branch 1 taken 74200 times.
371000 for (const auto& scvf : scvfs(fvGeometry))
55 {
56 296800 analyticalVelocityAtDofs_[scvf.dofIndex()][scvf.directionIndex()] = problem_->analyticalSolution(scvf.center(), time)[Indices::velocity(scvf.directionIndex())];
57 }
58
59 74200 analyticalPressure_[scv.dofIndex()] = problem_->analyticalSolution(scv.dofPosition(), time)[Indices::pressureIdx];
60
61
2/2
✓ Branch 0 taken 148400 times.
✓ Branch 1 taken 74200 times.
222600 for (int dirIdx = 0; dirIdx < dimWorld; ++dirIdx)
62 148400 analyticalVelocity_[scv.dofIndex()][dirIdx] = problem_->analyticalSolution(scv.center(), time)[Indices::velocity(dirIdx)];
63 }
64 }
65 9 }
66
67 /*!
68 * \brief Returns the analytical solution for the pressure
69 */
70
1/2
✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
9 const std::vector<Scalar>& getAnalyticalPressureSolution() const
71 {
72 return analyticalPressure_;
73 }
74
75 /*!
76 * \brief Returns the analytical solution for the velocity
77 */
78
1/2
✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
9 const std::vector<VelocityVector>& getAnalyticalVelocitySolution() const
79 {
80 return analyticalVelocity_;
81 }
82
83 /*!
84 * \brief Returns the analytical solution for the velocity at the faces
85 */
86
1/2
✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
9 const std::vector<VelocityVector>& getAnalyticalVelocitySolutionAtDofs() const
87 {
88 return analyticalVelocityAtDofs_;
89 }
90
91 private:
92 std::shared_ptr<const Problem> problem_;
93
94 std::vector<Scalar> analyticalPressure_;
95 std::vector<VelocityVector> analyticalVelocity_;
96 std::vector<VelocityVector> analyticalVelocityAtDofs_;
97 };
98
99 template<class Problem>
100 NavierStokesAnalyticalSolutionVectors(std::shared_ptr<Problem> p)
101 -> NavierStokesAnalyticalSolutionVectors<Problem>;
102
103 template<class Problem, class Scalar>
104 NavierStokesAnalyticalSolutionVectors(std::shared_ptr<Problem> p, Scalar t)
105 -> NavierStokesAnalyticalSolutionVectors<Problem, Scalar>;
106
107
108 namespace NavierStokesTest {
109
110 /*!
111 * \ingroup NavierStokesTests
112 * \brief Creates and provides the analytical solution in a form that can be written into the vtk output files
113 * \TODO make this the new NavierStokesAnalyticalSolutionVectors once all tests are ported to the new staggered
114 */
115 template<class MomentumProblem, class MassProblem, class Scalar = double>
116 class AnalyticalSolutionVectors
117 {
118 using MassGridGeometry = std::decay_t<decltype(std::declval<MassProblem>().gridGeometry())>;
119 using MomentumGridGeometry = std::decay_t<decltype(std::declval<MomentumProblem>().gridGeometry())>;
120 static constexpr int dimWorld = MassGridGeometry::GridView::dimensionworld;
121 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
122
123 using MassIndices = typename MassProblem::Indices;
124 using MomIndices = typename MomentumProblem::Indices;
125 public:
126 27 AnalyticalSolutionVectors(std::shared_ptr<const MomentumProblem> momentumProblem,
127 std::shared_ptr<const MassProblem> massProblem,
128 Scalar tInitial = 0.0)
129 27 : momentumProblem_(momentumProblem)
130
1/2
✓ Branch 2 taken 27 times.
✗ Branch 3 not taken.
27 , massProblem_(massProblem)
131
1/12
✓ Branch 1 taken 27 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
27 { update(tInitial); }
132
133 /*!
134 * \brief Creates the analytical solution in a form that can be written into the vtk output files
135 */
136 83 void update(Scalar time = 0.0)
137 {
138 83 analyticalPressure_.resize(massProblem_->gridGeometry().numDofs());
139 83 analyticalVelocity_.resize(massProblem_->gridGeometry().gridView().size(0));
140 83 analyticalVelocityAtDofs_.resize(momentumProblem_->gridGeometry().numDofs());
141
142 // cell-centers (pressure + velocity)
143 {
144
1/2
✓ Branch 1 taken 19 times.
✗ Branch 2 not taken.
83 auto fvGeometry = localView(massProblem_->gridGeometry());
145
1/2
✓ Branch 1 taken 19 times.
✗ Branch 2 not taken.
83 const auto gridView = massProblem_->gridGeometry().gridView();
146
6/10
✓ Branch 1 taken 505765 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 7 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 6736 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 6736 times.
✓ Branch 10 taken 7 times.
✓ Branch 13 taken 4 times.
✗ Branch 14 not taken.
1530626 for (const auto& element : elements(gridView))
147 {
148 // output velocity always on elements
149
2/4
✓ Branch 1 taken 6736 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 6736 times.
✗ Branch 5 not taken.
512418 const auto eIdx = massProblem_->gridGeometry().elementMapper().index(element);
150
1/2
✓ Branch 1 taken 6736 times.
✗ Branch 2 not taken.
512418 const auto center = element.geometry().center();
151
2/2
✓ Branch 0 taken 1024804 times.
✓ Branch 1 taken 512418 times.
1537222 for (int dirIdx = 0; dirIdx < dimWorld; ++dirIdx)
152 1024804 analyticalVelocity_[eIdx][dirIdx]
153
1/2
✓ Branch 1 taken 268800 times.
✗ Branch 2 not taken.
1024804 = momentumProblem_->analyticalSolution(center, time)[MomIndices::velocity(dirIdx)];
154
155 // this works for box and cc methods
156
1/2
✓ Branch 1 taken 4968 times.
✗ Branch 2 not taken.
514186 fvGeometry.bindElement(element);
157
3/3
✓ Branch 1 taken 512418 times.
✓ Branch 2 taken 505682 times.
✓ Branch 0 taken 13136 times.
1031236 for (const auto& scv : scvs(fvGeometry))
158 518818 analyticalPressure_[scv.dofIndex()]
159
1/2
✓ Branch 1 taken 134400 times.
✗ Branch 2 not taken.
518818 = massProblem_->analyticalSolution(scv.dofPosition(), time)[MassIndices::pressureIdx];
160 }
161 3 }
162
163 // dof positions (velocity)
164 {
165
1/2
✓ Branch 1 taken 16 times.
✗ Branch 2 not taken.
83 auto fvGeometry = localView(momentumProblem_->gridGeometry());
166
1/2
✓ Branch 1 taken 19 times.
✗ Branch 2 not taken.
83 const auto gridView = momentumProblem_->gridGeometry().gridView();
167
6/10
✓ Branch 1 taken 505765 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 7 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 6736 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 6736 times.
✓ Branch 10 taken 7 times.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
1530623 for (const auto& element : elements(gridView))
168 {
169
1/2
✓ Branch 1 taken 6736 times.
✗ Branch 2 not taken.
512418 fvGeometry.bindElement(element);
170
171 if constexpr (MomentumGridGeometry::discMethod == DiscretizationMethods::fcstaggered)
172
3/3
✓ Branch 0 taken 2029736 times.
✓ Branch 1 taken 514522 times.
✓ Branch 2 taken 1768 times.
2546026 for (const auto& scv : scvs(fvGeometry))
173 2036808 analyticalVelocityAtDofs_[scv.dofIndex()][scv.dofAxis()]
174
1/2
✓ Branch 1 taken 537600 times.
✗ Branch 2 not taken.
2036808 = momentumProblem_->analyticalSolution(scv.center(), time)[MomIndices::velocity(scv.dofAxis())];
175
176 else if constexpr (DiscretizationMethods::isCVFE<typename MomentumGridGeometry::DiscretizationMethod>)
177
2/2
✓ Branch 0 taken 9600 times.
✓ Branch 1 taken 3200 times.
12800 for (const auto& scv : scvs(fvGeometry))
178
2/2
✓ Branch 0 taken 19200 times.
✓ Branch 1 taken 9600 times.
28800 for (int dirIdx = 0; dirIdx < dimWorld; ++dirIdx)
179 19200 analyticalVelocityAtDofs_[scv.dofIndex()][dirIdx]
180 19200 = momentumProblem_->analyticalSolution(scv.dofPosition(), time)[MomIndices::velocity(dirIdx)];
181
182 else
183 DUNE_THROW(Dune::Exception, "Unknown discretization method: " << MomentumGridGeometry::discMethod);
184 }
185 1 }
186 83 }
187
188 /*!
189 * \brief Returns the analytical solution for the pressure
190 */
191
1/2
✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
27 const std::vector<Scalar>& analyticalPressureSolution() const
192 {
193
1/2
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
1 return analyticalPressure_;
194 }
195
196 /*!
197 * \brief Returns the analytical solution for the velocity
198 */
199
1/2
✓ Branch 1 taken 26 times.
✗ Branch 2 not taken.
27 const std::vector<VelocityVector>& analyticalVelocitySolution() const
200 {
201
1/2
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
1 return analyticalVelocity_;
202 }
203
204 /*!
205 * \brief Returns the analytical solution for the velocity at the faces
206 */
207 const std::vector<VelocityVector>& analyticalVelocitySolutionAtDofs() const
208 {
209 return analyticalVelocityAtDofs_;
210 }
211
212 private:
213 std::shared_ptr<const MomentumProblem> momentumProblem_;
214 std::shared_ptr<const MassProblem> massProblem_;
215
216 std::vector<Scalar> analyticalPressure_;
217 std::vector<VelocityVector> analyticalVelocity_;
218 std::vector<VelocityVector> analyticalVelocityAtDofs_;
219 };
220
221 template<class MomentumProblem, class MassProblem>
222 AnalyticalSolutionVectors(std::shared_ptr<MomentumProblem>, std::shared_ptr<MassProblem>)
223 -> AnalyticalSolutionVectors<MomentumProblem, MassProblem>;
224
225 template<class MomentumProblem, class MassProblem, class Scalar>
226 AnalyticalSolutionVectors(std::shared_ptr<MomentumProblem>, std::shared_ptr<MassProblem>, Scalar)
227 -> AnalyticalSolutionVectors<MomentumProblem, MassProblem, Scalar>;
228
229
230 } // end namespace NavierStokesTest
231 } // end namespace Dumux
232
233 #endif
234