ExaDG
Loading...
Searching...
No Matches
turbulence_model.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_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_TURBULENCE_MODEL_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_TURBULENCE_MODEL_H_
24
25// ExaDG
26#include <exadg/incompressible_navier_stokes/spatial_discretization/viscosity_model_base.h>
27
28namespace ExaDG
29{
30namespace IncNS
31{
32/*
33 * Turbulence model.
34 */
35template<int dim, typename Number>
36class TurbulenceModel : public ViscosityModelBase<dim, Number>
37{
38private:
41
42 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
43
44 typedef dealii::VectorizedArray<Number> scalar;
45 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
46
47 typedef std::pair<unsigned int, unsigned int> Range;
48
49 typedef CellIntegrator<dim, dim, Number> CellIntegratorU;
50 typedef FaceIntegrator<dim, dim, Number> FaceIntegratorU;
51
52public:
53 /*
54 * Constructor.
55 */
57
58 /*
59 * Destructor.
60 */
61 virtual ~TurbulenceModel();
62
63 /*
64 * Initialization function.
65 */
66 void
67 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
68 dealii::Mapping<dim> const & mapping_in,
69 std::shared_ptr<Operators::ViscousKernel<dim, Number>> viscous_kernel_in,
70 TurbulenceModelData const & turbulence_model_data_in,
71 unsigned int const dof_index_velocity_in);
72
77 void
78 set_viscosity(VectorType const & velocity) const final;
79
83 void
84 add_viscosity(VectorType const & velocity) const final;
85
89 void
90 calculate_filter_width(dealii::Mapping<dim> const & mapping);
91
92private:
93 void
94 cell_loop_set_coefficients(dealii::MatrixFree<dim, Number> const & data,
95 VectorType &,
96 VectorType const & src,
97 Range const & cell_range) const;
98
99 void
100 face_loop_set_coefficients(dealii::MatrixFree<dim, Number> const & data,
101 VectorType &,
102 VectorType const & src,
103 Range const & face_range) const;
104
105 void
106 boundary_face_loop_set_coefficients(dealii::MatrixFree<dim, Number> const & data,
107 VectorType &,
108 VectorType const & src,
109 Range const & face_range) const;
110
115 void
116 add_turbulent_viscosity(scalar & viscosity,
117 scalar const & filter_width,
118 tensor const & velocity_gradient,
119 double const & model_constant) const;
120
135 void
136 smagorinsky_model(scalar const & filter_width,
137 tensor const & velocity_gradient,
138 double const & C,
139 scalar & viscosity) const;
140
165 void
166 vreman_model(scalar const & filter_width,
167 tensor const & velocity_gradient,
168 double const & C,
169 scalar & viscosity) const;
170
196 void
197 wale_model(scalar const & filter_width,
198 tensor const & velocity_gradient,
199 double const & C,
200 scalar & viscosity) const;
201
219 void
220 sigma_model(scalar const & filter_width,
221 tensor const & velocity_gradient,
222 double const & C,
223 scalar & viscosity) const;
224
225 TurbulenceModelData turbulence_model_data;
226 dealii::AlignedVector<scalar> filter_width_vector;
227};
228
229} // namespace IncNS
230} // namespace ExaDG
231
232#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_TURBULENCE_MODEL_H_ */
Definition viscous_operator.h:49
Definition turbulence_model.h:37
void set_viscosity(VectorType const &velocity) const final
Definition turbulence_model.cpp:59
void add_viscosity(VectorType const &velocity) const final
Definition turbulence_model.cpp:68
void calculate_filter_width(dealii::Mapping< dim > const &mapping)
Definition turbulence_model.cpp:233
Definition viscosity_model_base.h:43
Definition driver.cpp:33
Definition viscosity_model_data.h:51