#include <deal.II/base/subscriptor.h>
#include <deal.II/base/utilities.h>
+#include <deal.II/lac/block_vector_base.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/solver.h>
#include <deal.II/lac/solver_control.h>
}
return *data[i];
}
+
+
+
+ template <typename VectorType,
+ std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
+ * = nullptr>
+ unsigned int
+ n_blocks(const VectorType &)
+ {
+ return 1;
+ }
+
+
+
+ template <typename VectorType,
+ std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
+ nullptr>
+ unsigned int
+ n_blocks(const VectorType &vector)
+ {
+ return vector.n_blocks();
+ }
+
+
+
+ template <typename VectorType,
+ std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
+ * = nullptr>
+ VectorType &
+ block(VectorType &vector, const unsigned int b)
+ {
+ AssertDimension(b, 0);
+ (void)b;
+ return vector;
+ }
+
+
+
+ template <typename VectorType,
+ std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
+ nullptr>
+ typename VectorType::BlockType &
+ block(VectorType &vector, const unsigned int b)
+ {
+ return vector.block(b);
+ }
+
} // namespace SolverIDRImplementation
} // namespace internal
VectorType &tmp_q = Q(i, x);
if (i != 0)
{
- for (auto indx : tmp_q.locally_owned_elements())
- tmp_q(indx) = normal_distribution(rng);
+ for (unsigned int b = 0;
+ b < internal::SolverIDRImplementation::n_blocks(tmp_q);
+ ++b)
+ for (auto indx : internal::SolverIDRImplementation::block(tmp_q, b)
+ .locally_owned_elements())
+ internal::SolverIDRImplementation::block(tmp_q, b)(indx) =
+ normal_distribution(rng);
tmp_q.compress(VectorOperation::insert);
}
else