From: Peter Munch Date: Fri, 18 Feb 2022 22:36:51 +0000 (+0100) Subject: FEEval::check_vector_compatibility: improve assert message X-Git-Tag: v9.4.0-rc1~457^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F13420%2Fhead;p=dealii.git FEEval::check_vector_compatibility: improve assert message --- diff --git a/include/deal.II/matrix_free/fe_evaluation.h b/include/deal.II/matrix_free/fe_evaluation.h index 0803f717a2..647edf9208 100644 --- a/include/deal.II/matrix_free/fe_evaluation.h +++ b/include/deal.II/matrix_free/fe_evaluation.h @@ -3279,11 +3279,15 @@ FEEvaluationBase:: "FEEvaluation. In that case, you must pass an " "std::vector or a BlockVector to " + "read_dof_values and distribute_local_to_global.")); - internal::check_vector_compatibility(*src[comp], *this->dof_info); + internal::check_vector_compatibility(*src[comp], + *this->matrix_free, + *this->dof_info); } else { - internal::check_vector_compatibility(*src[0], *this->dof_info); + internal::check_vector_compatibility(*src[0], + *this->matrix_free, + *this->dof_info); } // Case 2: contiguous indices which use reduced storage of indices and can diff --git a/include/deal.II/matrix_free/vector_access_internal.h b/include/deal.II/matrix_free/vector_access_internal.h index 54346511b3..1687a670bd 100644 --- a/include/deal.II/matrix_free/vector_access_internal.h +++ b/include/deal.II/matrix_free/vector_access_internal.h @@ -22,8 +22,12 @@ #include #include +#include #include +#if DEBUG +# include +#endif DEAL_II_NAMESPACE_OPEN @@ -169,15 +173,20 @@ namespace internal // 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, VectorType>::type * = nullptr> inline void check_vector_compatibility( - const VectorType & vec, - const internal::MatrixFreeFunctions::DoFInfo &dof_info) + const VectorType & vec, + const MatrixFree &matrix_free, + const internal::MatrixFreeFunctions::DoFInfo & dof_info) { (void)vec; + (void)matrix_free; (void)dof_info; AssertDimension(vec.size(), dof_info.vector_partitioner->size()); @@ -186,25 +195,80 @@ namespace internal // same as above for has_partitioners_are_compatible == true - template , VectorType>::type * = nullptr> inline void check_vector_compatibility( - const VectorType & vec, - const internal::MatrixFreeFunctions::DoFInfo &dof_info) + const VectorType & vec, + const MatrixFree &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 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 }