Refactor gmx_fatal()
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 13 Apr 2014 05:29:48 +0000 (08:29 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 16 Apr 2014 01:45:01 +0000 (03:45 +0200)
- Remove t_commrec and legacyheaders/ dependencies from
  utility/fatalerror.cpp by splitting gmx_fatal_collective() into a
  low-level helper (in fatalerror.cpp) and a high-level routine that
  knows t_commrec (in network.c).
- Move the declaration of gmx_fatal_collective to network.h to get rid
  of one more header in legacyheaders/.
- Refactor fatalerror.cpp such that the error handler is no longer
  responsible of terminating the program.  Split functionality into
  helper functions that can be called as appropriate, and remove a lot
  of duplication.
- Clean up gmx_abort() and code that calls it.
- Remove interactive prompt from the error handler exit path.

Change-Id: Ia0ab2e415a18bf6ed17bebd1fd8870946d16793e

src/gromacs/gmxlib/gmx_detect_hardware.c
src/gromacs/gmxlib/network.c
src/gromacs/legacyheaders/gmx_fatal_collective.h [deleted file]
src/gromacs/legacyheaders/network.h
src/gromacs/mdlib/domdec.c
src/gromacs/utility/basenetwork.cpp
src/gromacs/utility/basenetwork.h
src/gromacs/utility/exceptions.cpp
src/gromacs/utility/fatalerror.cpp
src/gromacs/utility/fatalerror.h
src/programs/mdrun/runner.c

index 09974276ea5641ecba35008d21b7782a8bf48899..62e67b8c58572bf46d3e5cd112f759229fc08c23 100644 (file)
@@ -48,7 +48,7 @@
 #include "types/enums.h"
 #include "types/hw_info.h"
 #include "types/commrec.h"
-#include "gmx_fatal_collective.h"
+#include "network.h"
 #include "md_logging.h"
 #include "gmx_cpuid.h"
 #include "gpu_utils.h"
index f1186ca8bf6d372803606459b5a495a88749139d..e46c7eea997b5137ca9b0f92bbcf57fac281d950 100644 (file)
  * To help us fund GROMACS development, we humbly ask that you cite
  * the research papers on the package. Check out http://www.gromacs.org.
  */
+#include "network.h"
+
 #ifdef HAVE_CONFIG_H
 #include <config.h>
 #endif
 
 #include <ctype.h>
+#include <stdarg.h>
 #include <string.h>
 
 #include "types/commrec.h"
-#include "network.h"
 #include "copyrite.h"
 #include "macros.h"
 
@@ -711,3 +713,32 @@ gmx_bool gmx_fexist_master(const char *fname, t_commrec *cr)
     }
     return bExist;
 }
