10#include <deal.II/lac/la_parallel_block_vector.h>
16 template <
int dim,
typename Number>
30 template <
typename Number>
31 using ScalarVector = dealii::LinearAlgebra::distributed::Vector<Number>;
36 template <
typename Number>
37 using BlockVector = dealii::LinearAlgebra::distributed::BlockVector<Number>;
47 template <
typename Number,
unsigned int problem_dim,
unsigned int prec_dim>
60 int problem_dimension = View::problem_dimension,
61 int prec_dimension = View::n_precomputed_values>
68 auto &[U, precomputed, V] = state_vector;
70 const unsigned int n_owned = offline_data.n_locally_owned();
71 const auto &partitioner = offline_data.scalar_partitioner();
73 for (
unsigned int i = 0; i < n_owned; ++i) {
74 if (!offline_data.affine_constraints().is_constrained(
75 partitioner->local_to_global(i)))
77 constexpr auto nan = std::numeric_limits<Number>::signaling_NaN();
78 U.write_tensor(dealii::Tensor<1, problem_dimension, Number>() * nan, i);
90 int problem_dimension = View::problem_dimension,
91 int prec_dimension = View::n_precomputed_values>
98 auto &[U, precomputed, V] = state_vector;
100 constexpr auto nan = std::numeric_limits<Number>::signaling_NaN();
101 const unsigned int n_owned = offline_data.n_locally_owned();
102 const auto block_size = offline_data.n_parabolic_state_vectors();
104 for (
unsigned int i = 0; i < n_owned; ++i) {
105 precomputed.write_tensor(
106 dealii::Tensor<1, prec_dimension, Number>() * nan, i);
107 for (
unsigned int b = 0; b < block_size; ++b) {
108 V.block(b).local_element(i) = nan;
125 int problem_dimension = View::problem_dimension,
126 int prec_dimension = View::n_precomputed_values>
131 auto &[U, precomputed, V] = state_vector;
136 V.reinit(block_size);
137 for (
unsigned int i = 0; i < block_size; ++i) {
143 using state_type =
typename View::state_type;
145 constexpr auto nan = std::numeric_limits<Number>::signaling_NaN();
148 for (
unsigned int i = 0; i < n_owned; ++i) {
149 U.write_tensor(state_type{} * nan, i);
150 precomputed.write_tensor(
151 dealii::Tensor<1, prec_dimension, Number>() * nan, i);
152 for (
unsigned int b = 0; b < block_size; ++b) {
153 V.block(b).local_element(i) = nan;
auto & n_locally_owned() const
const auto & scalar_partitioner() const
const auto & hyperbolic_vector_partitioner() const
const auto & precomputed_vector_partitioner() const
auto & n_parabolic_state_vectors() const
void reinit_state_vector(StateVector< Number, problem_dimension, prec_dimension > &state_vector, const OfflineData< dim, Number > &offline_data)
dealii::LinearAlgebra::distributed::Vector< Number > ScalarVector
void debug_poison_constrained_dofs(StateVector< Number, problem_dimension, prec_dimension > &state_vector, const OfflineData< dim, Number > &offline_data)
dealii::LinearAlgebra::distributed::BlockVector< Number > BlockVector
void debug_poison_precomputed_values(StateVector< Number, problem_dimension, prec_dimension > &state_vector, const OfflineData< dim, Number > &offline_data)
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Euler::Description Description