Move file opening out of init_pull()
authorBerk Hess <hess@kth.se>
Fri, 23 Mar 2018 11:14:23 +0000 (12:14 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 29 Mar 2018 13:02:03 +0000 (15:02 +0200)
This is a step towards making a constructor for t_pull.

Change-Id: I24583fa32066503ccfbe0abaef6a412241974b12

src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/gmxpreprocess/readir.h
src/gromacs/gmxpreprocess/readpull.cpp
src/gromacs/pulling/pull.cpp
src/gromacs/pulling/pull.h
src/programs/mdrun/runner.cpp

index 394bf80086ab539049d14c4da27261f5c3c48e28..c38b6825a3aca0dd3efa56a7fc6aaac74a6d8a17 100644 (file)
@@ -2252,7 +2252,7 @@ int gmx_grompp(int argc, char *argv[])
 
     if (ir->bPull)
     {
-        pull = set_pull_init(ir, sys, as_rvec_array(state.x.data()), state.box, state.lambda[efptMASS], oenv);
+        pull = set_pull_init(ir, sys, as_rvec_array(state.x.data()), state.box, state.lambda[efptMASS]);
     }
 
     /* Modules that supply external potential for pull coordinates
index 88f04363007803a57b4a15ebf8c1204c466f3b97..50999e15e972c359dc90120b4f9676eae94e99a9 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,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
@@ -147,8 +147,7 @@ void make_pull_coords(pull_params_t *pull);
 /* Process the pull coordinates after reading the pull groups */
 
 pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop,
-                      rvec *x, matrix box, real lambda,
-                      const gmx_output_env_t *oenv);
+                      rvec *x, matrix box, real lambda);
 /* Prints the initial pull group distances in x.
  * If requested, adds the current distance to the initial reference location.
  * Returns the pull_t pull work struct. This should be passed to finish_pull()
index f732e23395498580a08785ab0be07b156f51b708..f06509ddc38a3a01f5e74a0f9183acbe00effa07 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,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
@@ -503,8 +503,7 @@ void make_pull_coords(pull_params_t *pull)
 }
 
 pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop,
-                      rvec *x, matrix box, real lambda,
-                      const gmx_output_env_t *oenv)
+                      rvec *x, matrix box, real lambda)
 {
     pull_params_t *pull;
     pull_t        *pull_work;
@@ -513,7 +512,7 @@ pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop,
     double         t_start;
 
     pull      = ir->pull;
-    pull_work = init_pull(nullptr, pull, ir, 0, nullptr, mtop, nullptr, oenv, lambda, FALSE, ContinuationOptions());
+    pull_work = init_pull(nullptr, pull, ir, mtop, nullptr, lambda);
     auto mdAtoms = gmx::makeMDAtoms(nullptr, *mtop, *ir, false);
     auto md      = mdAtoms->mdatoms();
     atoms2md(mtop, ir, -1, nullptr, mtop->natoms, mdAtoms.get());
index 1625a7a1e83688c42082df10724cdb4a07d3a296..6ab60197335014b17a1540d3f57fb601fd7fcba4 100644 (file)
@@ -2117,11 +2117,8 @@ static void init_pull_group_index(FILE *fplog, const t_commrec *cr,
 
 struct pull_t *
 init_pull(FILE *fplog, const pull_params_t *pull_params, const t_inputrec *ir,
-          int nfile, const t_filenm fnm[],
           const gmx_mtop_t *mtop, const t_commrec *cr,
-          const gmx_output_env_t *oenv, real lambda,
-          gmx_bool bOutFile,
-          const ContinuationOptions &continuationOptions)
+          real lambda)
 {
     struct pull_t *pull;
     pull_comm_t   *comm;
@@ -2504,62 +2501,65 @@ init_pull(FILE *fplog, const pull_params_t *pull_params, const t_inputrec *ir,
     /* We still need to initialize the PBC reference coordinates */
     pull->bSetPBCatoms = TRUE;
 
-    /* Only do I/O when we are doing dynamics and if we are the MASTER */
     pull->out_x = nullptr;
     pull->out_f = nullptr;
-    if (bOutFile)
-    {
-        /* Check for px and pf filename collision, if we are writing
-           both files */
-        std::string px_filename, pf_filename;
-        std::string px_appended, pf_appended;
-        try
-        {
-            px_filename  = std::string(opt2fn("-px", nfile, fnm));
-            pf_filename  = std::string(opt2fn("-pf", nfile, fnm));
-        }
-        GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
 
+    return pull;
+}
 