+
+void gmx_fatal_collective(int f_errno, const char *file, int line,
+                          const t_commrec *cr, gmx_domdec_t *dd,
+                          const char *fmt, ...)
+{
+    va_list  ap;
+    gmx_bool bMaster, bFinalize;
+#ifdef GMX_MPI
+    int      result;
+    /* Check if we are calling on all processes in MPI_COMM_WORLD */
+    if (cr != NULL)
+    {
+        MPI_Comm_compare(cr->mpi_comm_mysim, MPI_COMM_WORLD, &result);
+    }
+    else
+    {
+        MPI_Comm_compare(dd->mpi_comm_all, MPI_COMM_WORLD, &result);
+    }
+    /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
+    bFinalize = (result != MPI_UNEQUAL);
+#else
+    bFinalize = TRUE;
+#endif
+    bMaster = (cr != NULL && MASTER(cr)) || (dd != NULL && DDMASTER(dd));
+
+    va_start(ap, fmt);
+    gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
+    va_end(ap);
+}
diff --git a/src/gromacs/legacyheaders/gmx_fatal_collective.h b/src/gromacs/legacyheaders/gmx_fatal_collective.h
deleted file mode 100644 (file)
index c7e5957..0000000
+++ /dev/null
@@ -1,67 +0,0 @@
-/*
- * This file is part of the GROMACS molecular simulation package.
- *
- * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
- * Copyright (c) 2001-2012, The GROMACS development team.
- * Copyright (c) 2012,2014, by the GROMACS development team, led by
- * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
- * and including many others, as listed in the AUTHORS file in the
- * top-level source directory and at http://www.gromacs.org.
- *
- * GROMACS is free software; you can 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.
- *
- * GROMACS is distributed in the hope that it will be useful,
- * but WITHOUT ANY WARRANTY; without even the implied warranty of
- * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
- * Lesser General Public License for more details.
- *
- * You should have received a copy of the GNU Lesser General Public
- * License along with GROMACS; if not, see
- * http://www.gnu.org/licenses, or write to the Free Software Foundation,
- * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
- *
- * If you want to redistribute modifications to GROMACS, please
- * consider that scientific software is very special. Version
- * control is crucial - bugs must be traceable. We will be happy to
- * consider code for inclusion in the official distribution, but
- * derived work must not be called official GROMACS. Details are found
- * in the README & COPYING files - if they are missing, get the
- * official version at http://www.gromacs.org.
- *
- * To help us fund GROMACS development, we humbly ask that you cite
- * the research papers on the package. Check out http://www.gromacs.org.
- */
-
-#ifndef _fatal_collective_h
-#define _fatal_collective_h
-
-#include "typedefs.h"
-
-#ifdef __cplusplus
-extern "C" {
-#endif
-
-
-void
-gmx_fatal_collective(int f_errno, const char *file, int line,
-                     const t_commrec *cr, gmx_domdec_t *dd,
-                     const char *fmt, ...);
-/* As gmx_fatal declared in gmx_fatal.h,
- * but only the master process prints the error message.
- * This should only be called one of the following two situations:
- * 1) On all nodes in cr->mpi_comm_mysim, with cr!=NULL,dd==NULL.
- * 2) On all nodes in dd->mpi_comm_all,   with cr==NULL,dd!=NULL.
- * This will call MPI_Finalize instead of MPI_Abort when possible,
- * This is useful for handling errors in code that is executed identically
- * for all processes.
- */
-
-
-#ifdef __cplusplus
-}
-#endif
-
-#endif  /* _fatal_collective_h */
index 90ae14c955470cdb466d9c358987b0b942062d42..42fe9ca74fdad8fd2f8e8ed5d30f63d06f6cabf1 100644 (file)
@@ -116,6 +116,20 @@ void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms);
 gmx_bool gmx_fexist_master(const char *fname, t_commrec *cr);
 /* Return TRUE when fname exists, FALSE otherwise, bcast from master to others */
 
+void
+gmx_fatal_collective(int f_errno, const char *file, int line,
+                     const t_commrec *cr, gmx_domdec_t *dd,
+                     const char *fmt, ...);
+/* As gmx_fatal declared in utility/fatalerror.h,
+ * but only the master process prints the error message.
+ * This should only be called one of the following two situations:
+ * 1) On all nodes in cr->mpi_comm_mysim, with cr!=NULL,dd==NULL.
+ * 2) On all nodes in dd->mpi_comm_all,   with cr==NULL,dd!=NULL.
+ * This will call MPI_Finalize instead of MPI_Abort when possible,
+ * This is useful for handling errors in code that is executed identically
+ * for all processes.
+ */
+
 /* This doesn't currently work if enabled (needs some header cleanup). */
 #ifdef DEBUG_GMX
 #define debug_gmx() do { FILE *fp = debug ? debug : stderr; \
index 98d1f541019ef940aec5a00660e8cf3fa2b54a83..0a10d70cefd5fc1c5ba431c7ba958c701f6f825e 100644 (file)
@@ -47,7 +47,7 @@
 #include "typedefs.h"
 #include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/fatalerror.h"
-#include "gmx_fatal_collective.h"
+#include "network.h"
 #include "vec.h"
 #include "domdec.h"
 #include "domdec_network.h"
index 00aa658db062badfcffefcb4306a3fdf7c715a34..d4adc345a8dfb72cca9202d837921d272b2d4c80 100644 (file)
@@ -242,7 +242,7 @@ int gmx_hostname_num()
 }
 
 #ifdef GMX_LIB_MPI
-void gmx_abort(int noderank, int nnodes, int errorno)
+void gmx_abort(int errorno)
 {
     const char *programName = "GROMACS";
     try
@@ -252,9 +252,11 @@ void gmx_abort(int noderank, int nnodes, int errorno)
     catch (const std::exception &)
     {
     }
+    const int nnodes   = gmx_node_num();
+    const int noderank = gmx_node_rank();
     if (nnodes > 1)
     {
-        std::fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
+        std::fprintf(stderr, "Halting parallel program %s on node %d out of %d\n",
                      programName, noderank, nnodes);
     }
     else
@@ -263,6 +265,6 @@ void gmx_abort(int noderank, int nnodes, int errorno)
     }
 
     MPI_Abort(MPI_COMM_WORLD, errorno);
