Support concurrent t_topology and gmx_mtop_t use
authorTeemu Murtola <teemu.murtola@gmail.com>
Mon, 20 Jun 2016 17:55:59 +0000 (20:55 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Tue, 11 Oct 2016 12:36:54 +0000 (14:36 +0200)
Make gmx_mtop_t_to_t_topology() take a boolean flag that specifies
whether to free memory from mtop or not.  This makes it possible to
convert an mtop to a t_topology and then use both in subsequent code.
This in turn makes a piecewise conversion from t_topology to gmx_mtop_t
possible.

Using these two concurrently currently always leads to memory leaks, but
these should mostly be temporary during the conversion, and can be
sorted out separately if necessary.

Part of #1862.

Change-Id: Iccd6039c2226d1dc617963878e5640bb53f09ad6

src/gromacs/fileio/confio.cpp
src/gromacs/fileio/tpxio.cpp
src/gromacs/gmxana/gmx_nmeig.cpp
src/gromacs/tools/check.cpp
src/gromacs/tools/convert_tpr.cpp
src/gromacs/tools/dump.cpp
src/gromacs/topology/mtop_util.cpp
src/gromacs/topology/mtop_util.h
src/gromacs/topology/topology.h

index 66cf9aadf82552fc6827c9895794280d67adbfd4..7735957498427c434613cdfc4ffe885eeab2c416 100644 (file)
@@ -360,27 +360,6 @@ static void read_stx_conf(const char *infile,
     }
 }
 
