ExaDG
Loading...
Searching...
No Matches
write_output.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_WRITE_OUTPUT_H_
23#define INCLUDE_EXADG_POSTPROCESSOR_WRITE_OUTPUT_H_
24
25// C/C++
26#include <fstream>
27
28// deal.II
29#include <deal.II/base/bounding_box.h>
30#include <deal.II/grid/grid_generator.h>
31#include <deal.II/grid/grid_out.h>
32#include <deal.II/numerics/data_out.h>
33#include <deal.II/numerics/data_out_faces.h>
34#include <deal.II/particles/data_out.h>
35#include <deal.II/particles/particle_handler.h>
36
37namespace ExaDG
38{
39template<int dim>
40void
41write_surface_mesh(dealii::Triangulation<dim> const & triangulation,
42 dealii::Mapping<dim> const & mapping,
43 unsigned int const n_subdivisions,
44 std::string const & folder,
45 std::string const & file,
46 unsigned int const counter,
47 MPI_Comm const & mpi_comm)
48{
49 // write surface mesh only
50 dealii::DataOutFaces<dim> data_out_surface(true /*surface only*/);
51 data_out_surface.attach_triangulation(triangulation);
52 data_out_surface.build_patches(mapping, n_subdivisions);
53 data_out_surface.write_vtu_with_pvtu_record(folder, file + "_surface", counter, mpi_comm, 4);
54}
55
56template<int dim>
57void
58write_boundary_IDs(dealii::Triangulation<dim> const & triangulation,
59 std::string const & folder,
60 std::string const & file,
61 MPI_Comm const & mpi_communicator)
62{
63 unsigned int const rank = dealii::Utilities::MPI::this_mpi_process(mpi_communicator);
64 unsigned int const n_ranks = dealii::Utilities::MPI::n_mpi_processes(mpi_communicator);
65
66 unsigned int const n_digits = static_cast<int>(std::ceil(std::log10(std::fabs(n_ranks))));
67
68 std::string filename = folder + file + "_boundary_IDs" + "." +
69 dealii::Utilities::int_to_string(rank, n_digits) + ".vtk";
70 std::ofstream output(filename);
71
72 dealii::GridOut grid_out;
73 dealii::GridOutFlags::Vtk flags;
74 flags.output_cells = false;
75 flags.output_faces = true;
76 flags.output_edges = false;
77 flags.output_only_relevant = false;
78 grid_out.set_flags(flags);
79 grid_out.write_vtk(triangulation, output);
80}
81
82template<int dim>
83void
84write_grid(dealii::Triangulation<dim> const & triangulation,
85 dealii::Mapping<dim> const & mapping,
86 unsigned int const n_subdivisions,
87 std::string const & folder,
88 std::string const & file,
89 unsigned int const & counter,
90 MPI_Comm const & mpi_comm)
91{
92 std::string filename = file + "_grid";
93
94 dealii::DataOut<dim> data_out;
95
96 dealii::DataOutBase::VtkFlags flags;
97 flags.write_higher_order_cells = n_subdivisions > 1;
98 data_out.set_flags(flags);
99
100 data_out.attach_triangulation(triangulation);
101 data_out.build_patches(mapping, n_subdivisions, dealii::DataOut<dim>::curved_inner_cells);
102 data_out.write_vtu_with_pvtu_record(folder, filename, counter, mpi_comm, 4);
103}
104
105template<int dim>
106void
107write_points(dealii::Triangulation<dim> const & triangulation,
108 dealii::Mapping<dim> const & mapping,
109 std::vector<dealii::Point<dim>> const & points,
110 std::string const & folder,
111 std::string const & file,
112 unsigned int const counter,
113 MPI_Comm const & mpi_comm)
114{
115 std::string filename = file + "_points";
116
117 dealii::Particles::ParticleHandler<dim, dim> particle_handler(triangulation, mapping);
118
119 particle_handler.insert_particles(points);
120
121 dealii::Particles::DataOut<dim, dim> particle_output;
122 particle_output.build_patches(particle_handler);
123 particle_output.write_vtu_with_pvtu_record(folder, filename, counter, mpi_comm);
124}
125
126template<int dim>
127void
128write_points_in_dummy_triangulation(std::vector<dealii::Point<dim>> const & points,
129 std::string const & folder,
130 std::string const & file,
131 unsigned int const counter,
132 MPI_Comm const & mpi_comm)
133{
134 dealii::BoundingBox<dim> bounding_box(points);
135 auto const boundary_points =
136 bounding_box.create_extended(1e-3 * std::pow(bounding_box.volume(), 1.0 / ((double)dim)))
137 .get_boundary_points();
138
139 dealii::Triangulation<dim> particle_dummy_tria;
140 dealii::GridGenerator::hyper_rectangle(particle_dummy_tria,
141 boundary_points.first,
142 boundary_points.second);
143
144 dealii::MappingQGeneric<dim> particle_dummy_mapping(1 /* mapping_degree */);
145
146 write_points(
147 particle_dummy_tria, particle_dummy_mapping, points, folder, file, counter, mpi_comm);
148}
149
150} // namespace ExaDG
151
152#endif /* INCLUDE_EXADG_POSTPROCESSOR_WRITE_OUTPUT_H_ */
Definition driver.cpp:33