Move code for managing mdrun -cpi and -append
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 20 Jan 2015 11:05:31 +0000 (12:05 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 16 Jun 2015 20:24:36 +0000 (22:24 +0200)
Prepares for refactoring to simplify and fix this handling. Kept logic
peculiar to handling mdrun command-line options near other such code.
Renamed functions to better suit their roles. Renamed bAppend
variables, in particular to separate the role of "the user asked for
appending" from "whether mdrun is actually going to do appending."

Change-Id: I063138a49d210e92cebe19e290dc6a0e7b72b1cf

src/gromacs/CMakeLists.txt
src/gromacs/gmxlib/checkpoint.cpp
src/gromacs/legacyheaders/checkpoint.h
src/gromacs/mdrunutility/CMakeLists.txt [new file with mode: 0644]
src/gromacs/mdrunutility/handlerestart.cpp [new file with mode: 0644]
src/gromacs/mdrunutility/handlerestart.h [new file with mode: 0644]
src/programs/mdrun/mdrun.cpp

index 3f426e02a3cd52f9fea61ee26c8b9ead4893e80f..1ef3f1e632fdc3837475aa8f01b0f11d8997e196 100644 (file)
@@ -98,6 +98,7 @@ add_subdirectory(ewald)
 add_subdirectory(fft)
 add_subdirectory(linearalgebra)
 add_subdirectory(math)
+add_subdirectory(mdrunutility)
 add_subdirectory(random)
 add_subdirectory(onlinehelp)
 add_subdirectory(options)
index dfe59a730939bc1958d6b40ac56cd453cd8038db..98897daaa66dabef938ac0a7ea9bd1540db568ce 100644 (file)
@@ -2517,145 +2517,24 @@ void list_checkpoint(const char *fn, FILE *out)
     done_state(&state);
 }
 
