ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
selected_components_extractor.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
10#include "state_vector.h"
11
12namespace ryujin
13{
14 template <typename Description, int dim, typename Number>
17
18 using View =
19 typename Description::template HyperbolicSystemView<dim, Number>;
20
22 using InitialPrecomputedVector = typename View::InitialPrecomputedVector;
24
26
27 static void check(const std::vector<std::string> &additional_names,
28 const std::vector<std::string> &selected)
29 {
30 const auto search = [&](const auto entry, const auto &names) {
31 const auto pos = std::find(std::begin(names), std::end(names), entry);
32 return pos != std::end(names);
33 };
34
35 for (const auto &entry : selected) {
36 const auto found = search(entry, View::component_names) ||
37 search(entry, View::primitive_component_names) ||
38 search(entry, View::precomputed_names) ||
39 search(entry, View::initial_precomputed_names) ||
40 search(entry, additional_names);
41 AssertThrow(
42 found,
43 dealii::ExcMessage(
44 "Invalid component name: \"" + entry +
45 "\" is not a valid conserved/primitive/precomputed/initial "
46 "component name."));
47 }
48 }
49
50 static std::vector<ScalarVector>
51 extract(const OfflineData<dim, Number> &offline_data,
52 const HyperbolicSystem &hyperbolic_system,
53 const StateVector &state_vector,
54 const InitialPrecomputedVector &initial_precomputed,
55 const std::vector<std::string> &additional_names,
56 const std::vector<std::reference_wrapper<const ScalarVector>>
57 &additional_vectors,
58 const std::vector<std::string> &selected)
59 {
60 /*
61 * Match the selected_components strings against conserved,
62 * primitive, and initial component names and record an index pair
63 * matching return vector position and component index:
64 */
65
66 std::vector<std::tuple<std::size_t, std::size_t>> conserved_indices;
67 std::vector<std::tuple<std::size_t, std::size_t>> primitive_indices;
68 std::vector<std::tuple<std::size_t, std::size_t>> precomputed_indices;
69 std::vector<std::tuple<std::size_t, std::size_t>> initial_indices;
70 std::vector<std::tuple<std::size_t, std::size_t>> additional_indices;
71
72 for (std::size_t i = 0; const auto &entry : selected) {
73 const auto search = [&](const auto &names, auto &indices) {
74 const auto pos = std::find(std::begin(names), std::end(names), entry);
75 if (pos != std::end(names)) {
76 const auto index = std::distance(std::begin(names), pos);
77 indices.push_back({i++, index});
78 return true;
79 }
80 return false;
81 };
82
83 if (search(View::component_names, conserved_indices))
84 ;
85 else if (search(View::primitive_component_names, primitive_indices))
86 ;
87 else if (search(View::precomputed_names, precomputed_indices))
88 ;
89 else if (search(View::initial_precomputed_names, initial_indices))
90 ;
91 else if (search(additional_names, additional_indices))
92 ;
93 else
94 AssertThrow(false, dealii::ExcInternalError());
95 }
96
97 std::vector<ScalarVector> extracted_components(selected.size());
98 const auto &scalar_partitioner = offline_data.scalar_partitioner();
99 for (auto &it : extracted_components)
100 it.reinit(scalar_partitioner);
101
102 for (const auto &[i, k] : conserved_indices) {
103 const auto &U = std::get<0>(state_vector);
104 U.extract_component(extracted_components[i], k);
105 }
106
107 if (!primitive_indices.empty()) {
108 const auto &U = std::get<0>(state_vector);
109 const unsigned int n_owned = scalar_partitioner->locally_owned_size();
110 const auto view = hyperbolic_system.template view<dim, Number>();
111 for (unsigned int i = 0; i < n_owned; ++i) {
112 const auto U_i = U.get_tensor(i);
113 const auto PU_i = view.to_primitive_state(U_i);
114 for (const auto &[j, k] : primitive_indices)
115 extracted_components[j].local_element(i) = PU_i[k];
116 }
117 }
118
119 for (const auto &[i, k] : precomputed_indices) {
120 const auto &prec = std::get<1>(state_vector);
121 prec.extract_component(extracted_components[i], k);
122 }
123
124 for (const auto &[i, k] : initial_indices) {
125 initial_precomputed.extract_component(extracted_components[i], k);
126 }
127
128 for (const auto &[i, k] : additional_indices) {
129 extracted_components[i] = additional_vectors[k];
130 }
131
132 return extracted_components;
133 }
134 };
135} // namespace ryujin
const auto & scalar_partitioner() const
Definition: offline_data.h:130
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
Definition: state_vector.h:33
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Definition: state_vector.h:53
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:34
typename Description::template HyperbolicSystemView< dim, Number > View
typename Description::HyperbolicSystem HyperbolicSystem
Vectors::ScalarVector< Number > ScalarVector
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)
typename View::InitialPrecomputedVector InitialPrecomputedVector
static void check(const std::vector< std::string > &additional_names, const std::vector< std::string > &selected)