Closed some .xvg files that were being left open
authorTeemu Murtola <teemu.murtola@gmail.com>
Thu, 13 Nov 2014 20:14:25 +0000 (22:14 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 8 Jan 2015 17:26:35 +0000 (18:26 +0100)
Closed donor properties file in gmx hbond after removing use of
function static file pointer.

Improved naming and scoping of some file pointers in
analyze_clusters().

Added done_ed() as a place to close the essential dynamics output
file. Later, destructor content for ED should go here (or preferably
RAII).

Removed pr_energy() from gmx_cluster since it was never called.

Change-Id: I4c288e924eb81e9e0f87b8e59ee64dacc0708567

src/gromacs/essentialdynamics/edsam.c
src/gromacs/essentialdynamics/edsam.h
src/gromacs/gmxana/gmx_cluster.c
src/gromacs/gmxana/gmx_hbond.c
src/programs/mdrun/runner.cpp

index 558eb13dad1be2f3acb3a83e655fc5428761b7ff..35855ffde15a2b9103dd24f1f297baa750d112b1 100644 (file)
@@ -3184,3 +3184,16 @@ void do_edsam(t_inputrec     *ir,
 
     ed->bFirst = FALSE;
 }
+
+void done_ed(gmx_edsam_t *ed)
+{
+    if (*ed)
+    {
+        /* ed->edo is opened sometimes with xvgropen, sometimes with
+         * gmx_fio_fopen, so we use the least common denominator for
+         * closing. */
+        gmx_fio_fclose((*ed)->edo);
+    }
+
+    /* TODO deallocate ed and set pointer to NULL */
+}
index 20cbaf442b0e3e2ed59abf8a06b282466389c6d1..e6ab9810cdc4dbe0c1e63d176d3da5b78e8874e1 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.
@@ -129,6 +129,11 @@ void dd_make_local_ed_indices(gmx_domdec_t *dd, gmx_edsam_t ed);
 void do_flood(t_commrec *cr, t_inputrec *ir, rvec x[], rvec force[], gmx_edsam_t ed,
               matrix box, gmx_int64_t step, gmx_bool bNS);
 
+/*! \brief Clean up
+ *
+ * \param ed                The essential dynamics data
+ */
+void done_ed(gmx_edsam_t *ed);
 
 #ifdef __cplusplus
 }
index 3d0b0bef8daf84cdf544fa2439f3bd55dd1420d0..034e10a2b05b03580b9ebb4575efdb2c44807c54 100644 (file)
@@ -136,11 +136,6 @@ typedef struct {
     int *nb;
 } t_nnb;
 
-void pr_energy(FILE *fp, real e)
-{
-    fprintf(fp, "Energy: %8.4f\n", e);
-}
-
 void cp_index(int nn, int from[], int to[])
 {
     int i;
@@ -1018,7 +1013,7 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
                              gmx_bool bFit, FILE *log, t_rgb rlo, t_rgb rhi,
                              const output_env_t oenv)
 {
-    FILE        *fp = NULL;
+    FILE        *size_fp = NULL;
     char         buf[STRLEN], buf1[40], buf2[40], buf3[40], *trxsfn;
     t_trxstatus *trxout  = NULL;
     t_trxstatus *trxsout = NULL;
@@ -1093,7 +1088,7 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
 
     if (clustidfn)
     {
-        fp = xvgropen(clustidfn, "Clusters", output_env_get_xvgr_tlabel(oenv), "Cluster #", oenv);
+        FILE *fp = xvgropen(clustidfn, "Clusters", output_env_get_xvgr_tlabel(oenv), "Cluster #", oenv);
         if (output_env_get_print_xvgr_codes(oenv))
         {
             fprintf(fp, "@    s0 symbol 2\n");
@@ -1108,11 +1103,10 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
     }
     if (sizefn)
     {
-        /* FIXME: This file is never closed. */
-        fp = xvgropen(sizefn, "Cluster Sizes", "Cluster #", "# Structures", oenv);
+        size_fp = xvgropen(sizefn, "Cluster Sizes", "Cluster #", "# Structures", oenv);
         if (output_env_get_print_xvgr_codes(oenv))
         {
-            fprintf(fp, "@g%d type %s\n", 0, "bar");
+            fprintf(size_fp, "@g%d type %s\n", 0, "bar");
         }
     }
     snew(structure, nf);
@@ -1161,7 +1155,7 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
         }
         if (sizefn)
         {
-            fprintf(fp, "%8d %8d\n", cl, nstr);
+            fprintf(size_fp, "%8d %8d\n", cl, nstr);
         }
         clrmsd  = 0;
         midstr  = 0;
@@ -1307,6 +1301,11 @@ static void analyze_clusters(int nf, t_clusters *clust, real **rmsd,
     {
         sfree(trxsfn);
     }
+
+    if (size_fp)
+    {
+        xvgrclose(size_fp);
+    }
 }
 
 static void convert_mat(t_matrix *mat, t_mat *rms)
index 2f9d4578f8d3c17699892fea7936929122c81722..8bf0e44054fbe92cdb36f5ef50e682a4efc4ece5 100644 (file)
@@ -3328,22 +3328,31 @@ static void init_hbframe(t_hbdata *hb, int nframes, real t)
     /*set_hb(hb->hbmap[i][j]->h[m],nframes-hb->hbmap[i][j]->n0,HB_NO);*/
 }
 
