ExaDG
Loading...
Searching...
No Matches
navier_stokes_calculators.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_OPERATORS_NAVIER_STOKES_CALCULATORS_H_
23#define INCLUDE_EXADG_OPERATORS_NAVIER_STOKES_CALCULATORS_H_
24
25// deal.II
26#include <deal.II/lac/la_parallel_vector.h>
27
28// ExaDG
29#include <exadg/incompressible_navier_stokes/spatial_discretization/operators/viscous_operator.h>
30#include <exadg/matrix_free/integrators.h>
31
32namespace ExaDG
33{
34template<int dim, typename Number>
35class DivergenceCalculator
36{
37public:
38 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
39
40 typedef DivergenceCalculator<dim, Number> This;
41
42 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
43 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
44
45 typedef dealii::VectorizedArray<Number> scalar;
46
47 DivergenceCalculator();
48
49 void
50 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
51 unsigned int const dof_index_vector_in,
52 unsigned int const dof_index_scalar_in,
53 unsigned int const quad_index_in);
54 void
55 compute_divergence(VectorType & dst, VectorType const & src) const;
56
57private:
58 void
59 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
60 VectorType & dst,
61 VectorType const & src,
62 std::pair<unsigned int, unsigned int> const & cell_range) const;
63
64 dealii::MatrixFree<dim, Number> const * matrix_free;
65
66 unsigned int dof_index_vector;
67 unsigned int dof_index_scalar;
68 unsigned int quad_index;
69};
70
71template<int dim, typename Number>
72class ShearRateCalculator
73{
74private:
75 typedef ShearRateCalculator<dim, Number> This;
76
77 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
78
79 typedef dealii::VectorizedArray<Number> scalar;
80 typedef dealii::SymmetricTensor<2, dim, dealii::VectorizedArray<Number>> symmetrictensor;
81
82 typedef std::pair<unsigned int, unsigned int> Range;
83
84 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
85 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
86
87public:
88 ShearRateCalculator();
89
90 void
91 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
92 unsigned int const dof_index_u_in,
93 unsigned int const dof_index_u_scalar_in,
94 unsigned int const quad_index_in);
95
96 void
97 compute_shear_rate(VectorType & dst, VectorType const & src) const;
98
99private:
100 void
101 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
102 VectorType & dst,
103 VectorType const & src,
104 Range const & cell_range) const;
105
106 dealii::MatrixFree<dim, Number> const * matrix_free;
107
108 unsigned int dof_index_u;
109 unsigned int dof_index_u_scalar;
110 unsigned int quad_index;
111};
112
113template<int dim, typename Number>
114class VorticityCalculator
115{
116public:
117 static unsigned int const number_vorticity_components = (dim == 2) ? 1 : dim;
118
119 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
120
121 typedef VorticityCalculator<dim, Number> This;
122
123 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
124
125 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
126
127 VorticityCalculator();
128
129 void
130 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
131 unsigned int const dof_index_in,
132 unsigned int const quad_index_in);
133
134 void
135 compute_vorticity(VectorType & dst, VectorType const & src) const;
136
137private:
138 void
139 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
140 VectorType & dst,
141 VectorType const & src,
142 std::pair<unsigned int, unsigned int> const & cell_range) const;
143
144 dealii::MatrixFree<dim, Number> const * matrix_free;
145 unsigned int dof_index;
146 unsigned int quad_index;
147};
148
149template<int dim, typename Number>
150class ViscosityCalculator
151{
152private:
153 typedef ViscosityCalculator<dim, Number> This;
154
155 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
156
157 typedef dealii::VectorizedArray<Number> scalar;
158
159 typedef std::pair<unsigned int, unsigned int> Range;
160
161 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
162
163public:
164 ViscosityCalculator();
165
166 void
167 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
168 unsigned int const dof_index_u_scalar_in,
169 unsigned int const quad_index_in,
170 IncNS::Operators::ViscousKernel<dim, Number> const & viscous_kernel_in);
171
172 void
173 get_viscosity(VectorType & dst, VectorType const & src) const;
174
175private:
176 void
177 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
178 VectorType & dst,
179 VectorType const & src,
180 Range const & cell_range) const;
181
182 dealii::MatrixFree<dim, Number> const * matrix_free;
183
184 unsigned int dof_index_u_scalar;
185 unsigned int quad_index;
186
188};
189
190template<int dim, typename Number>
191class MagnitudeCalculator
192{
193private:
194 typedef MagnitudeCalculator<dim, Number> This;
195
196 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
197
198 typedef dealii::VectorizedArray<Number> scalar;
199
200 typedef std::pair<unsigned int, unsigned int> Range;
201
202 typedef CellIntegrator<dim, dim, Number> IntegratorVector;
203 typedef CellIntegrator<dim, 1, Number> IntegratorScalar;
204
205public:
206 MagnitudeCalculator();
207
208 void
209 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
210 unsigned int const dof_index_u_in,
211 unsigned int const dof_index_u_scalar_in,
212 unsigned int const quad_index_in);
213
214 void
215 compute(VectorType & dst, VectorType const & src) const;
216
217private:
218 void
219 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
220 VectorType & dst,
221 VectorType const & src,
222 Range const & cell_range) const;
223
224 dealii::MatrixFree<dim, Number> const * matrix_free;
225
226 unsigned int dof_index_u;
227 unsigned int dof_index_u_scalar;
228 unsigned int quad_index;
229};
230
231template<int dim, typename Number>
232class QCriterionCalculator
233{
234private:
235 typedef QCriterionCalculator<dim, Number> This;
236
237 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
238
239 typedef dealii::VectorizedArray<Number> scalar;
240 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
241
242 typedef std::pair<unsigned int, unsigned int> Range;
243
244 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
245 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
246
247public:
248 QCriterionCalculator();
249
250 void
251 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
252 unsigned int const dof_index_u_in,
253 unsigned int const dof_index_u_scalar_in,
254 unsigned int const quad_index_in,
255 bool const compressible_flow);
256
257 void
258 compute(VectorType & dst, VectorType const & src) const;
259
260private:
261 void
262 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
263 VectorType & dst,
264 VectorType const & src,
265 Range const & cell_range) const;
266
267 dealii::MatrixFree<dim, Number> const * matrix_free;
268
269 unsigned int dof_index_u;
270 unsigned int dof_index_u_scalar;
271 unsigned int quad_index;
272 bool compressible_flow;
273};
274
275} // namespace ExaDG
276
277#endif /* INCLUDE_EXADG_OPERATORS_NAVIER_STOKES_CALCULATORS_H_ \
278 */
Definition viscous_operator.h:63
Definition driver.cpp:33