]> https://gitweb.dealii.org/ - dealii.git/commitdiff
FEEval::check_vector_compatibility: improve assert message 13420/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 18 Feb 2022 22:36:51 +0000 (23:36 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 21 Feb 2022 08:58:14 +0000 (09:58 +0100)
include/deal.II/matrix_free/fe_evaluation.h
include/deal.II/matrix_free/vector_access_internal.h

index 0803f717a2f9dbc07ea4de9fbc0c7ff8becad7cb..647edf920880123aab6362d3c79cb48045172ce8 100644 (file)
@@ -3279,11 +3279,15 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
                           "FEEvaluation. In that case, you must pass an "
                           "std::vector<VectorType> 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
index 54346511b39190ac6531f9339a8ed1fc8af69be4..1687a670bd7f3775bfc65287025bb2e10fd2f3bc 100644 (file)
 #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
 
@@ -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<!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());
@@ -186,25 +195,80 @@ namespace internal
 
 
   // 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
   }
 
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.