-    std::exit(1);
+    std::exit(errorno);
 }
 #endif
index 78aa66a782595194df4506b1f8c3fed67df64253..604f90a3a9d447494c696ab21b20a9261de89783 100644 (file)
@@ -101,7 +101,7 @@ int gmx_physicalnode_id_hash(void);
 int gmx_hostname_num(void);
 
 /** Abort the parallel run */
-void gmx_abort(int nodeid, int nnodes, int errorno);
+void gmx_abort(int errorno);
 
 #ifdef __cplusplus
 }
index 545d1547154d919e912ac575a4573a036442a91d..6dfc4330d83324b7c35ffc3d8d32beadb9d98f71 100644 (file)
@@ -516,7 +516,7 @@ int processExceptionAtExit(const std::exception & /*ex*/)
 #ifdef GMX_LIB_MPI
     // TODO: Consider moving the output done in gmx_abort() into the message
     // printing routine above, so that this could become a simple MPI_Abort().
-    gmx_abort(gmx_node_rank(), gmx_node_num(), returnCode);
+    gmx_abort(returnCode);
 #endif
     return returnCode;
 }
index f7544b3698b7582c9a2a7e6b62d265a20c5a8b16..7ca511886b10c5ad6e38c39e3d42fffe173c1f42 100644 (file)
  * the research papers on the package. Check out http://www.gromacs.org.
  */
 #include "fatalerror.h"
-#include "gromacs/legacyheaders/gmx_fatal_collective.h"
 
 #ifdef HAVE_CONFIG_H
 #include <config.h>
 #endif
 
-#include <ctype.h>
-#include <errno.h>
-#include <stdarg.h>
-#include <stdlib.h>
-#include <string.h>
+#include <cerrno>
+#include <cstdarg>
+#include <cstdlib>
+#include <cstring>
 
 #include <exception>
 
 #include "thread_mpi/threads.h"
 
-#include "gromacs/legacyheaders/types/commrec.h"
-
 #include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/baseversion.h"
+#include "gromacs/utility/common.h"
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/futil.h"
 #include "gromacs/utility/gmxmpi.h"
 static gmx_bool            bDebug         = FALSE;
 static FILE               *log_file       = NULL;
 
-static tMPI_Thread_mutex_t debug_mutex     = TMPI_THREAD_MUTEX_INITIALIZER;
-static tMPI_Thread_mutex_t where_mutex     = TMPI_THREAD_MUTEX_INITIALIZER;
+static tMPI_Thread_mutex_t debug_mutex    = TMPI_THREAD_MUTEX_INITIALIZER;
+static tMPI_Thread_mutex_t where_mutex    = TMPI_THREAD_MUTEX_INITIALIZER;
+
+static const char *const   gmxuser
+    = "Please report this to the mailing list (gmx-users@gromacs.org)";
 
 gmx_bool bDebugMode(void)
 {
@@ -121,7 +121,7 @@ void _where(const char *file, int line)
 
 static int fatal_errno = 0;
 
-static void quit_gmx(const char *msg)
+static void default_error_handler(const char *msg)
 {
     tMPI_Thread_mutex_lock(&debug_mutex);
     if (fatal_errno == 0)
@@ -143,174 +143,127 @@ static void quit_gmx(const char *msg)
         }
         perror(msg);
     }
+    tMPI_Thread_mutex_unlock(&debug_mutex);
+}
 