-
-static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
-{
-    int i;
-
-    /* Check if the output file name stored in the checkpoint file
-     * is one of the output file names of mdrun.
-     */
-    i = 0;
-    while (i < nfile &&
-           !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
-    {
-        i++;
-    }
-
-    return (i < nfile && gmx_fexist(fnm_cp));
-}
-
 /* This routine cannot print tons of data, since it is called before the log file is opened. */
-gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
-                                         t_commrec *cr,
-                                         gmx_bool bAppendReq,
-                                         int nfile, const t_filenm fnm[],
-                                         const char *part_suffix, gmx_bool *bAddPart)
+void
+read_checkpoint_simulation_part_and_filenames(t_fileio             *fp,
+                                              int                  *simulation_part,
+                                              int                  *nfiles,
+                                              gmx_file_position_t **outputfiles)
 {
-    t_fileio            *fp;
-    gmx_int64_t          step = 0;
-    double               t;
-    /* This next line is nasty because the sub-structures of t_state
-     * cannot be assumed to be zeroed (or even initialized in ways the
-     * rest of the code might assume). Using snew would be better, but
-     * this will all go away for 5.0. */
-    t_state              state;
-    int                  nfiles;
-    gmx_file_position_t *outputfiles;
-    int                  nexist, f;
-    gmx_bool             bAppend;
-    char                *fn, suf_up[STRLEN];
-
-    bAppend = FALSE;
+    gmx_int64_t step = 0;
+    double      t;
+    t_state     state;
 
-    if (SIMMASTER(cr))
-    {
-        if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
-        {
-            *simulation_part = 0;
-        }
-        else
-        {
-            init_state(&state, 0, 0, 0, 0, 0);
-
-            read_checkpoint_data(fp, simulation_part, &step, &t, &state,
-                                 &nfiles, &outputfiles);
-            if (gmx_fio_close(fp) != 0)
-            {
-                gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
-            }
-            done_state(&state);
-
-            if (bAppendReq)
-            {
-                nexist = 0;
-                for (f = 0; f < nfiles; f++)
-                {
-                    if (exist_output_file(outputfiles[f].filename, nfile, fnm))
-                    {
-                        nexist++;
-                    }
-                }
-                if (nexist == nfiles)
-                {
-                    bAppend = bAppendReq;
-                }
-                else if (nexist > 0)
-                {
-                    fprintf(stderr,
-                            "Output file appending has been requested,\n"
-                            "but some output files listed in the checkpoint file %s\n"
-                            "are not present or are named differently by the current program:\n",
-                            filename);
-                    fprintf(stderr, "output files present:");
-                    for (f = 0; f < nfiles; f++)
-                    {
-                        if (exist_output_file(outputfiles[f].filename,
-                                              nfile, fnm))
-                        {
-                            fprintf(stderr, " %s", outputfiles[f].filename);
-                        }
-                    }
-                    fprintf(stderr, "\n");
-                    fprintf(stderr, "output files not present or named differently:");
-                    for (f = 0; f < nfiles; f++)
-                    {
-                        if (!exist_output_file(outputfiles[f].filename,
-                                               nfile, fnm))
-                        {
-                            fprintf(stderr, " %s", outputfiles[f].filename);
-                        }
-                    }
-                    fprintf(stderr, "\n");
-
-                    gmx_fatal(FARGS, "File appending requested, but %d of the %d output files are not present or are named differently", nfiles-nexist, nfiles);
-                }
-            }
-
-            if (bAppend)
-            {
-                if (nfiles == 0)
-                {
-                    gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
-                }
-                fn = outputfiles[0].filename;
-                if (strlen(fn) < 4 ||
-                    gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
-                {
-                    gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
-                }
-                /* Set bAddPart to whether the suffix string '.part' is present
-                 * in the log file name.
-                 */
-                strcpy(suf_up, part_suffix);
-                upstring(suf_up);
-                *bAddPart = (strstr(fn, part_suffix) != NULL ||
-                             strstr(fn, suf_up) != NULL);
-            }
+    init_state(&state, 0, 0, 0, 0, 0);
 
-            sfree(outputfiles);
-        }
-    }
-    if (PAR(cr))
+    read_checkpoint_data(fp, simulation_part, &step, &t, &state,
+                         nfiles, outputfiles);
+    if (gmx_fio_close(fp) != 0)
     {
-        gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
-
-        if (*simulation_part > 0 && bAppendReq)
-        {
-            gmx_bcast(sizeof(bAppend), &bAppend, cr);
-            gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
-        }
+        gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
     }
-
-    return bAppend;
+    done_state(&state);
 }
index 52a8ff9ae227ba21fb50612ec315df1e99177498..c85a99fa2f3954303d804b8fc58892ffa340d3a7 100644 (file)
@@ -105,21 +105,20 @@ void read_checkpoint_part_and_step(const char  *filename,
                                    int         *simulation_part,
                                    gmx_int64_t *step);
 
-/* Read just the simulation 'generation' and with bAppendReq check files.
- * This is necessary already at the beginning of mdrun,
- * to be able to rename the logfile correctly.
- * When file appending is requested, checks which output files are present:
- * all present: return TRUE,
- * none present: return FALSE,
- * part present: fatal error.
- * When TRUE is returned, bAddPart will tell whether the simulation part
- * needs to be added to the output file name.
+/* ! \brief Read simulation part and output filenames from a checkpoint file
+ *
+ * Used by mdrun to handle restarts
+ *
+ * \param[in]  fp               Handle to open checkpoint file
+ * \param[out] simulation_part  The part of the simulation that wrote the checkpoint
+ * \param[out] nfiles           Number of output files from the previous run
+ * \param[out] outputfiles      Pointer to array of output file names from the previous run. Pointer is allocated in this function.
  */
-gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
-                                         t_commrec *cr,
-                                         gmx_bool bAppendReq,
-                                         int nfile, const t_filenm fnm[],
-                                         const char *part_suffix, gmx_bool *bAddPart);
+void
+read_checkpoint_simulation_part_and_filenames(t_fileio             *fp,
+                                              int                  *simulation_part,
+                                              int                  *nfiles,
+                                              gmx_file_position_t **outputfiles);
 
 #ifdef __cplusplus
 }
