ExaDG
Loading...
Searching...
No Matches
kinetic_energy_spectrum.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_POSTPROCESSOR_KINETIC_ENERGY_SPECTRUM_H_
23#define INCLUDE_EXADG_POSTPROCESSOR_KINETIC_ENERGY_SPECTRUM_H_
24
25// deal.II
26#include <deal.II/fe/fe_system.h>
27#include <deal.II/lac/la_parallel_vector.h>
28#include <deal.II/matrix_free/matrix_free.h>
29
30// ExaDG
31#include <exadg/postprocessor/time_control.h>
32#include <exadg/utilities/print_functions.h>
33
34namespace ExaDG
35{
36// forward declaration
37class DealSpectrumWrapper;
38
40{
42 : write_raw_data_to_files(false),
43 do_fftw(true),
44 directory("output/"),
45 filename("energy_spectrum"),
46 clear_file(true),
47 degree(0),
48 evaluation_points_per_cell(0),
49 exploit_symmetry(false),
50 n_cells_1d_coarse_grid(1),
51 refine_level(0),
52 length_symmetric_domain(dealii::numbers::PI)
53 {
54 }
55
56 void
57 print(dealii::ConditionalOStream & pcout) const
58 {
59 if(time_control_data.is_active)
60 {
61 // only implemented for unsteady problem
62 time_control_data.print(pcout, true /*unsteady*/);
63
64 pcout << std::endl << " Calculate kinetic energy spectrum:" << std::endl;
65 print_parameter(pcout, "Write raw data to files", write_raw_data_to_files);
66 print_parameter(pcout, "Do FFTW", do_fftw);
67 print_parameter(pcout, "Directory of output files", directory);
68 print_parameter(pcout, "Filename", filename);
69 print_parameter(pcout, "Clear file", clear_file);
70
71 print_parameter(pcout, "Evaluation points per cell", evaluation_points_per_cell);
72
73 print_parameter(pcout, "Exploit symmetry", exploit_symmetry);
74 if(exploit_symmetry)
75 {
76 print_parameter(pcout, "n_cells_1d_coarse_grid", n_cells_1d_coarse_grid);
77 print_parameter(pcout, "refine_level", refine_level);
78 print_parameter(pcout, "length_symmetric_domain", length_symmetric_domain);
79 }
80 }
81 }
82
83 TimeControlData time_control_data;
84
85 bool write_raw_data_to_files;
86 bool do_fftw;
87
88 // these parameters are only relevant if do_fftw = true
89 std::string directory;
90 std::string filename;
91 bool clear_file;
92
93 unsigned int degree;
94 unsigned int evaluation_points_per_cell;
95
96 // exploit symmetry for Navier-Stokes simulation and mirror dof-vector
97 // according to Taylor-Green symmetries for evaluation of energy spectrum.
98 bool exploit_symmetry;
99 unsigned int n_cells_1d_coarse_grid;
100 unsigned int refine_level;
101 double length_symmetric_domain;
102};
103
111template<int dim, typename Number>
113{
114public:
115 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
116
117 KineticEnergySpectrumCalculator(MPI_Comm const & mpi_comm);
118
119 void
120 setup(dealii::MatrixFree<dim, Number> const & matrix_free_data_in,
121 dealii::DoFHandler<dim> const & dof_handler_in,
122 KineticEnergySpectrumData const & data_in);
123
124 void
125 evaluate(VectorType const & velocity, double const time, bool const unsteady);
126
127 TimeControl time_control;
128
129private:
130 void
131 do_evaluate(VectorType const & velocity, double const time);
132
133 MPI_Comm const mpi_comm;
134
135 bool clear_files;
137 unsigned int const precision = 12;
138
139 std::shared_ptr<DealSpectrumWrapper> deal_spectrum_wrapper;
140
141 dealii::SmartPointer<dealii::DoFHandler<dim> const> dof_handler;
142
143 std::shared_ptr<VectorType> velocity_full;
144 std::shared_ptr<dealii::Triangulation<dim>> tria_full;
145 std::shared_ptr<dealii::FESystem<dim>> fe_full;
146 std::shared_ptr<dealii::DoFHandler<dim>> dof_handler_full;
147};
148} // namespace ExaDG
149
150#endif /* INCLUDE_EXADG_POSTPROCESSOR_KINETIC_ENERGY_SPECTRUM_H_ */
Definition kinetic_energy_spectrum.h:113
Definition time_control.h:64
Definition driver.cpp:33
Definition kinetic_energy_spectrum.h:40
Definition time_control.h:40