ExaDG
Loading...
Searching...
No Matches
statistics_manager.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_STATISTICS_MANAGER_H_
23#define INCLUDE_EXADG_POSTPROCESSOR_STATISTICS_MANAGER_H_
24
25// deal.II
26#include <deal.II/dofs/dof_handler.h>
27#include <deal.II/lac/la_parallel_vector.h>
28
29// ExaDG
30#include <exadg/postprocessor/time_control_statistics.h>
31#include <exadg/utilities/print_functions.h>
32
33namespace ExaDG
34{
35// turbulent channel data
36
38{
40 : cells_are_stretched(false),
41 viscosity(1.0),
42 density(1.0),
43 directory("output/"),
44 filename("channel")
45 {
46 }
47
48 void
49 print(dealii::ConditionalOStream & pcout)
50 {
51 if(time_control_data_statistics.time_control_data.is_active)
52 {
53 pcout << " Turbulent channel statistics:" << std::endl;
54
55 // only implemented for unsteady problem
56 pcout << " Time control:" << std::endl;
57 time_control_data_statistics.print(pcout, true /*unsteady*/);
58
59 print_parameter(pcout, "Cells are stretched", cells_are_stretched);
60 print_parameter(pcout, "Dynamic viscosity", viscosity);
61 print_parameter(pcout, "Density", density);
62 print_parameter(pcout, "Directory of output files", directory);
63 print_parameter(pcout, "Filename", filename);
64 }
65 }
66
67 TimeControlDataStatistics time_control_data_statistics;
68
69 // are cells stretched, i.e., is a volume manifold applied?
70 bool cells_are_stretched;
71
72 // dynamic viscosity
73 double viscosity;
74
75 // density
76 double density;
77
78 // directory and filename
79 std::string directory;
80 std::string filename;
81};
82
83template<int dim, typename Number>
85{
86public:
87 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
88
89 StatisticsManager(dealii::DoFHandler<dim> const & dof_handler_velocity,
90 dealii::Mapping<dim> const & mapping);
91
92 // The argument grid_transform indicates how the y-direction that is initially distributed from
93 // [0,1] is mapped to the actual grid. This must match the transformation applied to the
94 // triangulation, otherwise the identification of data will fail
95 void
96 setup(std::function<double(double const &)> const & grid_tranform,
97 TurbulentChannelData const & data);
98
99 void
100 evaluate(VectorType const & velocity, bool const unsteady);
101
102 void
103 write_output();
104
105 void
106 reset();
107
108 TimeControlStatistics time_control_statistics;
109
110private:
111 static unsigned int const n_points_y_per_cell_linear = 11;
112 unsigned int n_points_y_per_cell;
113
114 void
115 evaluate_statistics(VectorType const & velocity);
116
117 void
118 evaluate_statistics(const std::vector<VectorType> & velocity);
119
120 void
121 do_evaluate(const std::vector<VectorType const *> & velocity);
122
123 void
124 do_write_output(std::string const filename, double const dynamic_viscosity, double const density);
125
126 dealii::DoFHandler<dim> const & dof_handler;
127 dealii::Mapping<dim> const & mapping;
128 MPI_Comm mpi_comm;
129
130 // vector of y-coordinates at which statistical quantities are computed
131 std::vector<double> y_glob;
132
133 // mean velocity <u_i>, i=1,...,d (for all y-coordinates)
134 std::vector<std::vector<double>> vel_glob;
135
136 // square velocity <u_i²>, i=1,...,d (for all y-coordinates)
137 std::vector<std::vector<double>> velsq_glob;
138
139 // <u_1*u_2> = <u*v> (for all y-coordinates)
140 std::vector<double> veluv_glob;
141
142 // number of samples
143 int number_of_samples;
144
145 bool write_final_output;
146
148};
149
150} // namespace ExaDG
151
152#endif
Definition statistics_manager.h:85
Definition time_control_statistics.h:44
Definition driver.cpp:33
Definition time_control_statistics.h:31
Definition statistics_manager.h:38