diff --git a/src/gromacs/mdrunutility/CMakeLists.txt b/src/gromacs/mdrunutility/CMakeLists.txt
new file mode 100644 (file)
index 0000000..8f00f22
--- /dev/null
@@ -0,0 +1,40 @@
+#
+# This file is part of the GROMACS molecular simulation package.
+#
+# Copyright (c) 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.
+#
+# 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.
+
+file(GLOB MDRUNUTILITY_SOURCES *.cpp)
+set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${MDRUNUTILITY_SOURCES} PARENT_SCOPE)
+
+if (BUILD_TESTING)
+#    add_subdirectory(tests)
+endif()
diff --git a/src/gromacs/mdrunutility/handlerestart.cpp b/src/gromacs/mdrunutility/handlerestart.cpp
new file mode 100644 (file)
index 0000000..3038662
--- /dev/null
@@ -0,0 +1,286 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ *
+ * 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.
+ */
+/*! \internal \file
+ *
+ * \brief This file declares functions for mdrun to call to manage the
+ * details of doing a restart (ie. reading checkpoints, appending
+ * output files).
+ *
+ * \todo Clean up the error-prone logic here. Add doxygen.
+ *
+ * \author Berk Hess <hess@kth.se>
+ * \author Erik Lindahl <erik@kth.se>
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ *
+ * \ingroup module_mdrunutility
+ */
+
+#include "gmxpre.h"
+
+#include "handlerestart.h"
+
+#include <string.h>
+
+#include "gromacs/fileio/filenm.h"
+#include "gromacs/fileio/gmxfio.h"
+#include "gromacs/legacyheaders/checkpoint.h"
+#include "gromacs/legacyheaders/main.h"
+#include "gromacs/legacyheaders/types/commrec.h"
+#include "gromacs/utility/basedefinitions.h"
+#include "gromacs/utility/smalloc.h"
+
+/*! \brief Search for \p fnm_cp in fnm and return true iff found
+ *
+ * \todo This could be implemented sanely with a for loop. */
+static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm fnm[])
+{
+    int i;
+
+    /* Check if the output file name stored in the checkpoint file
+     * is one of the output file names of mdrun.
+     */
+    i = 0;
+    while (i < nfile &&
+           !(is_output(&fnm[i]) && strcmp(fnm_cp, fnm[i].fns[0]) == 0))
+    {
+        i++;
+    }
+
+    return (i < nfile && gmx_fexist(fnm_cp));
+}
+
+/*! \brief Support handling restarts
+ *
+ * \todo Clean this up (next patch)
+ *
+ * Read just the simulation 'generation' and with bTryToAppendFiles check files.
+ * This is is needed at the beginning of mdrun,
+ * to be able to rename the logfile correctly.
+ * When file appending is requested, checks which output files are present,
+ * and returns TRUE/FALSE in bDoAppendFiles if all or none are present.
+ * If only some output files are present, give a fatal error.
+ * When bDoAppendFiles is TRUE upon return, bAddPart will tell whether the simulation part
+ * needs to be added to the output file name.
+ *
+ * This routine cannot print tons of data, since it is called before
+ * the log file is opened. */
+static void
+read_checkpoint_data(const char *filename, int *simulation_part,
+                     t_commrec *cr,
+                     gmx_bool bTryToAppendFiles,
+                     int nfile, const t_filenm fnm[],
+                     const char *part_suffix,
+                     gmx_bool *bAddPart,
+                     gmx_bool *bDoAppendFiles)
+{
+    t_fileio            *fp;
+    int                  nfiles;
+    gmx_file_position_t *outputfiles;
+    int                  nexist, f;
+    char                *fn, suf_up[STRLEN];
+
+    *bDoAppendFiles = FALSE;
+
+    if (SIMMASTER(cr))
+    {
+        if (!gmx_fexist(filename) || (!(fp = gmx_fio_open(filename, "r")) ))
+        {
+            *simulation_part = 0;
+        }
+        else
+        {
+            read_checkpoint_simulation_part_and_filenames(fp,
+                                                          simulation_part,
+                                                          &nfiles,
+                                                          &outputfiles);
+
+            if (bTryToAppendFiles)
+            {
+                nexist = 0;
+                for (f = 0; f < nfiles; f++)
+                {
+                    if (exist_output_file(outputfiles[f].filename, nfile, fnm))
+                    {
+                        nexist++;
+                    }
+                }
+                if (nexist == nfiles)
+                {
+                    *bDoAppendFiles = bTryToAppendFiles;
+                }
+                else if (nexist > 0)
+                {
+                    fprintf(stderr,
+                            "Output file appending has been requested,\n"
+                            "but some output files listed in the checkpoint file %s\n"
+                            "are not present or are named differently by the current program:\n",
+                            filename);
+                    fprintf(stderr, "output files present:");
+                    for (f = 0; f < nfiles; f++)
+                    {
+                        if (exist_output_file(outputfiles[f].filename,
+                                              nfile, fnm))
+                        {
+                            fprintf(stderr, " %s", outputfiles[f].filename);
+                        }
+                    }
+                    fprintf(stderr, "\n");
+                    fprintf(stderr, "output files not present or named differently:");
+                    for (f = 0; f < nfiles; f++)
+                    {
+                        if (!exist_output_file(outputfiles[f].filename,
+                                               nfile, fnm))
+                        {
+                            fprintf(stderr, " %s", outputfiles[f].filename);
+                        }
+                    }
+                    fprintf(stderr, "\n");
+
+                    gmx_fatal(FARGS, "File appending requested, but %d of the %d output files are not present or are named differently", nfiles-nexist, nfiles);
+                }
+            }
+
+            if (*bDoAppendFiles)
+            {
+                if (nfiles == 0)
+                {
+                    gmx_fatal(FARGS, "File appending requested, but no output file information is stored in the checkpoint file");
+                }
+                fn = outputfiles[0].filename;
+                if (strlen(fn) < 4 ||
+                    gmx_strcasecmp(fn+strlen(fn)-4, ftp2ext(efLOG)) == 0)
+                {
+                    gmx_fatal(FARGS, "File appending requested, but the log file is not the first file listed in the checkpoint file");
+                }
+                /* Set bAddPart to whether the suffix string '.part' is present
+                 * in the log file name.
+                 */
+                strcpy(suf_up, part_suffix);
+                upstring(suf_up);
+                *bAddPart = (strstr(fn, part_suffix) != NULL ||
+                             strstr(fn, suf_up) != NULL);
+            }
+
+            sfree(outputfiles);
+        }
+    }
+    if (PAR(cr))
+    {
+        gmx_bcast(sizeof(*simulation_part), simulation_part, cr);
+
+        if (*simulation_part > 0 && bTryToAppendFiles)
+        {
+            gmx_bcast(sizeof(*bDoAppendFiles), bDoAppendFiles, cr);
+            gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
+        }
+    }
+}
+
+/* This routine cannot print tons of data, since it is called before the log file is opened. */
+void
+handleRestart(t_commrec *cr,
+              gmx_bool   bTryToAppendFiles,
+              const int  NFILE,
+              t_filenm   fnm[],
+              gmx_bool  *bDoAppendFiles,
+              gmx_bool  *bStartFromCpt)
+{
+    gmx_bool        bAddPart;
+    int             sim_part, sim_part_fn;
+    const char     *part_suffix = ".part";
+    FILE           *fpmulti;
+
+    bAddPart = !bTryToAppendFiles;
+
+    /* Check if there is ANY checkpoint file available */
+    sim_part    = 1;
+    sim_part_fn = sim_part;
+    if (opt2bSet("-cpi", NFILE, fnm))
+    {
+        read_checkpoint_data(opt2fn_master("-cpi", NFILE, fnm, cr),
+                             &sim_part_fn, cr,
+                             bTryToAppendFiles, NFILE, fnm,
+                             part_suffix, &bAddPart,
+                             bDoAppendFiles);
+        if (sim_part_fn == 0 && MULTIMASTER(cr))
+        {
+            fprintf(stdout, "No previous checkpoint file present, assuming this is a new run.\n");
+        }
+        else
+        {
+            sim_part = sim_part_fn + 1;
+        }
+
+        if (MULTISIM(cr) && MASTER(cr))
+        {
+            if (MULTIMASTER(cr))
+            {
+                /* Log file is not yet available, so if there's a
+                 * problem we can only write to stderr. */
+                fpmulti = stderr;
+            }
+            else
+            {
+                fpmulti = NULL;
+            }
+            check_multi_int(fpmulti, cr->ms, sim_part, "simulation part", TRUE);
+        }
+    }
+    else
+    {
+        *bDoAppendFiles = FALSE;
+    }
+    *bStartFromCpt = sim_part > 1;
+
+    if (!*bDoAppendFiles)
+    {
+        sim_part_fn = sim_part;
+    }
+
+    if (bAddPart)
+    {
+        char suffix[STRLEN];
+
+        /* Rename all output files (except checkpoint files) */
+        /* create new part name first (zero-filled) */
+        sprintf(suffix, "%s%04d", part_suffix, sim_part_fn);
+
+        add_suffix_to_output_names(fnm, NFILE, suffix);
+        if (MULTIMASTER(cr))
+        {
+            fprintf(stdout, "Checkpoint file is from part %d, new output files will be suffixed '%s'.\n", sim_part-1, suffix);
+        }
+    }
+}
diff --git a/src/gromacs/mdrunutility/handlerestart.h b/src/gromacs/mdrunutility/handlerestart.h
new file mode 100644 (file)
index 0000000..4d315cd
--- /dev/null
@@ -0,0 +1,91 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ *
+ * 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.
+ */
+/*! \defgroup module_mdrunutility Implementation of mdrun utility functionality
+ * \ingroup group_mdrun
+ *
+ * \brief This module contains code that implements general
+ * infrastructure for mdrun that does not suit any other module.
+ */
+/*! \libinternal \file
+ *
+ * \brief This file declares functions for mdrun to call to manage the
+ * details of doing a restart (ie. reading checkpoints, appending
+ * output files).
+ *
+ * \todo There may be other code in runner.cpp etc. that can usefully
+ * live here
+ *
+ * \author Berk Hess <hess@kth.se>
+ * \author Erik Lindahl <erik@kth.se>
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ *
+ * \inlibraryapi
+ * \ingroup module_mdrunutility
+ */
+
+#ifndef GMX_MDRUNUTILITY_HANDLERESTART_H
+#define GMX_MDRUNUTILITY_HANDLERESTART_H
+
+#include "gromacs/fileio/filenm.h"
+#include "gromacs/legacyheaders/types/commrec.h"
+#include "gromacs/utility/basedefinitions.h"
+
+/*! \brief Handle startup of mdrun, particularly regarding -cpi and -append
+ *
+ * If there is a checkpoint file, then prepare to start from that
+ * state. If there is also appending, then organize the file naming
+ * accordingly, if possible. Issues various fatal errors if the input
+ * conditions are inconsistent or broken. \p fnm is updated with
+ * suffix strings for part numbers if we are doing a restart from
+ * checkpoint and are not appending.
+ *
+ * Does communication to coordinate behaviour between all ranks of a
+ * simulation, and/or simulations.
+ *
+ * \param[in]    cr                 Communication structure
+ * \param[in]    bTryToAppendFiles  Whether mdrun -append was used
+ * \param[in]    NFILE              Size of fnm struct
+ * \param[inout] fnm                Filename parameters to mdrun
+ * \param[out]   bDoAppendFiles     Whether mdrun will append to files
+ * \param[out]   bStartFromCpt      Whether mdrun will start from the -cpi file
+ */
+void handleRestart(t_commrec *cr,
+                   gmx_bool   bTryToAppendFiles,
+                   const int  NFILE,
+                   t_filenm   fnm[],
+                   gmx_bool  *bDoAppendFiles,
+                   gmx_bool  *bStartFromCpt);
+
+#endif
index c584010d3d8f99f71cb6ab540391b0b22fa172ea..f924821d93c1c7cc20bb9ccacce4ec9c8641de03 100644 (file)
@@ -59,7 +59,6 @@
 
 #include "gromacs/commandline/pargs.h"
 #include "gromacs/fileio/filenm.h"
