Clarify distinct paths for mdrun log file handling
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 25 Sep 2018 11:48:22 +0000 (13:48 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 26 Sep 2018 19:12:53 +0000 (21:12 +0200)
Opening for appending is actually handled in the checkpointing code,
which is now explicit in mdrunner.

Also handled the error case when the log file cannot be opened.

Also removed unnecessary dependency on commrec.

Refs #2651

Change-Id: I0b8c7756a08c3d786571e7936e7f8f327fd17947

src/gromacs/mdrun/logging.cpp
src/gromacs/mdrun/logging.h
src/gromacs/mdrun/runner.cpp
src/programs/mdrun/mdrun.cpp

index b1707c0347879356e21db40f82e3c5b9e5d35ef9..8343e98d8a80d98a0dffd221d15fde4d9c4cb210 100644 (file)
@@ -44,7 +44,6 @@
 #include "logging.h"
 
 #include "gromacs/fileio/gmxfio.h"
-#include "gromacs/mdtypes/commrec.h"
 #include "gromacs/utility/binaryinformation.h"
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/exceptions.h"
 #include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/sysinfo.h"
 
-void gmx_log_open(const char *lognm, const t_commrec *cr,
-                  gmx_bool bAppendFiles, FILE** fplog)
+//! Implements aspects of logfile handling common to opening either for writing or appending.
+static void gmx_log_setup(gmx::BinaryInformationSettings settings,
+                          const int                      rankIndex,
+                          const int                      numRanks,
+                          FILE                         * fplog)
 {
     int    pid;
     char   host[256];
     char   timebuf[STRLEN];
-    FILE  *fp = *fplog;
 
-    if (!bAppendFiles)
-    {
-        fp = gmx_fio_fopen(lognm, bAppendFiles ? "a+" : "w+" );
-    }
-
-    gmx_fatal_set_log_file(fp);
+    GMX_RELEASE_ASSERT(fplog != nullptr, "Log file must be already open");
+    gmx_fatal_set_log_file(fplog);
 
     /* Get some machine parameters */
     gmx_gethostname(host, 256);
     pid = gmx_getpid();
     gmx_format_current_time(timebuf, STRLEN);
 
-    if (bAppendFiles)
-    {
-        fprintf(fp,
-                "\n"
-                "\n"
-                "-----------------------------------------------------------\n"
-                "Restarting from checkpoint, appending to previous log file.\n"
-                "\n"
-                );
-    }
-
-    fprintf(fp,
+    fprintf(fplog,
             "Log file opened on %s"
             "Host: %s  pid: %d  rank ID: %d  number of ranks:  %d\n",
-            timebuf, host, pid, cr->nodeid, cr->nnodes);
+            timebuf, host, pid, rankIndex, numRanks);
     try
     {
-        gmx::BinaryInformationSettings settings;
         settings.extendedInfo(true);
-        settings.copyright(!bAppendFiles);
-        gmx::printBinaryInformation(fp, gmx::getProgramContext(), settings);
+        gmx::printBinaryInformation(fplog, gmx::getProgramContext(), settings);
     }
     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
-    fprintf(fp, "\n");
+    fprintf(fplog, "\n");
 
-    fflush(fp);
+    fflush(fplog);
+}
 
-    *fplog = fp;
+void gmx_log_open(const char *lognm,
+                  const int   rankIndex,
+                  const int   numRanks,
+                  FILE     ** fplog)
+{
+    *fplog = gmx_fio_fopen(lognm, "w+");
+    if (*fplog == nullptr)
+    {
+        GMX_THROW(gmx::FileIOError("Could not open log file" + std::string(lognm)));
+    }
+    gmx::BinaryInformationSettings settings;
+    settings.copyright(true);
+    gmx_log_setup(settings, rankIndex, numRanks, *fplog);
+}
+
+void gmx_log_append(const int rankIndex,
+                    const int numRanks,
+                    FILE     *fplog)
+{
+    GMX_RELEASE_ASSERT(fplog != nullptr, "Log file must be already open");
+    fprintf(fplog,
+            "\n"
+            "\n"
+            "-----------------------------------------------------------\n"
+            "Restarting from checkpoint, appending to previous log file.\n"
+            "\n"
+            );
+    gmx::BinaryInformationSettings settings;
+    settings.copyright(false);
+    gmx_log_setup(settings, rankIndex, numRanks, fplog);
 }
 
 void gmx_log_close(FILE *fp)
index 9d9a2745f86ee26510fee973f6264845dcd4d2cc..bd8ed1dea08d6a5fcb06fbc5b710566c231b9f46 100644 (file)
 
 #include "gromacs/utility/basedefinitions.h"
 
-struct gmx_multisim_t;
-struct t_commrec;
+/*! \brief Open the log file for writing.
+ *
+ * \throws FileIOError when the log file cannot be opened.
+ * */
+void gmx_log_open(const char *lognm,
+                  int         rankIndex,
+                  int         numRanks,
+                  FILE      **fplog);
 
-/*! \brief Open the log file */
-void gmx_log_open(const char *fn, const t_commrec *cr,
-                  gmx_bool bAppendFiles, FILE** /*fplog*/);
+/*! \brief Use the open log file for appending.
+ *
+ * Does not throw.
+ */
+void gmx_log_append(int   rankIndex,
+                    int   numRanks,
+                    FILE *fplog);
 
-/*! \brief Close the log file */
+/*! \brief Close the log file.
+ *
+ * Does not throw.
+ */
 void gmx_log_close(FILE *fp);
 
 #endif
index ccd8365d1742f5d6899bff55916169f8d532669a..92d245e96af142b18f30d9b585a8fde8977643ff 100644 (file)
@@ -799,8 +799,7 @@ int Mdrunner::mdrunner()
 
     if (SIMMASTER(cr) && continuationOptions.appendFiles)
     {
-        gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
-                     continuationOptions.appendFiles, &fplog);
+        gmx_log_append(cr->nodeid, cr->nnodes, fplog);
         logOwner = buildLogger(fplog, nullptr);
         mdlog    = logOwner.logger();
     }
index 70238226404538da7e6620cb009fd203528d0b86..67d3aa357fc37a4d03ad06b18d09f7e60bc465d9 100644 (file)
@@ -490,8 +490,7 @@ int Mdrunner::mainFunction(int argc, char *argv[])
        there instead.  */
     if (MASTER(cr) && !continuationOptions.appendFiles)
     {
-        gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr,
-                     continuationOptions.appendFiles, &fplog);
+        gmx_log_open(ftp2fn(efLOG, nfile, fnm), cr->nodeid, cr->nnodes, &fplog);
     }
     else
     {