]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a specific exception for MPI failures.
authorDavid Wells <wellsd2@rpi.edu>
Wed, 9 Nov 2016 23:54:09 +0000 (18:54 -0500)
committerDavid Wells <wellsd2@rpi.edu>
Fri, 11 Nov 2016 20:03:51 +0000 (15:03 -0500)
include/deal.II/base/exceptions.h
source/base/exceptions.cc
tests/base/mpi_exceptions.cc [new file with mode: 0644]
tests/base/mpi_exceptions.mpirun=1.output [new file with mode: 0644]
tests/base/mpi_exceptions.mpirun=1.output.mpich [new file with mode: 0644]

index 8649ee67fe6f7819122b7b6b6ca2c612e718220e..ae3e10f9c2574e6fa29e36b3a6e284ab5a63b32d 100644 (file)
@@ -1082,6 +1082,40 @@ namespace StandardExceptions
                   << arg1);
 #endif
 //@}
+
+#ifdef DEAL_II_WITH_MPI
+  /**
+   * Exception for MPI errors. This exception is only defined if
+   * <code>deal.II</code> is compiled with MPI support. This exception should
+   * be used with <code>AssertThrow</code> to check error codes of MPI
+   * functions. For example:
+   * @code
+   * const int ierr = MPI_Isend(...);
+   * AssertThrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
+   * @endcode
+   * or, using the convenience macro <code>AssertThrowMPI</code>,
+   * @code
+   * const int ierr = MPI_Irecv(...);
+   * AssertThrowMPI(ierr);
+   * @endcode
+   *
+   * If the assertion fails then the error code will be used to print a helpful
+   * message to the screen by utilizing the <code>MPI_Error_string</code>
+   * function.
+   *
+   * @ingroup Exceptions
+   * @author David Wells, 2016
+   */
+  class ExcMPI : public dealii::ExceptionBase
+  {
+  public:
+    ExcMPI (const int error_code);
+
+    virtual void print_info (std::ostream &out) const;
+
+    const int error_code;
+  };
+#endif // DEAL_II_WITH_MPI
 } /*namespace StandardExceptions*/
 
 
@@ -1149,6 +1183,20 @@ namespace StandardExceptions
 #define AssertIsFinite(number) Assert(dealii::numbers::is_finite(number), \
                                       dealii::ExcNumberNotFinite(std::complex<double>(number)))
 
