ExaDG
Loading...
Searching...
No Matches
line_plot_calculation_statistics_homogeneous.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_POSTPROCESSOR_LINE_PLOT_CALCULATION_STATISTICS_HOMOGENEOUS_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_LINE_PLOT_CALCULATION_STATISTICS_HOMOGENEOUS_H_
24
25// C++
26#include <fstream>
27
28// deal.II
29#include <deal.II/dofs/dof_handler.h>
30#include <deal.II/fe/mapping_q.h>
31#include <deal.II/lac/la_parallel_vector.h>
32
33// ExaDG
34#include <exadg/incompressible_navier_stokes/postprocessor/line_plot_data.h>
35#include <exadg/postprocessor/time_control.h>
36
37namespace ExaDG
38{
39namespace IncNS
40{
41/*
42 * This function calculates statistics along lines over time
43 * and one spatial, homogeneous direction (averaging_direction = {0,1,2}), e.g.,
44 * in the x-direction with a line in the y-z plane.
45 *
46 * NOTE: This functionality can only be used for hypercube meshes and for geometries/meshes for
47 * which the cells are aligned with the coordinate axis.
48 */
49
50template<int dim, typename Number>
52{
53public:
54 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
55
56 typedef typename std::vector<
57 std::pair<typename dealii::DoFHandler<dim>::active_cell_iterator, dealii::Point<dim>>>
58 TYPE;
59
60 LinePlotCalculatorStatisticsHomogeneous(dealii::DoFHandler<dim> const & dof_handler_velocity_in,
61 dealii::DoFHandler<dim> const & dof_handler_pressure_in,
62 dealii::Mapping<dim> const & mapping_in,
63 MPI_Comm const & mpi_comm_in);
64
65 void
66 setup(LinePlotDataStatistics<dim> const & data_in);
67
68 void
69 evaluate(VectorType const & velocity, VectorType const & pressure);
70
71 void
72 write_output() const;
73
74 TimeControlStatistics time_control_statistics;
75
76private:
77 void
78 print_headline(std::ofstream & f, unsigned int const number_of_samples) const;
79
80 void
81 do_evaluate(VectorType const & velocity, VectorType const & pressure);
82
83 void
84 do_evaluate_velocity(VectorType const & velocity,
85 Line<dim> const & line,
86 unsigned int const line_iterator);
87
88 void
89 do_evaluate_pressure(VectorType const & pressure,
90 Line<dim> const & line,
91 unsigned int const line_iterator);
92
93 void
94 average_pressure_for_given_point(VectorType const & pressure,
95 TYPE const & vector_cells_and_ref_points,
96 double & length_local,
97 double & pressure_local);
98
99 void
100 find_points_and_weights(dealii::Point<dim> const & point_in_ref_coord,
101 std::vector<dealii::Point<dim>> & points,
102 std::vector<double> & weights,
103 unsigned int const averaging_direction,
104 dealii::QGauss<1> const & gauss_1d);
105
106 void
107 do_write_output() const;
108
109 mutable bool clear_files;
110
111 dealii::DoFHandler<dim> const & dof_handler_velocity;
112 dealii::DoFHandler<dim> const & dof_handler_pressure;
113 dealii::Mapping<dim> const & mapping;
114 MPI_Comm mpi_comm;
115
117
118 // Global points
119 std::vector<std::vector<dealii::Point<dim>>> global_points;
120
121 // For all lines: for all points along the line: list of all relevant cells and points in ref
122 // coordinates
123 std::vector<std::vector<std::vector<
124 std::pair<typename dealii::DoFHandler<dim>::active_cell_iterator, dealii::Point<dim>>>>>
125 cells_and_ref_points_velocity;
126
127 // For all lines: for all points along the line: list of all relevant cells and points in ref
128 // coordinates
129 std::vector<std::vector<std::vector<
130 std::pair<typename dealii::DoFHandler<dim>::active_cell_iterator, dealii::Point<dim>>>>>
131 cells_and_ref_points_pressure;
132
133 // For all lines: for pressure reference point: list of all relevant cells and points in ref
134 // coordinates
135 std::vector<std::vector<
136 std::pair<typename dealii::DoFHandler<dim>::active_cell_iterator, dealii::Point<dim>>>>
137 cells_and_ref_points_ref_pressure;
138
139 // number of samples for averaging in time
140 unsigned int number_of_samples;
141
142 // homogeneous direction for averaging in space
143 unsigned int averaging_direction;
144
145 // Velocity quantities
146 // For all lines: for all points along the line
147 std::vector<std::vector<dealii::Tensor<1, dim, double>>> velocity_global;
148
149 // Skin Friction quantities
150 // For all lines: for all points along the line
151 std::vector<std::vector<double>> wall_shear_global;
152
153 // Reynolds Stress quantities
154 // For all lines: for all points along the line
155 std::vector<std::vector<dealii::Tensor<2, dim, double>>> reynolds_global;
156
157 // Pressure quantities
158 // For all lines: for all points along the line
159 std::vector<std::vector<double>> pressure_global;
160 // For all lines
161 std::vector<double> reference_pressure_global;
162
163 // write final output
164 bool write_final_output;
165};
166
167} // namespace IncNS
168} // namespace ExaDG
169
170#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_LINE_PLOT_CALCULATION_STATISTICS_HOMOGENEOUS_H_ \
171 */
Definition line_plot_calculation_statistics_homogeneous.h:52
Definition time_control_statistics.h:44
Definition driver.cpp:33
Definition line_plot_data.h:227
Definition line_plot_data.h:93