#include <deal.II/base/vectorization.h>
#include <deal.II/matrix_free/dof_info.h>
+#include <deal.II/matrix_free/matrix_free.h>
#include <deal.II/matrix_free/type_traits.h>
+#if DEBUG
+# include <boost/algorithm/string/join.hpp>
+#endif
DEAL_II_NAMESPACE_OPEN
// version below is when has_partitioners_are_compatible == false
// FIXME: this is incorrect for PETSc/Trilinos MPI vectors
template <
+ int dim,
+ typename Number,
+ typename VectorizedArrayType,
typename VectorType,
typename std::enable_if<!has_partitioners_are_compatible<VectorType>,
VectorType>::type * = nullptr>
inline void
check_vector_compatibility(
- const VectorType & vec,
- const internal::MatrixFreeFunctions::DoFInfo &dof_info)
+ const VectorType & vec,
+ const MatrixFree<dim, Number, VectorizedArrayType> &matrix_free,
+ const internal::MatrixFreeFunctions::DoFInfo & dof_info)
{
(void)vec;
+ (void)matrix_free;
(void)dof_info;
AssertDimension(vec.size(), dof_info.vector_partitioner->size());
// same as above for has_partitioners_are_compatible == true
- template <typename VectorType,
+ template <int dim,
+ typename Number,
+ typename VectorizedArrayType,
+ typename VectorType,
typename std::enable_if<has_partitioners_are_compatible<VectorType>,
VectorType>::type * = nullptr>
inline void
check_vector_compatibility(
- const VectorType & vec,
- const internal::MatrixFreeFunctions::DoFInfo &dof_info)
+ const VectorType & vec,
+ const MatrixFree<dim, Number, VectorizedArrayType> &matrix_free,
+ const internal::MatrixFreeFunctions::DoFInfo & dof_info)
{
(void)vec;
+ (void)matrix_free;
(void)dof_info;
- Assert(vec.partitioners_are_compatible(*dof_info.vector_partitioner),
- ExcMessage(
- "The parallel layout of the given vector is not "
- "compatible with the parallel partitioning in MatrixFree. "
- "A potential reason is that you did not use "
- "MatrixFree::initialize_dof_vector() to get a "
- "compatible vector. If you did, the dof_handler_index "
- "used there and the one passed to the constructor of "
- "FEEvaluation do not match."));
+
+#if DEBUG
+ if (vec.partitioners_are_compatible(*dof_info.vector_partitioner) == false)
+ {
+ unsigned int dof_index = numbers::invalid_unsigned_int;
+
+ for (unsigned int i = 0; i < matrix_free.n_components(); ++i)
+ if (&matrix_free.get_dof_info(i) == &dof_info)
+ {
+ dof_index = i;
+ break;
+ }
+
+ Assert(dof_index != numbers::invalid_unsigned_int, ExcInternalError());
+
+ std::vector<std::string> dof_indices_with_compatible_partitioners;
+
+ for (unsigned int i = 0; i < matrix_free.n_components(); ++i)
+ if (vec.partitioners_are_compatible(
+ *matrix_free.get_dof_info(i).vector_partitioner))
+ dof_indices_with_compatible_partitioners.push_back(
+ std::to_string(i));
+
+ if (dof_indices_with_compatible_partitioners.empty())
+ {
+ Assert(false,
+ ExcMessage("The parallel layout of the given vector is "
+ "compatible neither with the Partitioner of the "
+ "current FEEvaluation with dof_handler_index=" +
+ std::to_string(dof_index) +
+ " nor with any Partitioner in MatrixFree. A "
+ "potential reason is that you did not use "
+ "MatrixFree::initialize_dof_vector() to get a "
+ "compatible vector."));
+ }
+ else
+ {
+ Assert(
+ false,
+ ExcMessage(
+ "The parallel layout of the given vector is "
+ "not compatible with the Partitioner of the "
+ "current FEEvaluation with dof_handler_index=" +
+ std::to_string(dof_index) +
+ ". However, the underlying "
+ "MatrixFree contains Partitioner objects that are compatible. "
+ "They have the following dof_handler_index values: " +
+ boost::algorithm::join(dof_indices_with_compatible_partitioners,
+ ", ") +
+ ". Did you want to pass any of these values to the "
+ "constructor of the current FEEvaluation object or "
+ "did you not use MatrixFree::initialize_dof_vector() "
+ "with dof_handler_index=" +
+ std::to_string(dof_index) +
+ " to get a "
+ "compatible vector?"));
+ }
+ }
+#endif
}