+#ifdef DEAL_II_WITH_MPI
+/**
+ * An assertion that checks whether or not an error code returned by an MPI
+ * function is equal to <code>MPI_SUCCESS</code>. If the check fails then an
+ * exception of type ExcMPI is thrown with the given error code as an
+ * argument.
+ *
+ * @ingroup Exceptions
+ * @author David Wells, 2016
+ */
+#define AssertThrowMPI(error_code) AssertThrow(error_code == MPI_SUCCESS, \
+                                               dealii::ExcMPI(error_code))
+#endif // DEAL_II_WITH_MPI
+
 #ifdef DEAL_II_WITH_CUDA
 /**
  * An assertion that checks that the error code produced by calling a CUDA
index 865b2c9e5e408b234fcaa9a3e45068710127baed..7f32809772b703ff75c5448d6f437d5a13b37de8 100644 (file)
 
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+#include <deal.II/base/mpi.h>
+
 #include <string>
 #include <cstring>
 #include <cstdlib>
 #include <iostream>
 #include <sstream>
 
+#ifdef DEAL_II_WITH_MPI
+#  include <mpi.h>
+#endif
+
 #ifdef DEAL_II_HAVE_GLIBC_STACKTRACE
 #  include <execinfo.h>
 #endif
@@ -336,6 +343,58 @@ void ExceptionBase::generate_message () const
 
 
 
+#ifdef DEAL_II_WITH_MPI
+namespace StandardExceptions
+{
+  ExcMPI::ExcMPI (const int error_code)
+    :
+    error_code (error_code)
+  {}
+
+  void ExcMPI::print_info (std::ostream &out) const
+  {
+    char error_name[MPI_MAX_ERROR_STRING];
+    error_name[0] = '\0';
+    int resulting_length = MPI_MAX_ERROR_STRING;
+
+    // get the string name of the error code by first converting it to an error
+    // class.
+    int error_class = 0;
+    int ierr = MPI_Error_class (error_code, &error_class);
+
+    // Check the output of the error printing functions. If either MPI function
+    // fails we should just print a less descriptive message.
+    bool error_name_known = ierr == MPI_SUCCESS;
+    if (error_name_known)
+      {
+        ierr = MPI_Error_string (error_class, error_name, &resulting_length);
+        error_name_known = error_name_known && (ierr == MPI_SUCCESS);
+      }
+
+    out << "deal.II encountered an error while calling an MPI function."
+        << std::endl;
+    if (error_name_known)
+      {
+        out << "The description of the error provided by MPI is \""
+            << error_name
+            << "\"."
+            << std::endl;
+      }
+    else
+      {
+        out << "This error code is not equal to any of the standard MPI error codes."
+            << std::endl;
+      }
+    out << "The numerical value of the original error code is "
+        << error_code
+        << "."
+        << std::endl;
+  }
+}
+#endif // DEAL_II_WITH_MPI
+
+
+
 namespace deal_II_exceptions
 {
   namespace internals
diff --git a/tests/base/mpi_exceptions.cc b/tests/base/mpi_exceptions.cc
new file mode 100644 (file)
index 0000000..61c1895
--- /dev/null
@@ -0,0 +1,68 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check that we can catch all MPI-1 error codes and print appropriate
+// exception messages to the screen.
+
+#include "../tests.h"
+#include <deal.II/base/exceptions.h>
+#include <fstream>
+
+#include <mpi.h>
+
+#include <vector>
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  std::vector<int> mpi_error_codes;
+  mpi_error_codes.push_back(MPI_ERR_BUFFER);
+  mpi_error_codes.push_back(MPI_ERR_COUNT);
+  mpi_error_codes.push_back(MPI_ERR_TYPE);
+  mpi_error_codes.push_back(MPI_ERR_TAG);
+  mpi_error_codes.push_back(MPI_ERR_COMM);
+  mpi_error_codes.push_back(MPI_ERR_RANK);
+  mpi_error_codes.push_back(MPI_ERR_REQUEST);
+  mpi_error_codes.push_back(MPI_ERR_ROOT);
+  mpi_error_codes.push_back(MPI_ERR_GROUP);
+  mpi_error_codes.push_back(MPI_ERR_OP);
+  mpi_error_codes.push_back(MPI_ERR_TOPOLOGY);
+  mpi_error_codes.push_back(MPI_ERR_UNKNOWN);
+  mpi_error_codes.push_back(MPI_ERR_OTHER);
+  mpi_error_codes.push_back(MPI_ERR_INTERN);
+  mpi_error_codes.push_back(MPI_ERR_LASTCODE);
+
+  for (unsigned int i = 0; i < mpi_error_codes.size(); ++i)
+    {
+      try
+        {
+          throw ExcMPI(mpi_error_codes[i]);
+        }
+      catch (const ExceptionBase &exc)
+        {
+          deallog << exc.get_exc_name() << std::endl;
+          exc.print_info(deallog.get_file_stream());
+        }
+    }
+
+
+  return 0;
+}
diff --git a/tests/base/mpi_exceptions.mpirun=1.output b/tests/base/mpi_exceptions.mpirun=1.output
new file mode 100644 (file)
index 0000000..dcd5592
--- /dev/null
@@ -0,0 +1,61 @@
+
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_BUFFER: invalid buffer pointer".
+The numerical value of the original error code is 1.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_COUNT: invalid count argument".
+The numerical value of the original error code is 2.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_TYPE: invalid datatype".
+The numerical value of the original error code is 3.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_TAG: invalid tag".
+The numerical value of the original error code is 4.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_COMM: invalid communicator".
+The numerical value of the original error code is 5.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_RANK: invalid rank".
+The numerical value of the original error code is 6.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_REQUEST: invalid request".
+The numerical value of the original error code is 7.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_ROOT: invalid root".
+The numerical value of the original error code is 8.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_GROUP: invalid group".
+The numerical value of the original error code is 9.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_OP: invalid reduce operation".
+The numerical value of the original error code is 10.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_TOPOLOGY: invalid communicator topology".
+The numerical value of the original error code is 11.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_UNKNOWN: unknown error".
+The numerical value of the original error code is 14.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_OTHER: known error not in list".
+The numerical value of the original error code is 16.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_INTERN: internal error".
+The numerical value of the original error code is 17.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "MPI_ERR_UNKNOWN: unknown error".
+The numerical value of the original error code is 92.
diff --git a/tests/base/mpi_exceptions.mpirun=1.output.mpich b/tests/base/mpi_exceptions.mpirun=1.output.mpich
new file mode 100644 (file)
index 0000000..8da187c
--- /dev/null
@@ -0,0 +1,61 @@
+
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid buffer pointer".
+The numerical value of the original error code was 1.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid count".
+The numerical value of the original error code was 2.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid datatype".
+The numerical value of the original error code was 3.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid tag".
+The numerical value of the original error code was 4.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid communicator".
+The numerical value of the original error code was 5.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid rank".
+The numerical value of the original error code was 6.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid MPI_Request".
+The numerical value of the original error code was 19.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid root".
+The numerical value of the original error code was 7.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid group".
+The numerical value of the original error code was 8.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid MPI_Op".
+The numerical value of the original error code was 9.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Invalid topology".
+The numerical value of the original error code was 10.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Unknown error.  Please file a bug report.".
+The numerical value of the original error code was 13.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Other MPI error".
+The numerical value of the original error code was 15.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Internal MPI error!".
+The numerical value of the original error code was 16.
+DEAL::
+deal.II encountered an error while calling an MPI function.
+The description of the error provided by MPI is "Unknown error class".
+The numerical value of the original error code was 1073741823.

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.