-#include "gromacs/legacyheaders/checkpoint.h"
 #include "gromacs/legacyheaders/macros.h"
 #include "gromacs/legacyheaders/main.h"
 #include "gromacs/legacyheaders/mdrun.h"
@@ -67,6 +66,7 @@
 #include "gromacs/legacyheaders/readinp.h"
 #include "gromacs/legacyheaders/typedefs.h"
 #include "gromacs/legacyheaders/types/commrec.h"
+#include "gromacs/mdrunutility/handlerestart.h"
 #include "gromacs/utility/fatalerror.h"
 
 #include "mdrun_main.h"
@@ -299,7 +299,7 @@ int gmx_mdrun(int argc, char *argv[])
     real            rdd                   = 0.0, rconstr = 0.0, dlb_scale = 0.8, pforce = -1;
     char           *ddcsx                 = NULL, *ddcsy = NULL, *ddcsz = NULL;
     real            cpt_period            = 15.0, max_hours = -1;
-    gmx_bool        bAppendFiles          = TRUE;
+    gmx_bool        bTryToAppendFiles     = TRUE;
     gmx_bool        bKeepAndNumCPT        = FALSE;
     gmx_bool        bResetCountersHalfWay = FALSE;
     output_env_t    oenv                  = NULL;
@@ -382,7 +382,7 @@ int gmx_mdrun(int argc, char *argv[])
           "Checkpoint interval (minutes)" },
         { "-cpnum",   FALSE, etBOOL, {&bKeepAndNumCPT},
           "Keep and number checkpoint files" },
