12#include <deal.II/base/function_parser.h>
13#include <deal.II/numerics/data_out.h>
14#include <deal.II/numerics/vector_tools.h>
18 template <
typename Description,
int dim,
typename Number>
20 const MPI_Comm &mpi_communicator,
24 const std::string &subsection )
25 : ParameterAcceptor(subsection)
26 , mpi_communicator_(mpi_communicator)
27 , offline_data_(&offline_data)
28 , hyperbolic_system_(&hyperbolic_system)
29 , parabolic_system_(¶bolic_system)
32 add_parameter(
"schlieren beta",
34 "Beta factor used in the exponential scale for the schlieren "
35 "and vorticity plots");
37 recompute_bounds_ =
true;
39 "schlieren recompute bounds",
41 "Recompute bounds for every output cycle. If set to false, bounds once "
42 "at the beginning and reused thereafter.");
44 static_assert(View::component_names.size() > 0,
45 "Need at least one scalar quantitity");
46 schlieren_quantities_.push_back(View::component_names[0]);
49 "schlieren quantities",
50 schlieren_quantities_,
51 "List of conserved quantities used for the schlieren postprocessor.");
53 if constexpr (dim > 1) {
55 "vorticity quantities",
56 vorticity_quantities_,
57 "List of conserved quantities used for the vorticity postprocessor.");
62 template <
typename Description,
int dim,
typename Number>
66 std::cout <<
"Postprocessor<dim, Number>::prepare()" << std::endl;
70 component_names_.clear();
71 schlieren_indices_.clear();
72 vorticity_indices_.clear();
74 const auto populate = [&](
const auto &strings,
77 const auto &cons = View::component_names;
78 const auto &prim = View::primitive_component_names;
79 for (
const auto &entry : strings) {
81 for (
const auto &[is_primitive, names] :
82 {std::make_pair(
false, cons), std::make_pair(
true, prim)}) {
83 const auto pos = std::find(std::begin(names), std::end(names), entry);
84 if (!found && pos != std::end(names)) {
85 const auto index = std::distance(std::begin(names), pos);
86 indices.push_back(std::make_pair(is_primitive, index));
87 component_names_.push_back(pre + entry);
93 dealii::ExcMessage(
"Invalid component name »" + entry +
"«"));
96 populate(schlieren_quantities_, schlieren_indices_,
"schlieren_");
97 populate(vorticity_quantities_, vorticity_indices_,
"vorticity_");
99 const auto &partitioner = offline_data_->scalar_partitioner();
101 quantities_.resize(component_names_.size());
102 for (
auto &it : quantities_)
103 it.reinit(partitioner);
107 template <
typename Description,
int dim,
typename Number>
112 std::cout <<
"Postprocessor<dim, Number>::compute()" << std::endl;
115 const auto &U = std::get<0>(state_vector);
117 using VA = dealii::VectorizedArray<Number>;
119 const auto &affine_constraints = offline_data_->affine_constraints();
121 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
122 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
123 const auto &cij_matrix = offline_data_->cij_matrix();
125 const unsigned int n_internal = offline_data_->n_locally_internal();
126 const unsigned int n_owned = offline_data_->n_locally_owned();
128 const unsigned int n_schlieren = schlieren_indices_.size();
129 Assert(n_schlieren == schlieren_quantities_.size(),
130 dealii::ExcInternalError());
131 const unsigned int n_vorticities = vorticity_indices_.size();
132 Assert(n_vorticities == vorticity_quantities_.size(),
133 dealii::ExcInternalError());
134 const unsigned int n_quantities = n_schlieren + n_vorticities;
135 Assert(n_quantities == quantities_.size(), dealii::ExcInternalError());
136 Assert(n_quantities == component_names_.size(), dealii::ExcInternalError());
145 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
146 using T =
decltype(sentinel);
147 unsigned int stride_size = get_stride_size<T>;
149 std::vector<grad_type<T>> local_schlieren_values(n_schlieren);
150 std::vector<curl_type<T>> local_vorticity_values(n_vorticities);
153 for (
unsigned int i = left; i < right; i += stride_size) {
155 for (
auto &it : local_schlieren_values)
157 for (
auto &it : local_vorticity_values)
161 const unsigned int row_length = sparsity_simd.row_length(i);
165 const unsigned int *js = sparsity_simd.columns(i);
166 for (
unsigned int col_idx = 0; col_idx < row_length;
167 ++col_idx, js += stride_size) {
169 const auto U_j = U.template get_tensor<T>(js);
170 const auto view = hyperbolic_system_->template view<dim, T>();
171 const auto prim_j = view.to_primitive_state(U_j);
173 const auto c_ij = cij_matrix.template get_tensor<T>(i, col_idx);
176 for (
const auto &[is_primitive, index] : schlieren_indices_) {
177 local_schlieren_values[k++] -=
178 c_ij * (is_primitive ? prim_j[index] : U_j[index]);
182 for (
const auto &[is_primitive, index] : vorticity_indices_) {
184 for (
unsigned int d = 0; d < dim; ++d)
185 q_j[d] = (is_primitive ? prim_j[index + d] : U_j[index + d]);
187 if constexpr (dim == 2) {
188 local_vorticity_values[k++][0] -= cross_product_2d(c_ij) * q_j;
189 }
else if constexpr (dim == 3) {
190 local_vorticity_values[k++] -= cross_product_3d(c_ij, q_j);
197 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
201 for (
const auto &schlieren : local_schlieren_values) {
202 const auto value_i = schlieren.norm() / m_i;
203 write_entry<T>(quantities_[k++], value_i, i);
206 for (
const auto &vorticity : local_vorticity_values) {
208 (dim == 2 ? vorticity[0] / m_i : vorticity.norm() / m_i);
209 write_entry<T>(quantities_[k++], value_i, i);
215 loop(Number(), n_internal, n_owned);
217 loop(VA(), 0, n_internal);
227 if (recompute_bounds_)
230 if (bounds_.size() != n_quantities) {
234 std::make_pair(Number(0.), std::numeric_limits<Number>::max()));
236 for (
unsigned int d = 0; d < n_quantities; ++d) {
237 auto &[q_max, q_min] = bounds_[d];
238 for (
unsigned int i = 0; i < n_owned; ++i) {
239 const auto q = quantities_[d].local_element(i);
240 q_max = std::max(q_max, std::abs(q));
241 q_min = std::min(q_min, std::abs(q));
243 q_max = dealii::Utilities::MPI::max(q_max, mpi_communicator_);
244 q_min = dealii::Utilities::MPI::min(q_min, mpi_communicator_);
245 Assert(q_max >= q_min, dealii::ExcInternalError());
254 constexpr Number eps = std::numeric_limits<Number>::epsilon();
255 constexpr Number floor = std::max(Number(1.0e-10), eps);
257 for (
unsigned int d = 0; d < n_quantities; ++d) {
258 auto &[q_max, q_min] = bounds_[d];
259 for (
unsigned int i = 0; i < n_owned; ++i) {
260 auto &q = quantities_[d].local_element(i);
262 const auto ratio = std::max(Number(0.), std::abs(q) - q_min - floor) /
263 std::max(q_max - q_min, eps);
265 const auto magnitude = Number(1.) - std::exp(-beta_ * ratio);
266 q = std::copysign(magnitude, q);
275 for (
auto &it : quantities_) {
276 affine_constraints.distribute(it);
277 it.update_ghost_values();
typename Description::HyperbolicSystem HyperbolicSystem
typename Description::ParabolicSystem ParabolicSystem
dealii::Tensor< 1, dim==2 ? 1 :dim, T > curl_type
dealii::Tensor< 1, dim, T > grad_type
void compute(const StateVector &state_vector) const
Postprocessor(const MPI_Comm &mpi_communicator, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system, const std::string &subsection="/Postprocessor")
typename View::StateVector StateVector
#define RYUJIN_PARALLEL_REGION_BEGIN
#define RYUJIN_PARALLEL_REGION_END