Clean up and improve fatal error handling
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 7 Oct 2015 19:11:02 +0000 (22:11 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 22 Oct 2015 12:14:04 +0000 (14:14 +0200)
- Hide gmx_strerror() within fatalerror.cpp, and clean up the
  implementation.
- Move more error formatting responsibility into the fatal error handler
  in fatalerror.cpp, and make it use the shared formatting code in
  errorformat.h.  This has a side effect of providing automatic line
  wrapping for gmx_fatal().
- Move all output from gmx_abort() to the fatal error output routines.
- Centralize the handling of aborting the program on various error
  conditions to gmx_exit_on_fatal_error().
- Flush various streams before calling std::abort() to try to ensure
  that there is no buffered data left (such as the just-printed fatal
  error).
- Call std::abort() instead of std::exit() in cases where the
  application cannot be guaranteed to be in a clean state: std::exit()
  calls destructors for all global and static singleton objects, which
  can cause issues if, e.g., other threads are still executing code.
- Remove essentially unused code for handling errno in gmx_fatal().
  There was a single instance passing anything else except 0, and the
  way perror() was used would not produce anything useful in this case,
  either.

Change-Id: I0a4c458308cde0e149f986b321caa7cd5cbb1eca

src/gromacs/fileio/tngio.cpp
src/gromacs/gmxana/powerspect.cpp
src/gromacs/utility/basenetwork.cpp
src/gromacs/utility/errorformat.cpp
src/gromacs/utility/exceptions.cpp
src/gromacs/utility/exceptions.h
src/gromacs/utility/fatalerror.cpp
src/gromacs/utility/fatalerror.h
src/gromacs/utility/gmxassert.cpp

index 55d04721118a938d5869d18f31aa28a51bf12b29..1324ecf554eff2b15d979ac5598f1b2dacc79661 100644 (file)
@@ -98,8 +98,7 @@ void gmx_tng_open(const char       *filename,
            no use case for GROMACS handling the non-fatal errors
            gracefully. */
         gmx_fatal(FARGS,
-                  "%s while opening %s for %s",
-                  gmx_strerror("file"),
+                  "File I/O error while opening %s for %s",
                   filename,
                   modeToVerb(mode));
     }
index 096052508d0d9eacb5ee7a08d34ab5d772dd6bd7..2c990608aaab5eb72160fa1d7060043200f4ac05 100644 (file)
@@ -67,12 +67,11 @@ void powerspectavg(real ***intftab, int tsteps, int xbins, int ybins, char **out
     int        n;                   /*time index*/
     int        fy  = ybins/2+1;     /* number of (symmetric) fourier y elements; */
     int        rfl = xbins*fy;      /*length of real - DFT == Symmetric 2D matrix*/
-    int        status;
 
 /*Prepare data structures for FFT, with time averaging of power spectrum*/
-    if ( (status = gmx_fft_init_2d_real(&fftp, xbins, ybins, GMX_FFT_FLAG_NONE) ) != 0)
+    if (gmx_fft_init_2d_real(&fftp, xbins, ybins, GMX_FFT_FLAG_NONE) != 0)
     {
-        gmx_fatal(status, __FILE__, __LINE__, "Error allocating FFT");
+        gmx_fatal(FARGS, "Error allocating FFT");
     }
 
 /*Initialize arrays*/
index ce1197c8cad4ba83d21927dd9117785692d0e0b2..57e7794bffc2812db6438057c84827536db065bd 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, 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.
 #include <cstdlib>
 #include <cstring>
 
-#include <exception>
-
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/gmxmpi.h"
-#include "gromacs/utility/programcontext.h"
 
 gmx_bool gmx_mpi_initialized(void)
 {
@@ -227,27 +224,7 @@ int gmx_physicalnode_id_hash(void)
 #ifdef GMX_LIB_MPI
 void gmx_abort(int errorno)
 {
-    const char *programName = "GROMACS";
-    try
-    {
-        programName = gmx::getProgramContext().displayName();
-    }
-    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 rank %d out of %d\n",
-                     programName, noderank, nnodes);
-    }
-    else
-    {
-        std::fprintf(stderr, "Halting program %s\n", programName);
-    }
-
     MPI_Abort(MPI_COMM_WORLD, errorno);
-    std::exit(errorno);
+    std::abort();
 }
 #endif
index 74cfeabb3261ead57745b56b0f824f418f84f5b1..a6fc1cf2a002f9ac2b81c28adb38e967cd0651a5 100644 (file)
@@ -48,6 +48,7 @@
 #include <cstring>
 
 #include "buildinfo.h"
+#include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/baseversion.h"
 #include "gromacs/utility/programcontext.h"
 #include "gromacs/utility/stringutil.h"
@@ -93,6 +94,10 @@ void printFatalErrorHeader(FILE *fp, const char *title,
     {
         std::fprintf(fp, "Function:    %s\n", func);
     }
+    if (gmx_node_num() > 1)
+    {
+        std::fprintf(fp, "MPI rank:    %d (out of %d)\n", gmx_node_rank(), gmx_node_num());
+    }
     std::fprintf(fp, "\n");
     std::fprintf(fp, "%s:\n", title);
 }
index 47d5ec2e42692793d7b8b631c49aaf24d501a95a..365a4e43b887dba03ecf56b09ba49e0e89ec13b4 100644 (file)
@@ -43,8 +43,6 @@
 
 #include "exceptions.h"
 
-#include "config.h"
-
 #include <cstring>
 
 #include <new>
@@ -58,6 +56,7 @@
 
 #include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/errorcodes.h"
+#include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/stringutil.h"
 #include "gromacs/utility/textwriter.h"
@@ -563,12 +562,22 @@ void formatExceptionMessageToWriter(TextWriter           *writer,
 int processExceptionAtExit(const std::exception & /*ex*/)
 {
     int returnCode = 1;
-#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(returnCode);
-#endif
+    // If we have more than one rank (whether real MPI or thread-MPI),
+    // we cannot currently know whether just one rank or all ranks
+    // actually threw the error, so we need to exit here.
+    // Returning would mean graceful cleanup, which is not possible if
+    // some code is still executing on other ranks/threads.
+    if (gmx_node_num() > 1)
+    {
+        gmx_exit_on_fatal_error(ExitType_Abort, returnCode);
+    }
     return returnCode;
 }
 
+void processExceptionAsFatalError(const std::exception &ex)
+{
+    printFatalErrorMessage(stderr, ex);
+    gmx_exit_on_fatal_error(ExitType_Abort, 1);
+}
+
 } // namespace gmx
index 67b81143b97d2fa4196d9ad81e14cb08bb3a7f93..fcb33a496d624e760e11dfaee8d4b5bb573ace85 100644 (file)
@@ -57,6 +57,8 @@
 #include <boost/exception/exception.hpp>
 #include <boost/exception/info.hpp>
 
+#include "gromacs/utility/basedefinitions.h"
+
 namespace gmx
 {
 
@@ -472,6 +474,15 @@ void formatExceptionMessageToWriter(TextWriter           *writer,
  */
 int processExceptionAtExit(const std::exception &ex);
 
+/*! \brief
+ * Helper function for terminating the program on an exception.
+ *
+ * \param[in] ex  Exception that is the cause for terminating the program.
+ *
+ * Does not throw, and does not return.
+ */
+gmx_noreturn void processExceptionAsFatalError(const std::exception &ex);
+
 /*! \brief
  * Macro for catching exceptions at C++ -> C boundary.
  *
@@ -496,8 +507,7 @@ int processExceptionAtExit(const std::exception &ex);
  */
 #define GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR \
     catch (const std::exception &ex) { \
-        ::gmx::printFatalErrorMessage(stderr, ex); \
-        ::std::exit(::gmx::processExceptionAtExit(ex)); \
+        ::gmx::processExceptionAsFatalError(ex); \
     }
 
 //! \}
index c2b040d51625b24d1b98041c4f9798cd647c10a8..75c2f6b5b6ca602a553dbcadebec032fce682924 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, 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.
@@ -52,6 +52,7 @@
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/baseversion.h"
 #include "gromacs/utility/cstringutil.h"
+#include "gromacs/utility/errorcodes.h"
 #include "gromacs/utility/futil.h"
 #include "gromacs/utility/programcontext.h"
 
@@ -60,6 +61,8 @@
 #include "gromacs/utility/gmxmpi.h"
 #endif
 
+#include "errorformat.h"
+
 static bool                bDebug         = false;
 static tMPI_Thread_mutex_t where_mutex    = TMPI_THREAD_MUTEX_INITIALIZER;
 
@@ -68,8 +71,6 @@ gmx_bool                   gmx_debug_at   = FALSE;
 
 static FILE               *log_file       = NULL;
 static tMPI_Thread_mutex_t error_mutex    = TMPI_THREAD_MUTEX_INITIALIZER;
-static const char *const   gmxuser
-    = "Please report this to the mailing list (gmx-users@gromacs.org)";
 
 void gmx_init_debug(const int dbglevel, const char *dbgfile)
 {
@@ -137,127 +138,120 @@ void gmx_fatal_set_log_file(FILE *fp)
     log_file = fp;
 }
 
-static int fatal_errno = 0;
-
-static void default_error_handler(const char *msg)
+static void default_error_handler(const char *title, const char *msg,
+                                  const char *file, int line)
 {
-    tMPI_Thread_mutex_lock(&error_mutex);
-    if (fatal_errno == 0)
+    if (log_file)
     {
-        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);
+        gmx::internal::printFatalErrorHeader(log_file, title, NULL, file, line);
+        gmx::internal::printFatalErrorMessageLine(log_file, msg, 0);
+        gmx::internal::printFatalErrorFooter(log_file);
     }
-    tMPI_Thread_mutex_unlock(&error_mutex);
+    gmx::internal::printFatalErrorHeader(stderr, title, NULL, file, line);
+    gmx::internal::printFatalErrorMessageLine(stderr, msg, 0);
+    gmx::internal::printFatalErrorFooter(stderr);
 }
 
-static void (*gmx_error_handler)(const char *msg) = default_error_handler;
+static gmx_error_handler_t gmx_error_handler = default_error_handler;
 
-void set_gmx_error_handler(void (*func)(const char *msg))
+void gmx_set_error_handler(gmx_error_handler_t func)
 {
-    // TODO: Either this is unnecessary, or also reads to the handler should be
-    // protected by a mutex.
     tMPI_Thread_mutex_lock(&error_mutex);
     gmx_error_handler = func;
     tMPI_Thread_mutex_unlock(&error_mutex);
 }
 
-static void call_error_handler(const char *key, const char *file, int line, const char *msg)
+static const char *gmx_strerror(const char *key)
 {
-    char        buf[10240], errerrbuf[1024];
-    const char *llines = "-------------------------------------------------------";
-    char       *strerr;
+    struct ErrorKeyEntry {
+        const char *key;
+        const char *msg;
+    };
+    ErrorKeyEntry map[] = {
+        { "call",   "Routine should not have been called" },
+        { "comm",   "Communication (parallel processing) problem" },
+        { "fatal",  "Fatal error" },
+        { "file",   "File input/output error" },
+        { "impl",   "Implementation restriction" },
+        { "incons", "Software inconsistency error" },
+        { "input",  "Input error or input inconsistency" },
+        { "mem",    "Memory allocation/freeing error" },
+        { "open",   "Cannot open file" },
+        { "range",  "Range checking error" }
+    };
 
-    if (msg == NULL)
+    if (key == NULL)
     {
-        sprintf(errerrbuf, "Empty fatal_error message. %s", gmxuser);
+        return "NULL error type (should not occur)";
     }
-    // In case ProgramInfo is not initialized and there is an issue with the
-    // initialization, fall back to "GROMACS".
-    const char *programName = "GROMACS";
-    try
+    for (const ErrorKeyEntry &entry : map)
     {
-        programName = gmx::getProgramContext().displayName();
+        if (std::strcmp(key, entry.key) == 0)
+        {
+            return entry.msg;
+        }
     }
-    catch (const std::exception &)
+    return gmx::getErrorCodeString(gmx::eeUnknownError);
+}
+
+static void call_error_handler(const char *key, const char *file, int line, const char *msg)
+{
+    if (msg == NULL)
     {
+        msg = "Empty gmx_fatal message (bug).";
     }
-
-    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);
+    tMPI_Thread_mutex_lock(&error_mutex);
+    gmx_error_handler(gmx_strerror(key), msg, file, line);
+    tMPI_Thread_mutex_unlock(&error_mutex);
 }
 
-gmx_noreturn static void do_exit(bool bMaster, bool bFinalize)
+void gmx_exit_on_fatal_error(ExitType exitType, int returnValue)
 {
+    if (log_file)
+    {
+        std::fflush(log_file);
+    }
     if (debug)
     {
-        fflush(debug);
+        std::fflush(debug);
     }
+    std::fflush(stdout);
+    std::fflush(stderr);
 
 #ifdef GMX_MPI
     if (gmx_mpi_initialized())
     {
-        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 if (bMaster)
+        switch (exitType)
         {
+            case ExitType_CleanExit:
+                MPI_Finalize();
+                break;
+            case ExitType_Abort:
 #ifdef GMX_LIB_MPI
-            gmx_abort(1);
+                gmx_abort(returnValue);
 #endif
-        }
-        else
-        {
-            /* Let all other processes wait till the master has printed
-             * the error message and issued MPI_Abort.
-             */
-            MPI_Barrier(MPI_COMM_WORLD);
+                break;
+            case ExitType_NonMasterAbort:
+                // Let all other processes wait till the master has printed
+                // the error message and issued MPI_Abort.
+                MPI_Barrier(MPI_COMM_WORLD);
+                break;
         }
     }
-#else
-    GMX_UNUSED_VALUE(bMaster);
-    GMX_UNUSED_VALUE(bFinalize);
 #endif
 
-    if (bDebugMode())
+    if (exitType == ExitType_CleanExit)
     {
-        std::abort();
+        std::exit(returnValue);
     }
-    std::exit(1);
+    // We could possibly use std::_Exit() or other C99/C++11 construct to still
+    // use the exit code, but we cannot use std::exit() if other threads may
+    // still be executing, since that would cause destructors to be called for
+    // global objects that may still be in use elsewhere.
+    std::abort();
 }
 
-void gmx_fatal_mpi_va(int f_errno, const char *file, int line,
+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)
 {
@@ -265,15 +259,15 @@ void gmx_fatal_mpi_va(int f_errno, const char *file, int line,
     {
         char msg[STRLEN];
         vsprintf(msg, fmt, ap);
-
-        tMPI_Thread_mutex_lock(&error_mutex);
-        fatal_errno = f_errno;
-        tMPI_Thread_mutex_unlock(&error_mutex);
-
         call_error_handler("fatal", file, line, msg);
     }
 
-    do_exit(bMaster, bFinalize);
+    ExitType exitType = ExitType_CleanExit;
+    if (!bFinalize)
+    {
+        exitType = bMaster ? ExitType_Abort : ExitType_NonMasterAbort;
+    }
+    gmx_exit_on_fatal_error(exitType, 1);
 }
 
 void gmx_fatal(int f_errno, const char *file, int line, const char *fmt, ...)
@@ -284,51 +278,10 @@ void gmx_fatal(int f_errno, const char *file, int line, const char *fmt, ...)
     va_end(ap);
 }
 
-char *gmx_strerror(const char *key)
-{
-    typedef struct {
-        const char *key, *msg;
-    } error_msg_t;
-    error_msg_t msg[] = {
-        { "bug",    "Possible bug" },
-        { "call",   "Routine should not have been called" },
-        { "comm",   "Communication (parallel processing) problem" },
-        { "fatal",  "Fatal error" },
-        { "cmd",    "Invalid command line argument" },
-        { "file",   "File input/output error" },
-        { "impl",   "Implementation restriction" },
-        { "incons", "Software inconsistency error" },
-        { "input",  "Input error or input inconsistency" },
-        { "mem",    "Memory allocation/freeing error" },
-        { "open",   "Can not open file" },
-        { "range",  "Range checking error" },
-        { NULL,     NULL}
-    };
-
-    if (key == NULL)
-    {
-        return strdup("Empty message");
-    }
-    else
-    {
-        for (size_t i = 0; msg[i].key != NULL; ++i)
-        {
-            if (strcmp(key, msg[i].key) == 0)
-            {
-                return strdup(msg[i].msg);
-            }
-        }
-        char buf[1024];
-        sprintf(buf, "No error message associated with key %s\n%s", key, gmxuser);
-        return strdup(buf);
-    }
-}
-
-
 void _gmx_error(const char *key, const char *msg, const char *file, int line)
 {
     call_error_handler(key, file, line, msg);
-    do_exit(true, false);
+    gmx_exit_on_fatal_error(ExitType_Abort, 1);
 }
 
 void _range_check(int n, int n_min, int n_max, const char *warn_str,
index 64cfe2dc6f14192a457bed397f58c245cd5dc96a..3dba1b8998d5e835050ade7bc35cc3a2d0ea63f4 100644 (file)
@@ -92,6 +92,9 @@ _where(const char *file, int line);
 void
 gmx_fatal_set_log_file(FILE *fp);
 
+/** Function pointer type for fatal error handler callback. */
+typedef void (*gmx_error_handler_t)(const char *title, const char *msg, const char *file, int line);
+
 /*! \brief
  * Sets an error handler for gmx_fatal() and other fatal error routines.
  *
@@ -99,12 +102,47 @@ gmx_fatal_set_log_file(FILE *fp);
  * \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.
+ * The message passed to the handler may be a multi-line string.
  *
  * \see gmx_fatal()
  */
-void
-set_gmx_error_handler(void (*func)(const char *msg));
+void gmx_set_error_handler(gmx_error_handler_t func);
+
+/** Identifies the state of the program on a fatal error. */
+enum ExitType
+{
+    /*! \brief
+     * Clean exit is possible.
+     *
+     * There should be no concurrently executing code that might be accessing
+     * global objects, and all MPI ranks should reach the same fatal error.
+     */
+    ExitType_CleanExit,
+    /*! \brief
+     * Program needs to be aborted.
+     *
+     * There are no preconditions for this state.
+     */
+    ExitType_Abort,
+    /*! \brief
+     * Program needs to be aborted, but some other rank is responsible of it.
+     *
+     * There should be some other MPI rank that reaches the same fatal error,
+     * but uses ExitType_Abort.  The other ranks can then use
+     * ExitType_NonMasterAbort to wait for that one rank to issue the abort.
+     */
+    ExitType_NonMasterAbort
+};
+
+/*! \brief
+ * Helper function to terminate the program on a fatal error.
+ *
+ * \param[in] exitType  Identifies the state of the program at the time of the
+ *    call, determining how the program can be terminated.
+ * \param[in] returnValue  Exit code for the program, for cases where it can be
+ *    used.
+ */
+gmx_noreturn void gmx_exit_on_fatal_error(enum ExitType exitType, int returnValue);
 
 /*! \brief
  * Low-level fatal error reporting routine for collective MPI errors.
@@ -119,8 +157,9 @@ set_gmx_error_handler(void (*func)(const char *msg));
  * here, since it would bring with it mdrun-specific dependencies).
  */
 gmx_noreturn 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_fatal_mpi_va(int fatal_errno, const char *file, int line,
+                 gmx_bool bMaster, gmx_bool bFinalize,
+                 const char *fmt, va_list ap);
 
 /*! \brief
  * Fatal error reporting routine for \Gromacs.
@@ -148,14 +187,6 @@ gmx_fatal(int fatal_errno, const char *file, int line, const char *fmt, ...);
 /** Helper macro to pass first three parameters to gmx_fatal(). */
 #define FARGS 0, __FILE__, __LINE__
 
-/*! \brief
- * Returns error message corresponding to a string key.
- *
- * This maps the strings used by gmx_error() to actual error messages.
- * Caller is responsible of freeing the returned string.
- */
-char *gmx_strerror(const char *key);
-
 /** Implementation for gmx_error(). */
 gmx_noreturn void _gmx_error(const char *key, const char *msg, const char *file, int line);
 /*! \brief
index ef61240a7635b0e648ac5dd56029894efcc91e2e..9ab6a11fc4cc7193f8a397ae133b2c16ca3117ef 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2011,2012,2014, by the GROMACS development team, led by
+ * Copyright (c) 2011,2012,2014,2015, 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.
@@ -46,6 +46,8 @@
 #include <cstdio>
 #include <cstdlib>
 
+#include "gromacs/utility/fatalerror.h"
+
 #include "errorformat.h"
 
 namespace gmx
@@ -62,7 +64,7 @@ void assertHandler(const char *condition, const char *msg,
     std::fprintf(stderr, "Condition: %s\n", condition);
     printFatalErrorMessageLine(stderr, msg, 0);
     printFatalErrorFooter(stderr);
-    std::abort();
+    gmx_exit_on_fatal_error(ExitType_Abort, 1);
 }
 
 }   // namespace internal