Fixed gmx rmsf -q -oq.
authorMark Abraham <mark.j.abraham@gmail.com>
Sun, 22 Jan 2017 13:56:19 +0000 (14:56 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Sun, 22 Jan 2017 14:07:37 +0000 (15:07 +0100)
xref was used instead of pdbx in that case, which led to a PDB file
containing B-factors to have the coordinates based on those from the
-s file, rather than -q file. The tool description (and the other
logic) suggests that using -q should be the behaviour, in particular
because the contents of xref were then unused. gmx rmsf -oq was
otherwise fine.

Also fixed where xref was being abused as a temporary vector, which
asks for more such trouble, so improved that (without changing
functionality).

Noted TODO that npdbatoms might not match top.atoms.nr leading to
possible segfaults.

Change-Id: Id2c14fa5ca93c2a8de87f9bfa144a6db38fcc97f

src/gromacs/gmxana/gmx_rmsf.c

index afca0061548c0f4ff860f1c8959b880dd47e6de4..7f90acadbdc7b987328c70fde4650524fced3d20 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, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2017, 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.
@@ -323,6 +323,7 @@ int gmx_rmsf(int argc, char *argv[])
         /* Read coordinates twice */
         read_stx_conf(opt2fn("-q", NFILE, fnm), title, pdbatoms, pdbx, NULL, NULL, pdbbox);
         read_stx_conf(opt2fn("-q", NFILE, fnm), title, refatoms, pdbx, NULL, NULL, pdbbox);
+        /* TODO Should this assert that npdbatoms == top.atoms.nr? */
     }
     else
     {
@@ -548,24 +549,26 @@ int gmx_rmsf(int argc, char *argv[])
         /* Write a .pdb file with B-factors and optionally anisou records */
         for (i = 0; i < isize; i++)
         {
-            rvec_inc(xref[index[i]], xcm);
+            rvec_inc(pdbx[index[i]], xcm);
         }
         write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
                                NULL, ePBC, pdbbox, isize, index);
     }
     if (opt2bSet("-ox", NFILE, fnm))
     {
-        /* Misuse xref as a temporary array */
+        rvec *bFactorX;
+        snew(bFactorX, top.atoms.nr);
         for (i = 0; i < isize; i++)
         {
             for (d = 0; d < DIM; d++)
             {
-                xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
+                bFactorX[index[i]][d] = xcm[d] + xav[i*DIM + d];
             }
         }
         /* Write a .pdb file with B-factors and optionally anisou records */
-        write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, NULL,
+        write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, bFactorX, NULL,
                                ePBC, pdbbox, isize, index);
+        sfree(bFactorX);
     }
     if (bAniso)
     {