11#include <deal.II/base/function_parser.h>
12#include <deal.II/numerics/data_out.h>
13#include <deal.II/numerics/vector_tools.h>
18 using namespace dealii;
21 template <
typename Description,
int dim,
typename Number>
31 const std::string &subsection )
32 : ParameterAcceptor(subsection)
33 , mpi_ensemble_(mpi_ensemble)
34 , offline_data_(&offline_data)
35 , hyperbolic_system_(&hyperbolic_system)
36 , parabolic_system_(¶bolic_system)
37 , postprocessor_(&postprocessor)
38 , initial_precomputed_(initial_precomputed)
43 add_parameter(
"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()");
49 add_parameter(
"manifolds",
51 "List of level set functions. The description is used to "
52 "only output cells that intersect the given level set.");
54 std::copy(std::begin(View::component_names),
55 std::end(View::component_names),
56 std::back_inserter(vtu_output_quantities_));
58 std::copy(std::begin(View::initial_precomputed_names),
59 std::end(View::initial_precomputed_names),
60 std::back_inserter(vtu_output_quantities_));
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.");
69 template <
typename Description,
int dim,
typename Number>
73 std::cout <<
"VTUOutput<dim, Number>::prepare()" << std::endl;
77 {
"alpha",
"smoothness_indicators"}, vtu_output_quantities_);
81 template <
typename Description,
int dim,
typename Number>
88 bool output_levelsets)
91 std::cout <<
"VTUOutput<dim, Number>::schedule_output()" << std::endl;
93 const auto &affine_constraints = offline_data_->affine_constraints();
100 auto selected_components =
105 initial_precomputed_,
106 {
"alpha",
"smoothness_indicators"},
107 {alpha_, smoothness_indicators_},
108 vtu_output_quantities_);
112 auto data_out = std::make_unique<dealii::DataOut<dim>>();
113 data_out->attach_dof_handler(offline_data_->dof_handler());
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);
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);
129 DataOutBase::VtkFlags flags(t,
132#
if DEAL_II_VERSION_GTE(9, 5, 0)
133 DataOutBase::CompressionLevel::best_speed);
135 DataOutBase::VtkFlags::best_speed);
137 data_out->set_flags(flags);
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;
147 data_out->build_patches(mapping, patch_order);
151 data_out->write_vtu_in_parallel(
152 name +
"_" + Utilities::to_string(cycle, 6) +
".vtu",
153 mpi_ensemble_.ensemble_communicator());
155 data_out->write_vtu_with_pvtu_record(
156 "", name, cycle, mpi_ensemble_.ensemble_communicator(), 6);
160 if (output_levelsets && manifolds_.size() != 0) {
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));
171 data_out->set_cell_selection([level_set_functions](
const auto &cell) {
172 if (!cell->is_active() || cell->is_artificial())
175 for (
const auto &function : level_set_functions) {
177 unsigned int above = 0;
178 unsigned int below = 0;
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)
185 if (function->value(vertex) <= 0. + 100. * eps)
187 if (above > 0 && below > 0)
194 data_out->build_patches(mapping, patch_order);
198 data_out->write_vtu_in_parallel(
199 name +
"-levelsets_" + Utilities::to_string(cycle, 6) +
".vtu",
200 mpi_ensemble_.ensemble_communicator());
202 data_out->write_vtu_with_pvtu_record(
206 mpi_ensemble_.ensemble_communicator(),
typename View::StateVector StateVector
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 ¶bolic_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
Vectors::ScalarVector< Number > ScalarVector
typename Description::ParabolicSystem ParabolicSystem
typename View::InitialPrecomputedVector InitialPrecomputedVector