-static void analyse_donor_props(const char *fn, t_hbdata *hb, int nframes, real t,
-                                const output_env_t oenv)
+static FILE *open_donor_properties_file(const char        *fn,
+                                        t_hbdata          *hb,
+                                        const output_env_t oenv)
 {
-    static FILE *fp    = NULL;
-    const char  *leg[] = { "Nbound", "Nfree" };
-    int          i, j, k, nbound, nb, nhtot;
+    FILE       *fp    = NULL;
+    const char *leg[] = { "Nbound", "Nfree" };
 
-    if (!fn)
+    if (!fn || !hb)
     {
-        return;
+        return NULL;
     }
-    if (!fp)
+
+    fp = xvgropen(fn, "Donor properties", output_env_get_xvgr_tlabel(oenv), "Number", oenv);
+    xvgr_legend(fp, asize(leg), leg, oenv);
+
+    return fp;
+}
+
+static void analyse_donor_properties(FILE *fp, t_hbdata *hb, int nframes, real t)
+{
+    int i, j, k, nbound, nb, nhtot;
+
+    if (!fp || !hb)
     {
-        // TODO: This file is never closed...
-        fp = xvgropen(fn, "Donor properties", output_env_get_xvgr_tlabel(oenv), "Number", oenv);
-        xvgr_legend(fp, asize(leg), leg, oenv);
+        return;
     }
     nbound = 0;
     nhtot  = 0;
@@ -3689,7 +3698,7 @@ int gmx_hbond(int argc, char *argv[])
     int                   grp, nabin, nrbin, bin, resdist, ihb;
     char                **leg;
     t_hbdata             *hb, *hbptr;
-    FILE                 *fp, *fpins = NULL, *fpnhb = NULL;
+    FILE                 *fp, *fpins = NULL, *fpnhb = NULL, *donor_properties = NULL;
     t_gridcell         ***grid;
     t_ncell              *icell, *jcell, *kcell;
     ivec                  ngrid;
@@ -3948,6 +3957,8 @@ int gmx_hbond(int argc, char *argv[])
     /*if (bSelected)
        snew(donors[gr0D], dons[gr0D].nrd);*/
 
+    donor_properties = open_donor_properties_file(opt2fn_null("-don", NFILE, fnm), hb, oenv);
+
     if (bHBmap)
     {
         printf("Making hbmap structure...");
@@ -4376,10 +4387,7 @@ int gmx_hbond(int argc, char *argv[])
 #pragma omp barrier
 #pragma omp single
                 {
-                    if (hb != NULL)
-                    {
-                        analyse_donor_props(opt2fn_null("-don", NFILE, fnm), hb, k, t, oenv);
-                    }
+                    analyse_donor_properties(donor_properties, hb, k, t);
                 }
 
 #pragma omp single
@@ -4449,6 +4457,12 @@ int gmx_hbond(int argc, char *argv[])
     free_grid(ngrid, &grid);
 
     close_trj(status);
+
+    if (donor_properties)
+    {
+        xvgrclose(donor_properties);
+    }
+
     if (fpnhb)
     {
         xvgrclose(fpnhb);
index 9093d6c05787392ef6994d98e494034668c432e8..229c25e3e1e14852a4407d6f770d625904654609 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) 2011,2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2011,2012,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.
@@ -1774,6 +1774,8 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
 
     rc = (int)gmx_get_stop_condition();
 
+    done_ed(&ed);
+
 #ifdef GMX_THREAD_MPI
     /* we need to join all threads. The sub-threads join when they
        exit this function, but the master thread needs to be told to