Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

FEEval::check_vector_compatibility: improve assert message #13420

Merged
merged 1 commit into from
Feb 21, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
8 changes: 6 additions & 2 deletions include/deal.II/matrix_free/fe_evaluation.h
Original file line number Diff line number Diff line change
Expand Up @@ -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
Expand Down
92 changes: 78 additions & 14 deletions include/deal.II/matrix_free/vector_access_internal.h
Original file line number Diff line number Diff line change
Expand Up @@ -22,8 +22,12 @@
#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

Expand Down Expand Up @@ -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());
Expand All @@ -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
}


Expand Down