-#ifdef GMX_LIB_MPI
-    if (gmx_mpi_initialized())
-    {
-        int  nnodes;
-        int  noderank;
+static void (*gmx_error_handler)(const char *msg) = default_error_handler;
 
-        nnodes   = gmx_node_num();
-        noderank = gmx_node_rank();
+void set_gmx_error_handler(void (*func)(const char *msg))
+{
+    // TODO: Either this is unnecessary, or also reads to the handler should be
+    // protected by a mutex.
+    tMPI_Thread_mutex_lock(&debug_mutex);
+    gmx_error_handler = func;
+    tMPI_Thread_mutex_unlock(&debug_mutex);
+}
 
-        if (nnodes > 1)
-        {
-            fprintf(stderr, "Error on node %d, will try to stop all the nodes\n",
-                    noderank);
-        }
-        gmx_abort(noderank, nnodes, -1);
-    }
-#endif
+static void call_error_handler(const char *key, const char *file, int line, const char *msg)
+{
+    char        buf[10240], errerrbuf[1024];
+    const char *llines = "-------------------------------------------------------";
+    char       *strerr;
 
-    if (debug)
+    if (msg == NULL)
     {
-        fflush(debug);
+        sprintf(errerrbuf, "Empty fatal_error message. %s", gmxuser);
     }
-    if (bDebugMode())
+    // In case ProgramInfo is not initialized and there is an issue with the
+    // initialization, fall back to "GROMACS".
+    const char *programName = "GROMACS";
+    try
+    {
+        programName = gmx::getProgramContext().displayName();
+    }
+    catch (const std::exception &)
     {
-        fprintf(stderr, "dump core (y/n):");
-        fflush(stderr);
-        if (toupper(getc(stdin)) != 'N')
-        {
-            (void) abort();
-        }
     }
 
-    exit(fatal_errno);
-    tMPI_Thread_mutex_unlock(&debug_mutex);
+    strerr = gmx_strerror(key);
+    sprintf(buf, "\n%s\nProgram %s, %s\n"
+            "Source code file: %s, line: %d\n\n"
+            "%s:\n%s\nFor more information and tips for troubleshooting, please check the GROMACS\n"
+            "website at http://www.gromacs.org/Documentation/Errors\n%s\n",
+            llines, programName, gmx_version(), file, line,
+            strerr, msg ? msg : errerrbuf, llines);
+    free(strerr);
+
+    gmx_error_handler(buf);
 }
 
-/* The function below should be identical to quit_gmx,
- * except that is does not actually quit and call gmx_abort.
- */
-static void quit_gmx_noquit(const char *msg)
+static void do_exit(bool bMaster, bool bFinalize)
 {
-    tMPI_Thread_mutex_lock(&debug_mutex);
-    if (!fatal_errno)
-    {
-        if (log_file)
-        {
-            fprintf(log_file, "%s\n", msg);
-        }
-        fprintf(stderr, "%s\n", msg);
-        /* we set it to no-zero because if this function is called, something
-           has gone wrong */
-        fatal_errno = 255;
-    }
-    else
-    {
-        if (fatal_errno != -1)
-        {
-            errno = fatal_errno;
-        }
-        perror(msg);
-    }
-
-#ifndef GMX_LIB_MPI
     if (debug)
     {
         fflush(debug);
     }
-    if (bDebugMode())
+
+#ifdef GMX_MPI
+    if (gmx_mpi_initialized())
     {
-        fprintf(stderr, "dump core (y/n):");
-        fflush(stderr);
-        if (toupper(getc(stdin)) != 'N')
+        if (bFinalize)
         {
-            (void) abort();
+            /* Broadcast the fatal error number possibly modified
+             * on the master process, in case the user would like
+             * to use the return status on a non-master process.
+             * The master process in cr and dd always has global rank 0.
+             */
+            MPI_Bcast(&fatal_errno, sizeof(fatal_errno), MPI_BYTE,
+                      0, MPI_COMM_WORLD);
+
+            /* Finalize nicely instead of aborting */
+            MPI_Finalize();
         }
-    }
+        else if (bMaster)
+        {
+#ifdef GMX_LIB_MPI
+            gmx_abort(1);
 #endif
-
-    tMPI_Thread_mutex_unlock(&debug_mutex);
-}
-
-void gmx_fatal(int f_errno, const char *file, int line, const char *fmt, ...)
-{
-    va_list ap;
-    char    msg[STRLEN];
-
-    va_start(ap, fmt);
-    vsprintf(msg, fmt, ap);
-    va_end(ap);
-
-    tMPI_Thread_mutex_lock(&debug_mutex);
-    fatal_errno = f_errno;
-    tMPI_Thread_mutex_unlock(&debug_mutex);
-
-    _gmx_error("fatal", msg, file, line);
-}
-
-void gmx_fatal_collective(int f_errno, const char *file, int line,
-                          const t_commrec *cr, gmx_domdec_t *dd,
-                          const char *fmt, ...)
-{
-    va_list     ap;
-    char        msg[STRLEN];
-#ifdef GMX_MPI
-    int         result;
+        }
+        else
+        {
+            /* Let all other processes wait till the master has printed
+             * the error message and issued MPI_Abort.
+             */
+            MPI_Barrier(MPI_COMM_WORLD);
+        }
+    }
+#else
+    GMX_UNUSED_VALUE(bMaster);
+    GMX_UNUSED_VALUE(bFinalize);
 #endif
 
