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 EXADG_OPERATORS_NAVIER_STOKES_CALCULATORS_H_
23#define 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{
34/*
35 * Calculator for the divergence of a vector field u(x) defined as the scalar quantity
36
37 * div(u) := sum_i=1,...,dim { d u(i) / d x(i) }.
38 */
39template<int dim, typename Number>
40class DivergenceCalculator
41{
42public:
43 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
44
45 typedef DivergenceCalculator<dim, Number> This;
46
47 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
48 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
49
50 typedef dealii::VectorizedArray<Number> scalar;
51
52 DivergenceCalculator();
53
54 void
55 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
56 unsigned int const dof_index_vector_in,
57 unsigned int const dof_index_scalar_in,
58 unsigned int const quad_index_in);
59
60 /*
61 * Compute the right-hand side of an L2 projection of the divergence of the vector field.
62 */
63 void
64 compute_projection_rhs(VectorType & dst_scalar_valued,
65 VectorType const & src_vector_valued) const;
66
67private:
68 void
69 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
70 VectorType & dst_scalar_valued,
71 VectorType const & src_vector_valued,
72 std::pair<unsigned int, unsigned int> const & cell_range) const;
73
74 dealii::MatrixFree<dim, Number> const * matrix_free;
75
76 unsigned int dof_index_vector;
77 unsigned int dof_index_scalar;
78 unsigned int quad_index;
79};
80
81/*
82 * Calculator for the shear rate according to
83 * [Galdi et al., 2008, Hemodynamical Flows: Modeling, Analysis and Simulation]
84 * of a vector field u(x) defined as the scalar quantity
85 *
86 * sqrt(2 * trace(sym_grad(u)^2)) = sqrt(2 * sym_grad(u) : sym_grad(u))
87 * where
88 * sym_grad(u) := 0.5 * (grad(u) + grad(u)^T)
89 *
90 * is the symmetric part of the vector field's gradient.
91 */
92template<int dim, typename Number>
93class ShearRateCalculator
94{
95private:
96 typedef ShearRateCalculator<dim, Number> This;
97
98 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
99
100 typedef dealii::VectorizedArray<Number> scalar;
101 typedef dealii::SymmetricTensor<2, dim, dealii::VectorizedArray<Number>> symmetric_tensor;
102
103 typedef std::pair<unsigned int, unsigned int> Range;
104
105 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
106 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
107
108public:
109 ShearRateCalculator();
110
111 void
112 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
113 unsigned int const dof_index_vector_in,
114 unsigned int const dof_index_scalar_in,
115 unsigned int const quad_index_in);
116
117 /*
118 * Compute the right-hand side of an L2 projection of the shear rate of the vector field.
119 */
120 void
121 compute_projection_rhs(VectorType & dst_scalar_valued,
122 VectorType const & src_vector_valued) const;
123
124private:
125 void
126 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
127 VectorType & dst_scalar_valued,
128 VectorType const & src_vector_valued,
129 Range const & cell_range) const;
130
131 dealii::MatrixFree<dim, Number> const * matrix_free;
132
133 unsigned int dof_index_vector;
134 unsigned int dof_index_scalar;
135 unsigned int quad_index;
136};
137
138/*
139 * Calculator for the vorticity of a vector field u(x) defined as
140 *
141 * omega := curl(u)
142 *
143 * 2D : omega = - d u(1) / d x(2) + d u(2) / d x(1)
144 *
145 * 3D : omega(1) = d u(3) / d x(2) - d u(2) / d x(3)
146 * omega(2) = d u(1) / d x(3) - d u(3) / d x(1)
147 * omega(3) = d u(2) / d x(1) - d u(1) / d x(2)
148 *
149 * which is scalar-valued in 2D and vector-valued in 3D. Therefore, we always write into a
150 * vector-valued destination vector with `dim` components, filling only the first component in 2D.
151 * The second component in the solution vector in the 2D case corresponds to projecting the zero
152 * function into the finite element space.
153 */
154template<int dim, typename Number>
155class VorticityCalculator
156{
157public:
158 static unsigned int const number_vorticity_components = (dim == 2) ? 1 : dim;
159
160 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
161
162 typedef VorticityCalculator<dim, Number> This;
163
164 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
165
166 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
167
168 VorticityCalculator();
169
170 void
171 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
172 unsigned int const dof_index_vector_in,
173 unsigned int const quad_index_in);
174
175 /*
176 * Compute the right-hand side of an L2 projection of the vorticity of the vector field.
177 */
178 void
179 compute_projection_rhs(VectorType & dst_vector_valued,
180 VectorType const & src_vector_valued) const;
181
182private:
183 void
184 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
185 VectorType & dst_vector_valued,
186 VectorType const & src_vector_valued,
187 std::pair<unsigned int, unsigned int> const & cell_range) const;
188
189 dealii::MatrixFree<dim, Number> const * matrix_free;
190 unsigned int dof_index_vector;
191 unsigned int quad_index;
192};
193
194/*
195 * Calculator to project the viscosity coefficients stored in the `ViscousKernel` into the finite
196 * element space. Note that the viscosity is *not updated* but only accessed.
197 */
198template<int dim, typename Number>
199class ViscosityCalculator
200{
201private:
202 typedef ViscosityCalculator<dim, Number> This;
203
204 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
205
206 typedef dealii::VectorizedArray<Number> scalar;
207
208 typedef std::pair<unsigned int, unsigned int> Range;
209
210 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
211
212public:
213 ViscosityCalculator();
214
215 void
216 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
217 unsigned int const dof_index_scalar_in,
218 unsigned int const quad_index_in,
219 IncNS::Operators::ViscousKernel<dim, Number> const & viscous_kernel_in);
220
221 /*
222 * Compute the right-hand side of an L2 projection of the stored viscosity.
223 */
224 void
225 compute_projection_rhs(VectorType & dst_scalar_valued) const;
226
227private:
228 void
229 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
230 VectorType & dst_scalar_valued,
231 VectorType const & src_scalar_valued,
232 Range const & cell_range) const;
233
234 dealii::MatrixFree<dim, Number> const * matrix_free;
235
236 unsigned int dof_index_scalar;
237 unsigned int quad_index;
238
240};
241
242template<int dim, typename Number>
243class MagnitudeCalculator
244{
245private:
246 typedef MagnitudeCalculator<dim, Number> This;
247
248 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
249
250 typedef dealii::VectorizedArray<Number> scalar;
251
252 typedef std::pair<unsigned int, unsigned int> Range;
253
254 typedef CellIntegrator<dim, dim, Number> IntegratorVector;
255 typedef CellIntegrator<dim, 1, Number> IntegratorScalar;
256
257public:
258 MagnitudeCalculator();
259
260 void
261 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
262 unsigned int const dof_index_vector_in,
263 unsigned int const dof_index_scalar_in,
264 unsigned int const quad_index_in);
265
266 /*
267 * Compute the right-hand side of an L2 projection of the magnitude of the vector field.
268 */
269 void
270 compute_projection_rhs(VectorType & dst_scalar_valued,
271 VectorType const & src_vector_valued) const;
272
273private:
274 void
275 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
276 VectorType & dst_scalar_valued,
277 VectorType const & src_vector_valued,
278 Range const & cell_range) const;
279
280 dealii::MatrixFree<dim, Number> const * matrix_free;
281
282 unsigned int dof_index_vector;
283 unsigned int dof_index_scalar;
284 unsigned int quad_index;
285};
286
287/*
288 * Calculator to project the Q-criterion into the finite element space.
289 */
290template<int dim, typename Number>
291class QCriterionCalculator
292{
293private:
294 typedef QCriterionCalculator<dim, Number> This;
295
296 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
297
298 typedef dealii::VectorizedArray<Number> scalar;
299 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
300
301 typedef std::pair<unsigned int, unsigned int> Range;
302
303 typedef CellIntegrator<dim, dim, Number> CellIntegratorVector;
304 typedef CellIntegrator<dim, 1, Number> CellIntegratorScalar;
305
306public:
307 QCriterionCalculator();
308
309 void
310 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
311 unsigned int const dof_index_vector_in,
312 unsigned int const dof_index_scalar_in,
313 unsigned int const quad_index_in,
314 bool const compressible_flow);
315
316 /*
317 * Compute the right-hand side of an L2 projection of the QCriterion of the vector field.
318 */
319 void
320 compute_projection_rhs(VectorType & dst_scalar_valued,
321 VectorType const & src_vector_valued) const;
322
323private:
324 void
325 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
326 VectorType & dst_scalar_valued,
327 VectorType const & src_vector_valued,
328 Range const & cell_range) const;
329
330 dealii::MatrixFree<dim, Number> const * matrix_free;
331
332 unsigned int dof_index_vector;
333 unsigned int dof_index_scalar;
334 unsigned int quad_index;
335 bool compressible_flow;
336};
337
338} // namespace ExaDG
339
340#endif /* EXADG_OPERATORS_NAVIER_STOKES_CALCULATORS_H_ */
Definition viscous_operator.h:64
Definition driver.cpp:33