-        { "-append",  FALSE, etBOOL, {&bAppendFiles},
+        { "-append",  FALSE, etBOOL, {&bTryToAppendFiles},
           "Append to previous output files when continuing from checkpoint instead of adding the simulation part number to all file names" },
         { "-nsteps",  FALSE, etINT64, {&nsteps},
           "Run this number of steps, overrides .mdp file option (-1 means infinite, -2 means use mdp option, smaller is invalid)" },
@@ -418,11 +418,8 @@ int gmx_mdrun(int argc, char *argv[])
     unsigned long   Flags;
     ivec            ddxyz;
     int             dd_node_order;
-    gmx_bool        bAddPart;
-    FILE           *fplog, *fpmulti;
-    int             sim_part, sim_part_fn;
-    const char     *part_suffix = ".part";
-    char            suffix[STRLEN];
+    gmx_bool        bDoAppendFiles, bStartFromCpt;
+    FILE           *fplog;
     int             rc;
     char          **multidir = NULL;
 
@@ -498,65 +495,8 @@ int gmx_mdrun(int argc, char *argv[])
 #endif
     }
 
-    bAddPart = !bAppendFiles;
-
-    /* Check if there is ANY checkpoint file available */
-    sim_part    = 1;
-    sim_part_fn = sim_part;
-    if (opt2bSet("-cpi", NFILE, fnm))
-    {
-        bAppendFiles =
-            read_checkpoint_simulation_part(opt2fn_master("-cpi", NFILE,
-                                                          fnm, cr),
-                                            &sim_part_fn, cr,
-                                            bAppendFiles, NFILE, fnm,
-                                            part_suffix, &bAddPart);
-        if (sim_part_fn == 0 && MULTIMASTER(cr))
-        {
-            fprintf(stdout, "No previous checkpoint file present, assuming this is a new run.\n");
-        }
-        else
-        {
-            sim_part = sim_part_fn + 1;
-        }
-
-        if (MULTISIM(cr) && MASTER(cr))
-        {
-            if (MULTIMASTER(cr))
-            {
-                /* Log file is not yet available, so if there's a
-                 * problem we can only write to stderr. */
-                fpmulti = stderr;
-            }
-            else
-            {
-                fpmulti = NULL;
-            }
-            check_multi_int(fpmulti, cr->ms, sim_part, "simulation part", TRUE);
-        }
-    }
-    else
-    {
-        bAppendFiles = FALSE;
-    }
-
-    if (!bAppendFiles)
-    {
-        sim_part_fn = sim_part;
-    }
-
-    if (bAddPart)
-    {
-        /* Rename all output files (except checkpoint files) */
-        /* create new part name first (zero-filled) */
-        sprintf(suffix, "%s%04d", part_suffix, sim_part_fn);
-
-        add_suffix_to_output_names(fnm, NFILE, suffix);
-        if (MULTIMASTER(cr))
-        {
-            fprintf(stdout, "Checkpoint file is from part %d, new output files will be suffixed '%s'.\n", sim_part-1, suffix);
-        }
-    }
+    handleRestart(cr, bTryToAppendFiles, NFILE, fnm,
+                  &bDoAppendFiles, &bStartFromCpt);
 
     Flags = opt2bSet("-rerun", NFILE, fnm) ? MD_RERUN : 0;
     Flags = Flags | (bDDBondCheck  ? MD_DDBONDCHECK  : 0);
