ExaDG
Loading...
Searching...
No Matches
pointwise_output_generator_base.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_COMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_POINTWISE_OUTPUT_GENERATOR_BASE_H_
23#define INCLUDE_COMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_POINTWISE_OUTPUT_GENERATOR_BASE_H_
24
25
26// deal.II
27#include <deal.II/base/conditional_ostream.h>
28#include <deal.II/base/mpi_remote_point_evaluation.h>
29#include <deal.II/base/point.h>
30
31#ifdef DEAL_II_WITH_HDF5
32# include <deal.II/base/hdf5.h>
33#endif
34
35#include <deal.II/distributed/tria.h>
36
37#include <deal.II/fe/mapping.h>
38
39#include <deal.II/numerics/vector_tools.h>
40
41// ExaDG
42#include <exadg/postprocessor/time_control.h>
43#include <exadg/utilities/tensor_utilities.h>
44
45
46namespace ExaDG
47{
48template<int dim>
50{
51 using point_value_type = typename dealii::Point<dim>::value_type;
52
54
55 TimeControlData time_control_data;
56
57 std::string directory;
58 std::string filename;
59
60 bool update_points_before_evaluation;
61
62 std::vector<dealii::Point<dim>> evaluation_points;
63
64 void
65 print(dealii::ConditionalOStream & pcout) const;
66};
67
68template<int dim, typename Number>
70{
71public:
72 using point_value_type = typename PointwiseOutputDataBase<dim>::point_value_type;
73
74 TimeControl time_control;
75
76protected:
77 void
78 do_evaluate(std::function<void()> const & write_solution, double const time, bool const unsteady);
79
80 PointwiseOutputGeneratorBase(MPI_Comm const & comm);
81
82 virtual ~PointwiseOutputGeneratorBase() = default;
83
84 void
85 setup_base(dealii::Triangulation<dim> const & triangulation_in,
86 dealii::Mapping<dim> const & mapping_in,
87 PointwiseOutputDataBase<dim> const & pointwise_output_data_in);
88
89 void
90 add_quantity(std::string const & name, unsigned int const n_components);
91
92 template<int n_components>
93 void
94 write_quantity(std::string const & name,
95 std::vector<dealii::Tensor<1, n_components, Number>> const & values,
96 unsigned int const first_component)
97 {
98 unsigned int const components = name_to_components.at(name);
99 for(unsigned int comp = 0; comp < components; ++comp)
100 {
101 extract_component_from_tensors(componentwise_result, values, first_component + comp);
102 write_component((components == 1) ? name : name + std::to_string(comp), componentwise_result);
103 }
104 }
105
106 void
107 write_quantity(std::string const & name, std::vector<Number> const & values)
108 {
109 write_component(name, values);
110 }
111
112 template<int n_components>
113 [[nodiscard]] std::vector<
114 typename dealii::FEPointEvaluation<n_components, dim, dim, Number>::value_type>
115 compute_point_values(dealii::LinearAlgebra::distributed::Vector<Number> const & solution,
116 dealii::DoFHandler<dim> const & dof_handler) const
117 {
118 return dealii::VectorTools::point_values<n_components>(*remote_evaluator,
119 dof_handler,
120 solution);
121 }
122
123private:
124 void
125 setup_remote_evaluator();
126
127 void
128 reinit_remote_evaluator();
129
130 void
131 create_hdf5_file();
132
133 void
134 write_evaluation_points(std::string const & name);
135
136#ifdef DEAL_II_WITH_HDF5
137 void
138 add_evaluation_points_dataset(dealii::HDF5::Group & group, std::string const & name);
139
140 void
141 add_time_dataset(dealii::HDF5::Group & group, std::string const & name);
142#endif
143
144 void
145 write_time(double time);
146
147 template<typename ComponentwiseContainerType>
148 void
149 write_component(std::string const & name, ComponentwiseContainerType const & componentwise_result)
150 {
151#ifdef DEAL_II_WITH_HDF5
152 auto dataset = hdf5_file->open_dataset("PhysicalInformation/" + name);
153
154 std::vector<hsize_t> hyperslab_offset = {0, time_control.get_counter()};
155 std::vector<hsize_t> hyperslab_dim = {pointwise_output_data.evaluation_points.size(), 1};
156
157 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
158 dataset.write_hyperslab(componentwise_result, hyperslab_offset, hyperslab_dim);
159 else
160 dataset.template write_none<Number>();
161#else
162 (void)name;
163 (void)componentwise_result;
164 AssertThrow(false, dealii::ExcMessage("deal.II is not compiled with HDF5!"));
165#endif
166 }
167
168 MPI_Comm const mpi_comm;
169
170 dealii::SmartPointer<dealii::Triangulation<dim> const> triangulation;
171 dealii::SmartPointer<dealii::Mapping<dim> const> mapping;
172 PointwiseOutputDataBase<dim> pointwise_output_data;
173 dealii::Vector<Number> componentwise_result;
174 unsigned int n_out_samples;
175 std::shared_ptr<dealii::Utilities::MPI::RemotePointEvaluation<dim>> remote_evaluator;
176 bool first_evaluation;
177
178 std::map<std::string, unsigned int> name_to_components;
179
180#ifdef DEAL_II_WITH_HDF5
181 std::unique_ptr<dealii::HDF5::File> hdf5_file;
182#endif
183};
184
185} // namespace ExaDG
186
187#endif /* INCLUDE_COMPRESSIBLE_NAVIER_STOKES_POSTPROCESSOR_POINTWISE_OUTPUT_GENERATOR_BASE_H_*/
Definition pointwise_output_generator_base.h:70
Definition time_control.h:64
Definition driver.cpp:33
void extract_component_from_tensors(dealii::Vector< Number > &dst, std::vector< dealii::Tensor< n_components1, n_components2, Number > > const &values, unsigned int const comp1, unsigned int const comp2)
Definition tensor_utilities.h:39
Definition pointwise_output_generator_base.h:50
Definition time_control.h:40