-static void done_gmx_groups_t(gmx_groups_t *g)
-{
-    int i;
-
-    for (i = 0; (i < egcNR); i++)
-    {
-        if (NULL != g->grps[i].nm_ind)
-        {
-            sfree(g->grps[i].nm_ind);
-            g->grps[i].nm_ind = NULL;
-        }
-        if (NULL != g->grpnr[i])
-        {
-            sfree(g->grpnr[i]);
-            g->grpnr[i] = NULL;
-        }
-    }
-    /* The contents of this array is in symtab, don't free it here */
-    sfree(g->grpname);
-}
-
 static void readConfAndAtoms(const char *infile,
                              t_symtab *symtab, char ***name, t_atoms *atoms,
                              int *ePBC,
@@ -471,9 +450,7 @@ gmx_bool read_tps_conf(const char *infile, t_topology *top, int *ePBC,
     snew(mtop, 1);
     readConfAndTopology(infile, &haveTopology, mtop, ePBC, x, v, box);
 
-    *top = gmx_mtop_t_to_t_topology(mtop);
-    /* In this case we need to throw away the group data too */
-    done_gmx_groups_t(&mtop->groups);
+    *top = gmx_mtop_t_to_t_topology(mtop, true);
     sfree(mtop);
 
     tpx_make_chain_identifiers(&top->atoms, &top->mols);
index b28ee11a7c1ffc531af8d3f13e30a82f5cf91b56..4cf763262654896a58f90ad9070a3fad7a0db211 100644 (file)
@@ -3357,7 +3357,7 @@ int read_tpx_top(const char *fn,
 
     ePBC = read_tpx(fn, ir, box, natoms, x, v, &mtop);
 
-    *top = gmx_mtop_t_to_t_topology(&mtop);
+    *top = gmx_mtop_t_to_t_topology(&mtop, true);
 
     return ePBC;
 }
index 0c363ed7c11c0c5c44c79a032e6fe3cd9b991ede..61f3bd68895b2b9182be697d2faebf07526b7c62 100644 (file)
@@ -402,7 +402,7 @@ int gmx_nmeig(int argc, char *argv[])
     }
     std::vector<size_t> atom_index = get_atom_index(&mtop);
 
-    top = gmx_mtop_t_to_t_topology(&mtop);
+    top = gmx_mtop_t_to_t_topology(&mtop, true);
 
     bM       = TRUE;
     int ndim = DIM*atom_index.size();
index 2238ae7d21dceadce34e8bbb056032b3111b1a83..d2fbd6998484feb45cb4c4d1e4616696946abd12 100644 (file)
@@ -110,8 +110,8 @@ static void comp_tpx(const char *fn1, const char *fn2,
          * We should implement direct mtop comparison,
          * but it might be useful to keep t_topology comparison as an option.
          */
-        top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
-        top[1] = gmx_mtop_t_to_t_topology(&mtop[1]);
+        top[0] = gmx_mtop_t_to_t_topology(&mtop[0], false);
+        top[1] = gmx_mtop_t_to_t_topology(&mtop[1], false);
         cmp_top(stdout, &top[0], &top[1], ftol, abstol);
         cmp_groups(stdout, &mtop[0].groups, &mtop[1].groups,
                    mtop[0].natoms, mtop[1].natoms);
@@ -133,7 +133,7 @@ static void comp_tpx(const char *fn1, const char *fn2,
              * We should implement direct mtop comparison,
              * but it might be useful to keep t_topology comparison as an option.
              */
-            top[0] = gmx_mtop_t_to_t_topology(&mtop[0]);
+            top[0] = gmx_mtop_t_to_t_topology(&mtop[0], true);
             cmp_top(stdout, &top[0], NULL, ftol, abstol);
         }
     }
index da7bff85de83ece3a2853e21a76f076b5dd5128b..21303d4339e861d03d90f75924cf42b02dcfefdd 100644 (file)
@@ -267,7 +267,7 @@ static void reduce_topology_x(int gnx, int index[],
     int         *invindex;
     int          i;
 
-    top      = gmx_mtop_t_to_t_topology(mtop);
+    top      = gmx_mtop_t_to_t_topology(mtop, false);
     bKeep    = bKeepIt(gnx, top.atoms.nr, index);
     invindex = invind(gnx, top.atoms.nr, index);
 
index c6fd2941d82cbcf4e416b76d408370da09a5bcdc..07d1c0fc8d7bb527f58009ffbd2deb0b2c87cad6 100644 (file)
@@ -105,7 +105,7 @@ static void list_tpx(const char *fn, gmx_bool bShowNumbers, const char *mdpfn,
     {
         if (bSysTop)
         {
-            top = gmx_mtop_t_to_t_topology(&mtop);
+            top = gmx_mtop_t_to_t_topology(&mtop, false);
         }
 
         if (available(stdout, &tpx, 0, fn))
index 0f140660761e9be2188045d02c0b3987c8e8edcd..f83312fb27d1484ab68883cfe68ff5709e2e0c0b 100644 (file)
@@ -956,7 +956,28 @@ gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,
     return top;
 }
 
-t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
+static void done_gmx_groups_t(gmx_groups_t *g)
+{
+    int i;
+
+    for (i = 0; (i < egcNR); i++)
+    {
+        if (NULL != g->grps[i].nm_ind)
+        {
+            sfree(g->grps[i].nm_ind);
+            g->grps[i].nm_ind = NULL;
+        }
+        if (NULL != g->grpnr[i])
+        {
+            sfree(g->grpnr[i]);
+            g->grpnr[i] = NULL;
+        }
+    }
+    /* The contents of this array is in symtab, don't free it here */
+    sfree(g->grpname);
+}
+
+t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop)
 {
     int            mt, mb;
     gmx_localtop_t ltop;
@@ -975,23 +996,23 @@ t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop)
     top.bIntermolecularInteractions = mtop->bIntermolecularInteractions;
     top.symtab                      = mtop->symtab;
 
-    /* We only need to free the moltype and molblock data,
-     * all other pointers have been copied to top.
-     *
-     * Well, except for the group data, but we can't free those, because they
-     * are used somewhere even after a call to this function.
-     */
-    for (mt = 0; mt < mtop->nmoltype; mt++)
+    if (freeMTop)
     {
-        done_moltype(&mtop->moltype[mt]);
-    }
-    sfree(mtop->moltype);
+        // Free pointers that have not been copied to top.
+        for (mt = 0; mt < mtop->nmoltype; mt++)
+        {
+            done_moltype(&mtop->moltype[mt]);
+        }
+        sfree(mtop->moltype);
 
-    for (mb = 0; mb < mtop->nmolblock; mb++)
-    {
-        done_molblock(&mtop->molblock[mb]);
+        for (mb = 0; mb < mtop->nmolblock; mb++)
+        {
+            done_molblock(&mtop->molblock[mb]);
+        }
+        sfree(mtop->molblock);
+
+        done_gmx_groups_t(&mtop->groups);
     }
-    sfree(mtop->molblock);
 
     return top;
 }
index 622f4cd2c88b7f6adaf8037e163d0a255f6dcd19..8196d71ee55f056b896061b3905a814ac10b8499 100644 (file)
@@ -200,7 +200,7 @@ gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
 
 
 /* Generate a 'local' topology for the whole system.
- * When feeEnergyInteractionsAtEnd == true, the free energy interactions will
+ * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
  * be sorted to the end.
  */
 gmx_localtop_t *
@@ -208,10 +208,14 @@ gmx_mtop_generate_local_top(const gmx_mtop_t *mtop, bool freeEnergyInteractionsA
 
 
 /* Converts a gmx_mtop_t struct to t_topology.
- * All memory relating only to mtop will be freed.
+ *
+ * If freeMTop == true, memory related to mtop will be freed so that done_top()
+ * on the result value will free all memory.
+ * If freeMTop == false, mtop and the return value will share some of their
+ * memory, and there is currently no way to consistently free all the memory.
  */
 t_topology
-gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop);
+gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop);
 
 /*! \brief Get vector of atoms indices from topology
  *
index 3599ac893a90921dfc01f231825b2efbd964aac7..6b36d1f66ba8a6f5cfcb481b75ba813b152237e5 100644 (file)
@@ -150,9 +150,6 @@ void done_molblock(gmx_molblock_t *molb);
 void done_mtop(gmx_mtop_t *mtop, gmx_bool bDoneSymtab);
 void done_top(t_topology *top);
 
-t_atoms *mtop2atoms(gmx_mtop_t *mtop);
-/* generate a t_atoms struct for the system from gmx_mtop_t */
-
 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
              gmx_bool bShowNumbers);
 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top, gmx_bool bShowNumbers);