Re-implement gmx tune_pme -cpi
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 20 Jan 2015 11:09:21 +0000 (12:09 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 3 Mar 2015 16:07:57 +0000 (17:07 +0100)
The previous implementation just copied mdrun, but that implementation
does things that are more complicated than tune_pme needs, and in turn
slightly complicated the implementation for mdrun.

Change-Id: Idab32148bc04f3a1c62ab78de8f77179e43317aa

src/gromacs/gmxana/gmx_tune_pme.c
src/gromacs/gmxlib/checkpoint.cpp
src/gromacs/legacyheaders/checkpoint.h
src/programs/mdrun/mdrun.cpp

index e2c3dcc09cfedfa1e97f4f76ace856bee0afc665..d5b4c667edbb1ccd7cf05b2e37c92d5ab9b1481f 100644 (file)
@@ -2430,22 +2430,20 @@ int gmx_tune_pme(int argc, char *argv[])
     create_command_line_snippets(bAppendFiles, bKeepAndNumCPT, bResetCountersHalfWay, presteps,
                                  NFILE, fnm, &cmd_args_bench, &cmd_args_launch, ExtraArgs);
 
-    /* Read in checkpoint file if requested */
+    /* Prepare to use checkpoint file if requested */
     sim_part = 1;
     if (opt2bSet("-cpi", NFILE, fnm))
     {
-        snew(cr, 1);
-        cr->duty = DUTY_PP; /* makes the following routine happy */
-        read_checkpoint_simulation_part(opt2fn("-cpi", NFILE, fnm),
-                                        &sim_part, &cpt_steps, cr,
-                                        FALSE, NFILE, fnm, NULL, NULL);
-        sfree(cr);
-        sim_part++;
-        /* sim_part will now be 1 if no checkpoint file was found */
-        if (sim_part <= 1)
+        const char *filename = opt2fn("-cpi", NFILE, fnm);
+        int         cpt_sim_part;
+        read_checkpoint_part_and_step(filename,
+                                      &cpt_sim_part, &cpt_steps);
+        if (cpt_sim_part == 0)
         {
-            gmx_fatal(FARGS, "Checkpoint file %s not found!", opt2fn("-cpi", NFILE, fnm));
+            gmx_fatal(FARGS, "Checkpoint file %s could not be read!", filename);
         }
+        /* tune_pme will run the next part of the simulation */
+        sim_part = cpt_sim_part + 1;
     }
 
     /* Open performance output file and write header info */
index f53a55d2ae49a1ea664c5ac1ef88c1ca6e114cfb..ec964b3bd9f1966ae60530d60fd909b0bee8116c 100644 (file)
@@ -2258,6 +2258,44 @@ void load_checkpoint(const char *fn, FILE **fplog,
     ir->simulation_part += 1;
 }
 
+void read_checkpoint_part_and_step(const char  *filename,
+                                   int         *simulation_part,
+                                   gmx_int64_t *step)
+{
+    int       file_version;
+    char     *version, *btime, *buser, *bhost, *fprog, *ftime;
+    int       double_prec;
+    int       eIntegrator;
+    int       nppnodes, npme;
+    ivec      dd_nc;
+    int       flags_eks, flags_enh, flags_dfh;
+    double    t;
+    t_state   state;
+    t_fileio *fp;
+
+    if (filename == NULL ||
+        !gmx_fexist(filename) ||
+        (!(fp = gmx_fio_open(filename, "r"))))
+    {
+        *simulation_part = 0;
+        *step            = 0;
+        return;
+    }
+
+    /* Not calling initializing state before use is nasty, but all we
+       do is read into its member variables and throw the struct away
+       again immediately. */
+
+    do_cpt_header(gmx_fio_getxdr(fp), TRUE, &file_version,
+                  &version, &btime, &buser, &bhost, &double_prec, &fprog, &ftime,
+                  &eIntegrator, simulation_part, step, &t, &nppnodes, dd_nc, &npme,
+                  &state.natoms, &state.ngtc, &state.nnhpres, &state.nhchainlength,
+                  &(state.dfhist.nlambda), &state.flags, &flags_eks, &flags_enh, &flags_dfh,
+                  &state.edsamstate.nED, &state.swapstate.eSwapCoords, NULL);
+
+    gmx_fio_close(fp);
+}
+
 static void read_checkpoint_data(t_fileio *fp, int *simulation_part,
                                  gmx_int64_t *step, double *t, t_state *state,
                                  int *nfiles, gmx_file_position_t **outputfiles)
@@ -2499,7 +2537,7 @@ static gmx_bool exist_output_file(const char *fnm_cp, int nfile, const t_filenm
 
 /* 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,
-                                         gmx_int64_t *cpt_step, t_commrec *cr,
+                                         t_commrec *cr,
                                          gmx_bool bAppendReq,
                                          int nfile, const t_filenm fnm[],
                                          const char *part_suffix, gmx_bool *bAddPart)
@@ -2618,10 +2656,6 @@ gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_p
             gmx_bcast(sizeof(*bAddPart), bAddPart, cr);
         }
     }
-    if (NULL != cpt_step)
-    {
-        *cpt_step = step;
-    }
 
     return bAppend;
 }
index 78f329dc232d728fbccc0fea5a932ab7b773f766..52a8ff9ae227ba21fb50612ec315df1e99177498 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.
@@ -91,6 +91,20 @@ void read_checkpoint_trxframe(t_fileio *fp, t_trxframe *fr);
 /* Print the complete contents of checkpoint file fn to out */
 void list_checkpoint(const char *fn, FILE *out);
 
+/* ! \brief Read simulation step and part from a checkpoint file
+ *
+ * Used by tune_pme to handle tuning with a checkpoint file as part of the input.
+ *
+ * \param[in]  filename         Name of checkpoint file
+ * \param[out] simulation_part  The part of the simulation that wrote the checkpoint
+ * \param[out] step             The final step number of the simulation that wrote the checkpoint
+ *
+ * The output variables will both contain 0 if filename is NULL, the file
+ * does not exist, or is not readable. */
+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.
@@ -102,7 +116,7 @@ void list_checkpoint(const char *fn, FILE *out);
  * needs to be added to the output file name.
  */
 gmx_bool read_checkpoint_simulation_part(const char *filename, int *simulation_part,
-                                         gmx_int64_t *step, t_commrec *cr,
+                                         t_commrec *cr,
                                          gmx_bool bAppendReq,
                                          int nfile, const t_filenm fnm[],
                                          const char *part_suffix, gmx_bool *bAddPart);
index 1462f0f4267a4e073b03339a49d106127ab8b861..989c8b4ddcad9f9f37065237a1b092fe80725367 100644 (file)
@@ -491,7 +491,7 @@ int gmx_mdrun(int argc, char *argv[])
         bAppendFiles =
             read_checkpoint_simulation_part(opt2fn_master("-cpi", NFILE,
                                                           fnm, cr),
-                                            &sim_part_fn, NULL, cr,
+                                            &sim_part_fn, cr,
                                             bAppendFiles, NFILE, fnm,
                                             part_suffix, &bAddPart);
         if (sim_part_fn == 0 && MULTIMASTER(cr))