@@ -565,10 +505,10 @@ int gmx_mdrun(int argc, char *argv[])
     Flags = Flags | (bConfout      ? MD_CONFOUT      : 0);
     Flags = Flags | (bRerunVSite   ? MD_RERUN_VSITE  : 0);
     Flags = Flags | (bReproducible ? MD_REPRODUCIBLE : 0);
-    Flags = Flags | (bAppendFiles  ? MD_APPENDFILES  : 0);
+    Flags = Flags | (bDoAppendFiles ? MD_APPENDFILES  : 0);
     Flags = Flags | (opt2parg_bSet("-append", asize(pa), pa) ? MD_APPENDFILESSET : 0);
     Flags = Flags | (bKeepAndNumCPT ? MD_KEEPANDNUMCPT : 0);
-    Flags = Flags | (sim_part > 1    ? MD_STARTFROMCPT : 0);
+    Flags = Flags | (bStartFromCpt ? MD_STARTFROMCPT : 0);
     Flags = Flags | (bResetCountersHalfWay ? MD_RESETCOUNTERSHALFWAY : 0);
     Flags = Flags | (bIMDwait      ? MD_IMDWAIT      : 0);
     Flags = Flags | (bIMDterm      ? MD_IMDTERM      : 0);
@@ -577,7 +517,7 @@ int gmx_mdrun(int argc, char *argv[])
     /* We postpone opening the log file if we are appending, so we can
        first truncate the old log file and append to the correct position
        there instead.  */
-    if (MASTER(cr) && !bAppendFiles)
+    if (MASTER(cr) && !bDoAppendFiles)
     {
         gmx_log_open(ftp2fn(efLOG, NFILE, fnm), cr,
                      Flags & MD_APPENDFILES, &fplog);
@@ -601,7 +541,7 @@ int gmx_mdrun(int argc, char *argv[])
 
     /* Log file has to be closed in mdrunner if we are appending to it
        (fplog not set here) */
-    if (MASTER(cr) && !bAppendFiles)
+    if (MASTER(cr) && !bDoAppendFiles)
     {
         gmx_log_close(fplog);
     }