Fix trjconv when running without a structure file
authorPaul Bauer <paul.bauer.q@gmail.com>
Wed, 22 Aug 2018 09:58:58 +0000 (11:58 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 22 Aug 2018 16:18:43 +0000 (18:18 +0200)
When running trjconv without a reference structure file,
the topology information would be left uninitialized and
cause subsequent crashes during the deallocation.

Added tests to cover this behaviour.

Fixes #2619

Change-Id: I501c560ee0f8afca00bc78fd6f81b6841b5a5b57

docs/release-notes/2018/2018.3.rst
src/gromacs/gmxana/gmx_trjconv.cpp
src/gromacs/gmxana/tests/gmx_trjconv.cpp

index c8126a3da033ea714a4e213ae266dd3df2ea17ad..412d6539e8e416374d624e91969c755cf9c5a2f6 100644 (file)
@@ -111,6 +111,15 @@ the option. Now it only outputs values that are requested to be processed.
 
 :issue:`2565`
 
+Fix trjconv when not providing structure file
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+trjconv would fail with a segmentation violation when running without any structure
+file due to incomplete initialization of the topology data structure. This fix adds
+the missing checks that prevents the failure.
+
+:issue:`2619`
+
 Fixes to improve portability
 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^
 
index 3bd115e7fd982e8d852f30facc811213a2ec8c2a..4556135d889fe30fbd1b1804fb239a443bffab12 100644 (file)
@@ -870,7 +870,7 @@ int gmx_trjconv(int argc, char *argv[])
     real             *w_rls = nullptr;
     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
 #define SKIP 10
-    t_topology        top;
+    t_topology       *top   = nullptr;
     gmx_mtop_t       *mtop  = nullptr;
     gmx_conect        gc    = nullptr;
     int               ePBC  = -1;
@@ -1093,13 +1093,14 @@ int gmx_trjconv(int argc, char *argv[])
 
         if (bTPS)
         {
-            read_tps_conf(top_file, &top, &ePBC, &xp, nullptr, top_box,
+            snew(top, 1);
+            read_tps_conf(top_file, top, &ePBC, &xp, nullptr, top_box,
                           bReset || bPBCcomRes);
-            std::strncpy(top_title, *top.name, 255);
+            std::strncpy(top_title, *top->name, 255);
             top_title[255] = '\0';
-            atoms          = &top.atoms;
+            atoms          = &top->atoms;
 
-            if (0 == top.mols.nr && (bCluster || bPBCcomMol))
+            if (0 == top->mols.nr && (bCluster || bPBCcomMol))
             {
                 gmx_fatal(FARGS, "Option -pbc %s requires a .tpr file for the -s option", pbc_opt[pbc_enum]);
             }
@@ -1122,11 +1123,11 @@ int gmx_trjconv(int argc, char *argv[])
 
             if (bCONECT)
             {
-                gc = gmx_conect_generate(&top);
+                gc = gmx_conect_generate(top);
             }
             if (bRmPBC)
             {
-                gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
+                gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
             }
         }
 
@@ -1219,7 +1220,7 @@ int gmx_trjconv(int argc, char *argv[])
                store original location (to put structure back) */
             if (bRmPBC)
             {
-                gmx_rmpbc(gpbc, top.atoms.nr, top_box, xp);
+                gmx_rmpbc(gpbc, top->atoms.nr, top_box, xp);
             }
             copy_rvec(xp[index[0]], x_shift);
             reset_x_ndim(nfitdim, ifit, ind_fit, atoms->nr, nullptr, xp, w_rls);
@@ -1530,7 +1531,7 @@ int gmx_trjconv(int argc, char *argv[])
                 }
                 else if (bCluster)
                 {
-                    calc_pbc_cluster(ecenter, ifit, &top, ePBC, fr.x, ind_fit, fr.box);
+                    calc_pbc_cluster(ecenter, ifit, top, ePBC, fr.x, ind_fit, fr.box);
                 }
 
                 if (bPFit)
@@ -1709,7 +1710,7 @@ int gmx_trjconv(int argc, char *argv[])
                         if (bPBCcomMol)
                         {
                             put_molecule_com_in_box(unitcell_enum, ecenter,
-                                                    &top.mols,
+                                                    &top->mols,
                                                     natoms, atoms->atom, ePBC, fr.box, fr.x);
                         }
                         /* Copy the input trxframe struct to the output trxframe struct */
@@ -1978,7 +1979,11 @@ int gmx_trjconv(int argc, char *argv[])
     }
 
     sfree(mtop);
-    done_top(&top);
+    if (top)
+    {
+        done_top(top);
+        sfree(top);
+    }
     sfree(xp);
     sfree(xmem);
     sfree(vmem);
index 795fcdb3ad0fd393097b7846d105b84a41908eb3..49a8a1dad696826997646a6f59a98723a7e9a958 100644 (file)
@@ -104,4 +104,34 @@ INSTANTIATE_TEST_CASE_P(NoFatalErrorWhenWritingFrom,
                         TrjconvWithIndexGroupSubset,
                             ::testing::ValuesIn(trajectoryFileNames));
 
+class TrjconvWithoutTopologyFile : public gmx::test::CommandLineTestBase,
+                                   public ::testing::WithParamInterface<const char *>
+{
+    public:
+        void runTest(const char *fileName)
+        {
+            auto &cmdline = commandLine();
+
+            setInputFile("-f", fileName);
+            setInputFile("-n", "spc2.ndx");
+            setOutputFile("-o", "spc-traj.trr", gmx::test::NoTextMatch());
+
+            gmx::test::StdioTestHelper stdioHelper(&fileManager());
+            stdioHelper.redirectStringToStdin("SecondWaterMolecule\n");
+
+            /* As mentioned above, the tests don't check much besides
+             * that trjconv does not crash.
+             */
+            ASSERT_EQ(0, gmx_trjconv(cmdline.argc(), cmdline.argv()));
+        }
+};
+
+TEST_P(TrjconvWithoutTopologyFile, WithDifferentInputFormats)
+{
+    runTest(GetParam());
+}
+
+INSTANTIATE_TEST_CASE_P(NoFatalErrorWhenWritingFrom,
+                        TrjconvWithoutTopologyFile,
+                            ::testing::ValuesIn(trajectoryFileNames));
 } // namespace