-#ifdef GMX_MPI
-    /* Check if we are calling on all processes in MPI_COMM_WORLD */
-    if (cr != NULL)
-    {
-        MPI_Comm_compare(cr->mpi_comm_mysim, MPI_COMM_WORLD, &result);
-    }
-    else
+    if (bDebugMode())
     {
-        MPI_Comm_compare(dd->mpi_comm_all, MPI_COMM_WORLD, &result);
+        std::abort();
     }
-    /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
-    const bool bFinalize = (result != MPI_UNEQUAL);
-#else
-    const bool bFinalize = true;
-#endif
+    std::exit(1);
+}
 
-    if ((cr != NULL && MASTER(cr)  ) ||
-        (dd != NULL && DDMASTER(dd)))
+void gmx_fatal_mpi_va(int f_errno, const char *file, int line,
+                      gmx_bool bMaster, gmx_bool bFinalize,
+                      const char *fmt, va_list ap)
+{
+    if (bMaster)
     {
-        va_start(ap, fmt);
+        char msg[STRLEN];
         vsprintf(msg, fmt, ap);
-        va_end(ap);
 
         tMPI_Thread_mutex_lock(&debug_mutex);
         fatal_errno = f_errno;
         tMPI_Thread_mutex_unlock(&debug_mutex);
 
-        if (bFinalize)
-        {
-            /* Use an error handler that does not quit */
-            set_gmx_error_handler(quit_gmx_noquit);
-        }
-
-        _gmx_error("fatal", msg, file, line);
+        call_error_handler("fatal", file, line, msg);
     }
 
-#ifdef GMX_MPI
-    if (bFinalize)
-    {
-        /* Broadcast the fatal error number possibly modified
-         * on the master process, in case the user would like
-         * to use the return status on a non-master process.
-         * The master process in cr and dd always has global rank 0.
-         */
-        MPI_Bcast(&fatal_errno, sizeof(fatal_errno), MPI_BYTE,
-                  0, MPI_COMM_WORLD);
-
-        /* Finalize nicely instead of aborting */
-        MPI_Finalize();
-    }
-    else
-    {
-        /* Let all other processes wait till the master has printed
-         * the error message and issued MPI_Abort.
-         */
-        MPI_Barrier(MPI_COMM_WORLD);
-    }
-#endif
+    do_exit(bMaster, bFinalize);
+}
 
