ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
vtu_output.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2024 by the ryujin authors
4//
5
6#pragma once
7
9#include "vtu_output.h"
10
11#include <deal.II/base/function_parser.h>
12#include <deal.II/numerics/data_out.h>
13#include <deal.II/numerics/vector_tools.h>
14
15
16namespace ryujin
17{
18 using namespace dealii;
19
20
21 template <typename Description, int dim, typename Number>
23 const MPIEnsemble &mpi_ensemble,
24 const OfflineData<dim, Number> &offline_data,
25 const HyperbolicSystem &hyperbolic_system,
26 const ParabolicSystem &parabolic_system,
27 const Postprocessor<Description, dim, Number> &postprocessor,
28 const InitialPrecomputedVector &initial_precomputed,
29 const ScalarVector &alpha,
30 const ScalarVector &smoothness_indicators,
31 const std::string &subsection /*= "VTUOutput"*/)
32 : ParameterAcceptor(subsection)
33 , mpi_ensemble_(mpi_ensemble)
34 , offline_data_(&offline_data)
35 , hyperbolic_system_(&hyperbolic_system)
36 , parabolic_system_(&parabolic_system)
37 , postprocessor_(&postprocessor)
38 , initial_precomputed_(initial_precomputed)
39 , alpha_(alpha)
40 , smoothness_indicators_(smoothness_indicators)
41 {
42 use_mpi_io_ = true;
43 add_parameter("use mpi io",
44 use_mpi_io_,
45 "If enabled write out one vtu file via MPI IO using "
46 "write_vtu_in_parallel() instead of independent output files "
47 "via write_vtu_with_pvtu_record()");
48
49 add_parameter("manifolds",
50 manifolds_,
51 "List of level set functions. The description is used to "
52 "only output cells that intersect the given level set.");
53
54 std::copy(std::begin(View::component_names),
55 std::end(View::component_names),
56 std::back_inserter(vtu_output_quantities_));
57
58 std::copy(std::begin(View::initial_precomputed_names),
59 std::end(View::initial_precomputed_names),
60 std::back_inserter(vtu_output_quantities_));
61
62 add_parameter("vtu output quantities",
63 vtu_output_quantities_,
64 "List of conserved, primitive, precomputed, or postprocessed "
65 "quantities that will be written to the vtu files.");
66 }
67
68
69 template <typename Description, int dim, typename Number>
71 {
72#ifdef DEBUG_OUTPUT
73 std::cout << "VTUOutput<dim, Number>::prepare()" << std::endl;
74#endif
75
77 {"alpha", "smoothness_indicators"}, vtu_output_quantities_);
78 }
79
80
81 template <typename Description, int dim, typename Number>
83 const StateVector &state_vector,
84 std::string name,
85 Number t,
86 unsigned int cycle,
87 bool output_full,
88 bool output_levelsets)
89 {
90#ifdef DEBUG_OUTPUT
91 std::cout << "VTUOutput<dim, Number>::schedule_output()" << std::endl;
92#endif
93 const auto &affine_constraints = offline_data_->affine_constraints();
94
95 /*
96 * Extract quantities and store in ScalarVectors so that we can call
97 * DataOut::add_data_vector()
98 */
99
100 auto selected_components =
102 *offline_data_,
103 *hyperbolic_system_,
104 state_vector,
105 initial_precomputed_,
106 {"alpha", "smoothness_indicators"},
107 {alpha_, smoothness_indicators_},
108 vtu_output_quantities_);
109
110 /* prepare DataOut: */
111
112 auto data_out = std::make_unique<dealii::DataOut<dim>>();
113 data_out->attach_dof_handler(offline_data_->dof_handler());
114
115 for (unsigned int d = 0; d < selected_components.size(); ++d) {
116 affine_constraints.distribute(selected_components[d]);
117 selected_components[d].update_ghost_values();
118 data_out->add_data_vector(selected_components[d],
119 vtu_output_quantities_[d],
120 DataOut<dim>::type_dof_data);
121 }
122
123 const auto n_quantities = postprocessor_->n_quantities();
124 for (unsigned int i = 0; i < n_quantities; ++i)
125 data_out->add_data_vector(postprocessor_->quantities()[i],
126 postprocessor_->component_names()[i],
127 DataOut<dim>::type_dof_data);
128
129 DataOutBase::VtkFlags flags(t,
130 cycle,
131 true,
132#if DEAL_II_VERSION_GTE(9, 5, 0)
133 DataOutBase::CompressionLevel::best_speed);
134#else
135 DataOutBase::VtkFlags::best_speed);
136#endif
137 data_out->set_flags(flags);
138
139 const auto &discretization = offline_data_->discretization();
140 const auto &mapping = discretization.mapping();
141 const auto patch_order =
142 std::max(1u, discretization.polynomial_degree()) - 1u;
143
144 /* Perform output: */
145
146 if (output_full) {
147 data_out->build_patches(mapping, patch_order);
148
149 if (use_mpi_io_) {
150 /* MPI-based synchronous IO */
151 data_out->write_vtu_in_parallel(
152 name + "_" + Utilities::to_string(cycle, 6) + ".vtu",
153 mpi_ensemble_.ensemble_communicator());
154 } else {
155 data_out->write_vtu_with_pvtu_record(
156 "", name, cycle, mpi_ensemble_.ensemble_communicator(), 6);
157 }
158 }
159
160 if (output_levelsets && manifolds_.size() != 0) {
161 /*
162 * Specify an output filter that selects only cells for output that are
163 * in the viscinity of a specified set of output planes:
164 */
165
166 std::vector<std::shared_ptr<FunctionParser<dim>>> level_set_functions;
167 for (const auto &expression : manifolds_)
168 level_set_functions.emplace_back(
169 std::make_shared<FunctionParser<dim>>(expression));
170
171 data_out->set_cell_selection([level_set_functions](const auto &cell) {
172 if (!cell->is_active() || cell->is_artificial())
173 return false;
174
175 for (const auto &function : level_set_functions) {
176
177 unsigned int above = 0;
178 unsigned int below = 0;
179
180 for (unsigned int v : cell->vertex_indices()) {
181 const auto vertex = cell->vertex(v);
182 constexpr auto eps = std::numeric_limits<Number>::epsilon();
183 if (function->value(vertex) >= 0. - 100. * eps)
184 above++;
185 if (function->value(vertex) <= 0. + 100. * eps)
186 below++;
187 if (above > 0 && below > 0)
188 return true;
189 }
190 }
191 return false;
192 });
193
194 data_out->build_patches(mapping, patch_order);
195
196 if (use_mpi_io_) {
197 /* MPI-based synchronous IO */
198 data_out->write_vtu_in_parallel(
199 name + "-levelsets_" + Utilities::to_string(cycle, 6) + ".vtu",
200 mpi_ensemble_.ensemble_communicator());
201 } else {
202 data_out->write_vtu_with_pvtu_record(
203 "",
204 name + "-levelsets",
205 cycle,
206 mpi_ensemble_.ensemble_communicator(),
207 6);
208 }
209 }
210
211 /* Explicitly delete pointer to free up memory early: */
212 data_out.reset();
213 }
214
215} /* namespace ryujin */
void schedule_output(const StateVector &state_vector, std::string name, Number t, unsigned int cycle, bool output_full=true, bool output_cutplanes=true)
VTUOutput(const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem &parabolic_system, const Postprocessor< Description, dim, Number > &postprocessor, const InitialPrecomputedVector &initial_precomputed, const ScalarVector &alpha, const ScalarVector &smoothness_indicators, const std::string &subsection="/VTUOutput")
typename Description::HyperbolicSystem HyperbolicSystem
Definition: vtu_output.h:38
Vectors::ScalarVector< Number > ScalarVector
Definition: vtu_output.h:54
typename Description::ParabolicSystem ParabolicSystem
Definition: vtu_output.h:39
typename View::InitialPrecomputedVector InitialPrecomputedVector
Definition: vtu_output.h:53
static std::vector< ScalarVector > extract(const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const StateVector &state_vector, const InitialPrecomputedVector &initial_precomputed, const std::vector< std::string > &additional_names, const std::vector< std::reference_wrapper< const ScalarVector > > &additional_vectors, const std::vector< std::string > &selected)
static void check(const std::vector< std::string > &additional_names, const std::vector< std::string > &selected)