ExaDG
Loading...
Searching...
No Matches
output_generator.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_OUTPUT_GENERATOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_OUTPUT_GENERATOR_H_
24
25#include <exadg/postprocessor/output_data_base.h>
26#include <exadg/postprocessor/solution_field.h>
27#include <exadg/postprocessor/time_control.h>
28
29namespace ExaDG
30{
31namespace IncNS
32{
34{
36 : write_vorticity(false),
37 write_divergence(false),
38 write_shear_rate(false),
39 write_velocity_magnitude(false),
40 write_vorticity_magnitude(false),
41 write_streamfunction(false),
42 write_q_criterion(false),
43 mean_velocity(TimeControlData()),
44 write_cfl(false),
45 write_aspect_ratio(false)
46 {
47 }
48
49 void
50 print(dealii::ConditionalOStream & pcout, bool unsteady)
51 {
52 OutputDataBase::print(pcout, unsteady);
53
54 print_parameter(pcout, "Write vorticity", write_vorticity);
55 print_parameter(pcout, "Write divergence", write_divergence);
56 print_parameter(pcout, "Write shear rate", write_shear_rate);
57 print_parameter(pcout, "Write velocity magnitude", write_velocity_magnitude);
58 print_parameter(pcout, "Write vorticity magnitude", write_vorticity_magnitude);
59 print_parameter(pcout, "Write streamfunction", write_streamfunction);
60 print_parameter(pcout, "Write Q criterion", write_q_criterion);
61
62 mean_velocity.print(pcout, unsteady);
63 }
64
65 // write vorticity of velocity field
66 bool write_vorticity;
67
68 // write divergence of velocity field
69 bool write_divergence;
70
71 // write shear rate in velocity field
72 bool write_shear_rate;
73
74 // write velocity magnitude
75 bool write_velocity_magnitude;
76
77 // write vorticity magnitude
78 bool write_vorticity_magnitude;
79
80 // Calculate streamfunction in order to visualize streamlines!
81 // Note that this option is only available in 2D!
82 // To calculate the streamfunction Psi, a Poisson equation is solved
83 // with homogeneous Dirichlet BC's. Accordingly, this approach can only be
84 // used for flow problems where the whole boundary is one streamline, e.g.,
85 // cavity-type flow problems where the velocity is tangential to the boundary
86 // on one part of the boundary and 0 (no-slip) on the rest of the boundary.
87 bool write_streamfunction;
88
89 // write Q criterion
90 bool write_q_criterion;
91
92 // Average velocity field over time for statistically steady, turbulent
93 // flow problems in order to visualize the time-averaged velocity field.
94 // Of course, this module can also be used for statistically unsteady problems,
95 // but in this case the mean velocity field is probably not meaningful.
96 TimeControlData mean_velocity;
97
98
99 // write cfl
100 bool write_cfl;
101
102 // write aspect ratio
103 bool write_aspect_ratio;
104};
105
106template<int dim, typename Number>
108
109template<int dim, typename Number>
111{
112public:
113 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
114
115 OutputGenerator(MPI_Comm const & comm);
116
117 void
118 setup(dealii::DoFHandler<dim> const & dof_handler_velocity_in,
119 dealii::DoFHandler<dim> const & dof_handler_pressure_in,
120 dealii::Mapping<dim> const & mapping_in,
121 OutputData const & output_data_in);
122
123 void
124 evaluate(VectorType const & velocity,
125 VectorType const & pressure,
126 std::vector<dealii::SmartPointer<SolutionField<dim, Number>>> const & additional_fields,
127 double const time,
128 bool const unsteady) const;
129
130 TimeControl time_control;
131
132private:
133 MPI_Comm const mpi_comm;
134
135 OutputData output_data;
136
137 dealii::SmartPointer<dealii::DoFHandler<dim> const> dof_handler_velocity;
138 dealii::SmartPointer<dealii::DoFHandler<dim> const> dof_handler_pressure;
139 dealii::SmartPointer<dealii::Mapping<dim> const> mapping;
140};
141
142} // namespace IncNS
143} // namespace ExaDG
144
145#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_OUTPUT_GENERATOR_H_ */
Definition output_generator.h:111
Definition spatial_operator_base.h:68
Definition solution_field.h:39
Definition time_control.h:64
Definition driver.cpp:33
Definition output_generator.h:34
Definition output_data_base.h:31
Definition time_control.h:40