-    exit(fatal_errno);
+void gmx_fatal(int f_errno, const char *file, int line, const char *fmt, ...)
+{
+    va_list ap;
+    va_start(ap, fmt);
+    gmx_fatal_mpi_va(f_errno, file, line, TRUE, FALSE, fmt, ap);
+    va_end(ap);
 }
 
 /*
@@ -339,17 +292,6 @@ void init_debug(const int dbglevel, const char *dbgfile)
     tMPI_Thread_mutex_unlock(&debug_mutex);
 }
 
-static const char *gmxuser = "Please report this to the mailing list (gmx-users@gromacs.org)";
-
-static void        (*gmx_error_handler)(const char *msg) = quit_gmx;
-
-void set_gmx_error_handler(void (*func)(const char *msg))
-{
-    tMPI_Thread_mutex_lock(&debug_mutex);
-    gmx_error_handler = func;
-    tMPI_Thread_mutex_unlock(&debug_mutex);
-}
-
 char *gmx_strerror(const char *key)
 {
     typedef struct {
@@ -393,37 +335,8 @@ char *gmx_strerror(const char *key)
 
 void _gmx_error(const char *key, const char *msg, const char *file, int line)
 {
-    char        buf[10240], errerrbuf[1024];
-    const char *llines = "-------------------------------------------------------";
-    char       *strerr;
-
-    /* protect the audience from suggestive discussions */
-
-    if (msg == NULL)
-    {
-        sprintf(errerrbuf, "Empty fatal_error message. %s", gmxuser);
-    }
-    // In case ProgramInfo is not initialized and there is an issue with the
-    // initialization, fall back to "GROMACS".
-    const char *programName = "GROMACS";
-    try
-    {
-        programName = gmx::getProgramContext().displayName();
-    }
-    catch (const std::exception &)
-    {
-    }
-
-    strerr = gmx_strerror(key);
-    sprintf(buf, "\n%s\nProgram %s, %s\n"
-            "Source code file: %s, line: %d\n\n"
-            "%s:\n%s\nFor more information and tips for troubleshooting, please check the GROMACS\n"
-            "website at http://www.gromacs.org/Documentation/Errors\n%s\n",
-            llines, programName, gmx_version(), file, line,
-            strerr, msg ? msg : errerrbuf, llines);
-    free(strerr);
-
-    gmx_error_handler(buf);
+    call_error_handler(key, file, line, msg);
+    do_exit(true, false);
 }
 
 void _range_check(int n, int n_min, int n_max, const char *warn_str,
index 90808bd0fab29869700d61946f74f23732f7db61..a30f99e49fabb72cb2b15f7e1df3145aab766975 100644 (file)
@@ -44,6 +44,7 @@
 #ifndef GMX_UTILITY_FATALERROR_H
 #define GMX_UTILITY_FATALERROR_H
 
+#include <stdarg.h>
 #include <stdio.h>
 
 #include "basedefinitions.h"
@@ -80,6 +81,22 @@ _where(const char *file, int line);
 /** Prints filename and line to stdlog. */
 #define where() _where(__FILE__, __LINE__)
 
+/*! |brief
+ * Low-level fatal error reporting routine for collective MPI errors.
+ *
+ * This function works as gmx_fatal(), but provides additional control for
+ * cases where it is known that the same error occurs on multiple MPI ranks.
+ * The error handler is called only if \p bMaster is `TRUE`, and MPI_Finalize()
+ * is called instead of MPI_Abort() in MPI-enabled \Gromacs if \p bFinalize is
+ * `TRUE`.
+ *
+ * This is used to implement gmx_fatal_collective() (which cannot be declared
+ * here, since it would bring with it mdrun-specific dependencies).
+ */
+void
+gmx_fatal_mpi_va(int fatal_errno, const char *file, int line, gmx_bool bMaster,
+                 gmx_bool bFinalize, const char *fmt, va_list ap) GMX_ATTRIBUTE_NORETURN;
+
 /*! \brief
  * Fatal error reporting routine for \Gromacs.
  *
@@ -93,7 +110,7 @@ _where(const char *file, int line);
  * The format of \p fmt uses printf()-like formatting.
  *
  * In case all MPI processes want to stop with the same fatal error,
- * use gmx_fatal_collective(), declared in gmx_fatal_collective.h,
+ * use gmx_fatal_collective(), declared in network.h,
  * to avoid having as many error messages as processes.
  *
  * The first three parameters can be provided through ::FARGS:
@@ -197,10 +214,10 @@ void _gmx_error(const char *key, const char *msg, const char *file, int line) GM
 /*! \brief
  * Sets an error handler for gmx_fatal() and other fatal error routines.
  *
- * The default handler prints the message and aborts the program.
- * If you set a custom handler, it must also abort the program, otherwise
- * \Gromacs will behave unpredictably (most likely, it crashes shortly after
- * the fatal error).
+ * The default handler prints the message.
+ * \Gromacs will terminate the program after the error handler returns.
+ * To make gmx_fatal_collective() work, the error handler should not terminate
+ * the program, as it cannot know what is the desired way of termination.
  * The string passed to the handler may be a multi-line string.
  *
  * \see gmx_fatal()
index 97fb51791bdf0b7aa90023c4f4b3c455a36a34ca..845b6665ebc29a43ea86091e695037e82741cbcf 100644 (file)
@@ -73,7 +73,6 @@
 #include "gmx_detect_hardware.h"
 #include "gmx_omp_nthreads.h"
 #include "gromacs/gmxpreprocess/calc_verletbuf.h"
-#include "gmx_fatal_collective.h"
 #include "membed.h"
 #include "macros.h"
 #include "gmx_thread_affinity.h"