-        if ((pull->params.nstxout != 0) &&
-            (pull->params.nstfout != 0) &&
-            (px_filename == pf_filename))
+void init_pull_output_files(pull_t                    *pull,
+                            int                        nfile,
+                            const t_filenm             fnm[],
+                            const gmx_output_env_t    *oenv,
+                            const ContinuationOptions &continuationOptions)
+{
+    /* Check for px and pf filename collision, if we are writing
+       both files */
+    std::string px_filename, pf_filename;
+    std::string px_appended, pf_appended;
+    try
+    {
+        px_filename  = std::string(opt2fn("-px", nfile, fnm));
+        pf_filename  = std::string(opt2fn("-pf", nfile, fnm));
+    }
+    GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
+
+    if ((pull->params.nstxout != 0) &&
+        (pull->params.nstfout != 0) &&
+        (px_filename == pf_filename))
+    {
+        if (!opt2bSet("-px", nfile, fnm) && !opt2bSet("-pf", nfile, fnm))
         {
-            if (!opt2bSet("-px", nfile, fnm) && !opt2bSet("-pf", nfile, fnm))
-            {
-                /* We are writing both pull files but neither set directly. */
-                try
-                {
-                    px_appended   = append_before_extension(px_filename, "_pullx");
-                    pf_appended   = append_before_extension(pf_filename, "_pullf");
-                }
-                GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
-                pull->out_x   = open_pull_out(px_appended.c_str(), pull, oenv,
-                                              TRUE, continuationOptions);
-                pull->out_f = open_pull_out(pf_appended.c_str(), pull, oenv,
-                                            FALSE, continuationOptions);
-                return pull;
-            }
-            else
+            /* We are writing both pull files but neither set directly. */
+            try
             {
-                /* If one of -px and -pf is set but the filenames are identical: */
-                gmx_fatal(FARGS, "Identical pull_x and pull_f output filenames %s",
-                          px_filename.c_str());
+                px_appended   = append_before_extension(px_filename, "_pullx");
+                pf_appended   = append_before_extension(pf_filename, "_pullf");
             }
-        }
-        if (pull->params.nstxout != 0)
-        {
-            pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv,
+            GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
+            pull->out_x = open_pull_out(px_appended.c_str(), pull, oenv,
                                         TRUE, continuationOptions);
+            pull->out_f = open_pull_out(pf_appended.c_str(), pull, oenv,
+                                        FALSE, continuationOptions);
+            return;
         }
-        if (pull->params.nstfout != 0)
+        else
         {
-            pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
-                                        FALSE, continuationOptions);
+            /* If at least one of -px and -pf is set but the filenames are identical: */
+            gmx_fatal(FARGS, "Identical pull_x and pull_f output filenames %s",
+                      px_filename.c_str());
         }
     }
-
-    return pull;
+    if (pull->params.nstxout != 0)
+    {
+        pull->out_x = open_pull_out(opt2fn("-px", nfile, fnm), pull, oenv,
+                                    TRUE, continuationOptions);
+    }
+    if (pull->params.nstfout != 0)
+    {
+        pull->out_f = open_pull_out(opt2fn("-pf", nfile, fnm), pull, oenv,
+                                    FALSE, continuationOptions);
+    }
 }
 
 static void destroy_pull_group(pull_group_work_t *pgrp)
index 73bb40586289b59d68963521ad4227e44a9db5dd..7f45fec6f875b217ce520ca85f268fb4c64e1e4d 100644 (file)
@@ -224,27 +224,33 @@ void dd_make_local_pull_groups(const t_commrec *cr,
  * \param fplog       General output file, normally md.log.
  * \param pull_params The pull input parameters containing all pull settings.
  * \param ir          The inputrec.
- * \param nfile       Number of files.
- * \param fnm         Standard filename struct.
  * \param mtop        The topology of the whole system.
  * \param cr          Struct for communication info.
- * \param oenv        Output options.
  * \param lambda      FEP lambda.
- * \param bOutFile    Open output files?
- * \param continuationOptions  Options for continuing from checkpoint file
  */
 struct pull_t *init_pull(FILE                      *fplog,
                          const pull_params_t       *pull_params,
                          const t_inputrec          *ir,
-                         int                        nfile,
-                         const t_filenm             fnm[],
                          const gmx_mtop_t          *mtop,
                          const t_commrec           *cr,
-                         const gmx_output_env_t    *oenv,
-                         real                       lambda,
-                         gmx_bool                   bOutFile,
-                         const ContinuationOptions &continuationOptions);
+                         real                       lambda);
 
+/*! \brief Set up and open the pull output files, when requested.
+ *
+ * NOTE: This should only be called on the master rank and only when
+ *       doing dynamics (e.g. not with energy minimization).
+ *
+ * \param pull        The pull work data struct
+ * \param nfile       Number of files.
+ * \param fnm         Standard filename struct.
+ * \param oenv        Output options.
+ * \param continuationOptions  Options for continuing from checkpoint file
+ */
+void init_pull_output_files(pull_t                    *pull,
+                            int                        nfile,
+                            const t_filenm             fnm[],
+                            const gmx_output_env_t    *oenv,
+                            const ContinuationOptions &continuationOptions);
 
 /*! \brief Close the pull output files.
  *
index 306b21777e267e1d6908f3c2ee7abed74fd583ce..9e8fe71c58d63619af8c1e384fced0983c16dc16 100644 (file)
@@ -1305,10 +1305,14 @@ int Mdrunner::mdrunner()
         {
             /* Initialize pull code */
             inputrec->pull_work =
-                init_pull(fplog, inputrec->pull, inputrec, nfile, fnm,
-                          mtop, cr, oenv, inputrec->fepvals->init_lambda,
-                          EI_DYNAMICS(inputrec->eI) && MASTER(cr),
-                          continuationOptions);
+                init_pull(fplog, inputrec->pull, inputrec,
+                          mtop, cr, inputrec->fepvals->init_lambda);
+            if (EI_DYNAMICS(inputrec->eI) && MASTER(cr))
+            {
+                init_pull_output_files(inputrec->pull_work,
+                                       nfile, fnm, oenv,
+                                       continuationOptions);
+            }
         }
 
         if (inputrec->bRot)