ExaDG
Loading...
Searching...
No Matches
operator.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_COMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_H_
23#define INCLUDE_EXADG_COMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATOR_H_
24
25// deal.II
26#include <deal.II/fe/fe_dgq.h>
27#include <deal.II/fe/fe_system.h>
28#include <deal.II/lac/la_parallel_vector.h>
29
30// ExaDG
31#include <exadg/compressible_navier_stokes/spatial_discretization/calculators.h>
32#include <exadg/compressible_navier_stokes/spatial_discretization/interface.h>
33#include <exadg/compressible_navier_stokes/spatial_discretization/kernels_and_operators.h>
34#include <exadg/compressible_navier_stokes/user_interface/boundary_descriptor.h>
35#include <exadg/compressible_navier_stokes/user_interface/field_functions.h>
36#include <exadg/compressible_navier_stokes/user_interface/parameters.h>
37#include <exadg/grid/grid.h>
38#include <exadg/matrix_free/matrix_free_data.h>
39#include <exadg/operators/inverse_mass_operator.h>
40#include <exadg/operators/navier_stokes_calculators.h>
41
42namespace ExaDG
43{
44namespace CompNS
45{
46template<int dim, typename Number>
47class Operator : public dealii::Subscriptor, public Interface::Operator<Number>
48{
49private:
50 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
51
52public:
53 Operator(std::shared_ptr<Grid<dim> const> grid,
54 std::shared_ptr<dealii::Mapping<dim> const> mapping,
55 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
56 std::shared_ptr<FieldFunctions<dim> const> field_functions,
57 Parameters const & param,
58 std::string const & field,
59 MPI_Comm const & mpi_comm);
60
61 void
62 fill_matrix_free_data(MatrixFreeData<dim, Number> & matrix_free_data) const;
63
67 void
68 setup();
69
75 void
76 setup(std::shared_ptr<dealii::MatrixFree<dim, Number> const> matrix_free,
77 std::shared_ptr<MatrixFreeData<dim, Number> const> matrix_free_data);
78
79 dealii::types::global_dof_index
80 get_number_of_dofs() const;
81
82 // initialization of DoF vectors
83 void
84 initialize_dof_vector(VectorType & src) const final;
85
86 void
87 initialize_dof_vector_scalar(VectorType & src) const;
88
89 void
90 initialize_dof_vector_dim_components(VectorType & src) const;
91
92 // set initial conditions
93 void
94 prescribe_initial_conditions(VectorType & src, double const time) const final;
95
96 /*
97 * This function is used in case of explicit time integration:
98 * This function evaluates the right-hand side operator, the
99 * convective and viscous terms (subsequently multiplied by -1.0 in order
100 * to shift these terms to the right-hand side of the equations)
101 * and finally applies the inverse mass operator.
102 */
103 void
104 evaluate(VectorType & dst, VectorType const & src, Number const time) const final;
105
106 void
107 evaluate_convective(VectorType & dst, VectorType const & src, Number const time) const;
108
109 void
110 evaluate_viscous(VectorType & dst, VectorType const & src, Number const time) const;
111
112 void
113 evaluate_convective_and_viscous(VectorType & dst,
114 VectorType const & src,
115 Number const time) const;
116
117 void
118 apply_inverse_mass(VectorType & dst, VectorType const & src) const;
119
120 // getters
121 dealii::MatrixFree<dim, Number> const &
122 get_matrix_free() const;
123
124 dealii::Mapping<dim> const &
125 get_mapping() const;
126
127 dealii::FiniteElement<dim> const &
128 get_fe() const;
129
130 dealii::DoFHandler<dim> const &
131 get_dof_handler() const;
132
133 dealii::DoFHandler<dim> const &
134 get_dof_handler_scalar() const;
135
136 dealii::DoFHandler<dim> const &
137 get_dof_handler_vector() const;
138
139 unsigned int
140 get_dof_index_vector() const;
141
142 unsigned int
143 get_dof_index_scalar() const;
144
145 unsigned int
146 get_quad_index_standard() const;
147
148 // pressure
149 void
150 compute_pressure(VectorType & dst, VectorType const & src) const;
151
152 // velocity
153 void
154 compute_velocity(VectorType & dst, VectorType const & src) const;
155
156 // temperature
157 void
158 compute_temperature(VectorType & dst, VectorType const & src) const;
159
160 // vorticity
161 void
162 compute_vorticity(VectorType & dst, VectorType const & src) const;
163
164 // divergence
165 void
166 compute_divergence(VectorType & dst, VectorType const & src) const;
167
168 // shear rate
169 void
170 compute_shear_rate(VectorType & dst, VectorType const & src) const;
171
172 double
173 get_wall_time_operator_evaluation() const final;
174
175 // global CFL criterion: calculates the time step size for a given global maximum velocity
176 double
177 calculate_time_step_cfl_global() const final;
178
179 // Calculate time step size according to diffusion term
180 double
181 calculate_time_step_diffusion() const final;
182
183private:
184 void
185 initialize_dof_handler_and_constraints();
186
187 void
188 setup_operators();
189
190 unsigned int
191 get_dof_index_all() const;
192
193 unsigned int
194 get_quad_index_overintegration_conv() const;
195
196 unsigned int
197 get_quad_index_overintegration_vis() const;
198
199 unsigned int
200 get_quad_index_l2_projections() const;
201
202 /*
203 * Grid
204 */
205 std::shared_ptr<Grid<dim> const> grid;
206
207 /*
208 * Mapping
209 */
210 std::shared_ptr<dealii::Mapping<dim> const> mapping;
211
212 /*
213 * User interface: Boundary conditions and field functions.
214 */
215 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor;
216 std::shared_ptr<FieldFunctions<dim> const> field_functions;
217
218 /*
219 * List of parameters.
220 */
221 Parameters const & param;
222
223 std::string const field;
224
225 /*
226 * Basic finite element ingredients.
227 */
228
229 std::shared_ptr<dealii::FiniteElement<dim>> fe; // all (dim+2) components: (rho, rho u, rho E)
230 std::shared_ptr<dealii::FiniteElement<dim>> fe_vector; // e.g. velocity
231 std::shared_ptr<dealii::FiniteElement<dim>> fe_scalar; // scalar quantity, e.g, pressure
232
233 // dealii::DoFHandler
234 dealii::DoFHandler<dim> dof_handler; // all (dim+2) components: (rho, rho u, rho E)
235 dealii::DoFHandler<dim> dof_handler_vector; // e.g. velocity
236 dealii::DoFHandler<dim> dof_handler_scalar; // scalar quantity, e.g, pressure
237
238 std::string const dof_index_all = "all_fields";
239 std::string const dof_index_vector = "vector";
240 std::string const dof_index_scalar = "scalar";
241
242 std::string const quad_index_standard = "standard";
243 std::string const quad_index_overintegration_conv = "overintegration_conv";
244 std::string const quad_index_overintegration_vis = "overintegration_vis";
245
246 std::string const quad_index_l2_projections = quad_index_standard;
247 // alternative: use more accurate over-integration strategy
248 // std::string const quad_index_l2_projections = quad_index_overintegration_conv;
249
250 /*
251 * Constraints.
252 */
253 dealii::AffineConstraints<Number> constraint;
254
255 /*
256 * Matrix-free operator evaluation.
257 */
258 std::shared_ptr<MatrixFreeData<dim, Number> const> matrix_free_data;
259 std::shared_ptr<dealii::MatrixFree<dim, Number> const> matrix_free;
260
261 /*
262 * Basic operators.
263 */
264 MassOperator<dim, Number> mass_operator;
265 BodyForceOperator<dim, Number> body_force_operator;
266 ConvectiveOperator<dim, Number> convective_operator;
267 ViscousOperator<dim, Number> viscous_operator;
268
269 /*
270 * Merged operators.
271 */
272 CombinedOperator<dim, Number> combined_operator;
273
275 InverseMassOperator<dim, dim, Number> inverse_mass_vector;
276 InverseMassOperator<dim, 1, Number> inverse_mass_scalar;
277
278 // L2 projections to calculate derived quantities
279 p_u_T_Calculator<dim, Number> p_u_T_calculator;
280 VorticityCalculator<dim, Number> vorticity_calculator;
281 DivergenceCalculator<dim, Number> divergence_calculator;
282 ShearRateCalculator<dim, Number> shear_rate_calculator;
283
284 /*
285 * MPI
286 */
287 MPI_Comm const mpi_comm;
288
289 /*
290 * Output to screen.
291 */
292 dealii::ConditionalOStream pcout;
293
294 // wall time for operator evaluation
295 mutable double wall_time_operator_evaluation;
296};
297
298} // namespace CompNS
299} // namespace ExaDG
300
301#endif /* INCLUDE_CONVECTION_DIFFUSION_DG_CONVECTION_DIFFUSION_OPERATION_H_ */
Definition kernels_and_operators.h:301
Definition kernels_and_operators.h:1689
Definition kernels_and_operators.h:510
Definition interface.h:35
Definition kernels_and_operators.h:417
Definition operator.h:48
void setup()
Definition operator.cpp:257
Definition parameters.h:36
Definition kernels_and_operators.h:949
Definition calculators.h:43
Definition navier_stokes_calculators.h:35
Definition grid.h:40
Definition inverse_mass_operator.h:53
Definition navier_stokes_calculators.h:72
Definition navier_stokes_calculators.h:114
Definition driver.cpp:33
Definition boundary_descriptor.h:109
Definition field_functions.h:31
Definition matrix_free_data.h:40