ExaDG
Loading...
Searching...
No Matches
line_plot_calculation_statistics.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_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_LINE_PLOT_CALCULATION_STATISTICS_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 by averaging over time.
43 *
44 * Additionally, averaging in circumferential direction can be performed if desired.
45 *
46 * General assumptions:
47 *
48 * - we assume straight lines and an equidistant distribution of the evaluation points along each
49 * line.
50 *
51 * Assumptions for averaging in circumferential direction:
52 *
53 * - to define the plane in which we want to perform the averaging in circumferential direction,
54 * a normal vector has to be specified that has to be oriented normal to the straight line and
55 * normal to the averaging plane. To construct the sample points for averaging in circumferential
56 * direction, we assume that the first point of the line (line.begin) defines the center of the
57 * circle, and that the other points in circumferential direction can be constructed by rotating
58 * the vector from the center of the circle (line.begin) to the current point along the line
59 * around the normal vector.
60 */
61
62template<int dim, typename Number>
64{
65public:
66 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
67
68 typedef typename std::vector<
69 std::pair<typename dealii::DoFHandler<dim>::active_cell_iterator, dealii::Point<dim>>>
70 TYPE;
71
72 LinePlotCalculatorStatistics(dealii::DoFHandler<dim> const & dof_handler_velocity_in,
73 dealii::DoFHandler<dim> const & dof_handler_pressure_in,
74 dealii::Mapping<dim> const & mapping_in,
75 MPI_Comm const & mpi_comm_in);
76
77 void
78 setup(LinePlotDataStatistics<dim> const & data_in);
79
80 void
81 evaluate(VectorType const & velocity, VectorType const & pressure);
82
83 void
84 write_output() const;
85
86 TimeControlStatistics time_control_statistics;
87
88private:
89 void
90 print_headline(std::ofstream & f, unsigned int const number_of_samples) const
91 {
92 f << "number of samples: N = " << number_of_samples << std::endl;
93 }
94
95 void
96 initialize_cell_data(VectorType const & velocity, VectorType const & pressure);
97
98 void
99 do_evaluate(VectorType const & velocity, VectorType const & pressure);
100
101 void
102 do_evaluate_velocity(VectorType const & velocity,
103 Line<dim> const & line,
104 unsigned int const line_iterator);
105
106 void
107 do_evaluate_pressure(VectorType const & pressure,
108 Line<dim> const & line,
109 unsigned int const line_iterator);
110
111 void
112 do_write_output() const;
113
114 mutable bool clear_files;
115
116 dealii::DoFHandler<dim> const & dof_handler_velocity;
117 dealii::DoFHandler<dim> const & dof_handler_pressure;
118 dealii::Mapping<dim> const & mapping;
119 MPI_Comm mpi_comm;
120
122
123 // Global points
124 std::vector<std::vector<dealii::Point<dim>>> global_points;
125
126 bool cell_data_has_been_initialized;
127
128 // For all lines: for all points along the line: for all relevant cells: dof index of first dof of
129 // current cell and all shape function values
130 std::vector<std::vector<
131 std::vector<std::pair<std::vector<dealii::types::global_dof_index>, std::vector<Number>>>>>
132 cells_global_velocity;
133
134 // For all lines: for all points along the line: for all relevant cells: dof index of first dof of
135 // current cell and all shape function values
136 std::vector<std::vector<
137 std::vector<std::pair<std::vector<dealii::types::global_dof_index>, std::vector<Number>>>>>
138 cells_global_pressure;
139
140 // number of samples for averaging in time
141 unsigned int number_of_samples;
142
143 // Velocity quantities
144 // For all lines: for all points along the line
145 std::vector<std::vector<dealii::Tensor<1, dim, Number>>> velocity_global;
146
147 // Pressure quantities
148 // For all lines: for all points along the line
149 std::vector<std::vector<Number>> pressure_global;
150
151 bool write_final_output;
152};
153
154} // namespace IncNS
155} // namespace ExaDG
156
157#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_LINE_PLOT_CALCULATION_STATISTICS_H_ \
158 */
Definition line_plot_calculation_statistics.h:64
Definition time_control_statistics.h:44
Definition driver.cpp:33
Definition line_plot_data.h:227
Definition line_plot_data.h:93