Remove ineffective code for writing more TER fields to pdb files
authorMark Abraham <mark.j.abraham@gmail.com>
Sun, 10 Feb 2019 13:28:49 +0000 (14:28 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 7 Mar 2019 12:13:12 +0000 (13:13 +0100)
A weakly typed string containing a residue type name was passed to a
function expecting a residue name. This has led to fewer TER records
being output than expected by the author of ff9a9692e1f4. The intent
then was presumably to refer to the residue type of the previous
chain, but in fact the residue type of that residue type was used,
which is always "Other". Thus no extra TER fields were ever written.

The intended feature for PDB writing is not clearly vital, and has
never worked, so is now removed. These days, it would only be active
if pdb2gmx wrote .pdb output, since in no other case do we build chain
numbers.

Change-Id: I4e5ed09171772e17fd0fc11aff239bf7aaa966eb

src/gromacs/fileio/confio.cpp
src/gromacs/fileio/pdbio.cpp
src/gromacs/fileio/pdbio.h
src/gromacs/fileio/trxio.cpp
src/gromacs/gmxana/gmx_chi.cpp
src/gromacs/gmxana/gmx_confrms.cpp
src/gromacs/gmxana/gmx_do_dssp.cpp
src/gromacs/gmxana/gmx_trjconv.cpp
src/gromacs/gmxpreprocess/editconf.cpp

index 1c66a2495662be1a9ae902eea4bae01a4fba1101..9def609e43eb390aa181e3e78f54003be6bf21ea 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,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -101,7 +101,7 @@ void write_sto_conf_indexed(const char *outfile, const char *title,
         case efENT:
         case efPQR:
             out = gmx_fio_fopen(outfile, "w");
-            write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, nullptr, TRUE, ftp == efPQR);
+            write_pdbfile_indexed(out, title, atoms, x, ePBC, box, ' ', -1, nindex, index, nullptr, ftp == efPQR);
             gmx_fio_fclose(out);
             break;
         case efESP:
@@ -151,7 +151,7 @@ void write_sto_conf(const char *outfile, const char *title, const t_atoms *atoms
         case efBRK:
         case efENT:
             out = gmx_fio_fopen(outfile, "w");
-            write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, nullptr, TRUE);
+            write_pdbfile(out, title, atoms, x, ePBC, box, ' ', -1, nullptr);
             gmx_fio_fclose(out);
             break;
         case efESP:
index b6882696e6c618c97c0b7a0336b86799fd8137e1..14efea694f678d1fe27bc1bea435b48b0340d3fc 100644 (file)
@@ -305,7 +305,7 @@ void write_pdbfile_indexed(FILE *out, const char *title,
                            const t_atoms *atoms, const rvec x[],
                            int ePBC, const matrix box, char chainid,
                            int model_nr, int nindex, const int index[],
-                           gmx_conect conect, gmx_bool bTerSepChains,
+                           gmx_conect conect,
                            bool usePqrFormat)
 {
     gmx_conect_t   *gc = static_cast<gmx_conect_t *>(conect);
@@ -339,34 +339,13 @@ void write_pdbfile_indexed(FILE *out, const char *title,
 
     fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
 
-    int         lastchainnum      = -1;
-    std::string prevRestype;
-    std::string lastRestype;
-
     ResidueType rt;
     for (int ii = 0; ii < nindex; ii++)
     {
-        int i             = index[ii];
-        int resind        = atoms->atom[i].resind;
-        int chainnum      = atoms->resinfo[resind].chainnum;
-        lastRestype = prevRestype;
-        prevRestype = rt.typeNameForIndexedResidue(*atoms->resinfo[resind].name);
-
-        /* Add a TER record if we changed chain, and if either the previous or this chain is protein/DNA/RNA. */
-        if (bTerSepChains && ii > 0 && chainnum != lastchainnum)
-        {
-            /* Only add TER if the previous chain contained protein/DNA/RNA. */
-            if (rt.namedResidueHasType(lastRestype, "Protein") ||
-                rt.namedResidueHasType(lastRestype, "DNA") ||
-                rt.namedResidueHasType(lastRestype, "RNA"))
-            {
-                fprintf(out, "TER\n");
-            }
-            lastchainnum    = chainnum;
-        }
-
-        std::string resnm = *atoms->resinfo[resind].name;
-        std::string nm    = *atoms->atomname[i];
+        int         i             = index[ii];
+        int         resind        = atoms->atom[i].resind;
+        std::string resnm         = *atoms->resinfo[resind].name;
+        std::string nm            = *atoms->atomname[i];
 
         /* rename HG12 to 2HG1, etc. */
         nm    = xlate_atomname_gmx2pdb(nm);
@@ -462,7 +441,7 @@ void write_pdbfile_indexed(FILE *out, const char *title,
 }
 
 void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms, const rvec x[],
-                   int ePBC, const matrix box, char chainid, int model_nr, gmx_conect conect, gmx_bool bTerSepChains)
+                   int ePBC, const matrix box, char chainid, int model_nr, gmx_conect conect)
 {
     int i, *index;
 
@@ -472,7 +451,7 @@ void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms, const rve
         index[i] = i;
     }
     write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
-                          atoms->nr, index, conect, bTerSepChains, false);
+                          atoms->nr, index, conect, false);
     sfree(index);
 }
 
index 46a3ca40d7d8f21389b996311d1aa3abc23ca3cd..98badfc6a790b718b898e655bbd3c6c228233178 100644 (file)
@@ -91,13 +91,13 @@ void gmx_write_pdb_box(FILE *out, int ePBC, const matrix box);
 void write_pdbfile_indexed(FILE *out, const char *title, const t_atoms *atoms,
                            const rvec x[], int ePBC, const matrix box, char chain,
                            int model_nr, int nindex, const int index[],
-                           gmx_conect conect, gmx_bool bTerSepChains,
+                           gmx_conect conect,
                            bool usePqrFormat);
 /* REALLY low level */
 
 void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms,
                    const rvec x[], int ePBC, const matrix box, char chain,
-                   int model_nr, gmx_conect conect, gmx_bool bTerSepChains);
+                   int model_nr, gmx_conect conect);
 /* Low level pdb file writing routine.
  *
  *          ONLY FOR SPECIAL PURPOSES,
index f9a33fefac166d668280adb6786602c12c3328d8..356358b5378a7d912f77143ce129028af11132eb 100644 (file)
@@ -428,7 +428,7 @@ int write_trxframe_indexed(t_trxstatus *status, const t_trxframe *fr, int nind,
             else
             {
                 write_pdbfile_indexed(gmx_fio_getfp(status->fio), title, fr->atoms,
-                                      fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, TRUE, FALSE);
+                                      fr->x, -1, fr->box, ' ', fr->step, nind, ind, gc, FALSE);
             }
             break;
         case efG96:
@@ -592,7 +592,7 @@ int write_trxframe(t_trxstatus *status, t_trxframe *fr, gmx_conect gc)
             {
                 write_pdbfile(gmx_fio_getfp(status->fio), title,
                               fr->atoms, fr->x, fr->bPBC ? fr->ePBC : -1, fr->box,
-                              ' ', fr->step, gc, TRUE);
+                              ' ', fr->step, gc);
             }
             break;
         case efG96:
index 4931615de89e916d1218abdf4dcb7dbda95d8549..bda1e768f9f927259e4717d3564959e363b8cafa 100644 (file)
@@ -1167,7 +1167,7 @@ static void order_params(FILE *log,
         fprintf(fp, "REMARK generated by g_chi\n");
         fprintf(fp, "REMARK "
                 "B-factor field contains negative of dihedral order parameters\n");
-        write_pdbfile(fp, nullptr, atoms, x, ePBC, box, ' ', 0, nullptr, TRUE);
+        write_pdbfile(fp, nullptr, atoms, x, ePBC, box, ' ', 0, nullptr);
         x0 = y0 = z0 = 1000.0;
         for (i = 0; (i < atoms->nr); i++)
         {
index 4b43f2a77537b0f5a50de9110b8d911d5a01b3c7..0ad9d7c3faffe9eef39d7a4d3fddd483bc716cd5 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,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -801,9 +801,9 @@ int gmx_confrms(int argc, char *argv[])
             fp = gmx_ffopen(outfile, "w");
             if (!bOne)
             {
-                write_pdbfile(fp, *top1->name, atoms1, x1, ePBC1, box1, ' ', 1, nullptr, TRUE);
+                write_pdbfile(fp, *top1->name, atoms1, x1, ePBC1, box1, ' ', 1, nullptr);
             }
-            write_pdbfile(fp, *top2->name, atoms2, x2, ePBC2, box2, ' ', bOne ? -1 : 2, nullptr, TRUE);
+            write_pdbfile(fp, *top2->name, atoms2, x2, ePBC2, box2, ' ', bOne ? -1 : 2, nullptr);
             gmx_ffclose(fp);
             break;
         case efGRO:
index b70cba7deb4a64b6ffca9ebbfbc34fca234c9637..a466b222c2989d9a63c0a4cb7227eb53c1d971dc 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) 2012,2013,2014,2015,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, 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.
@@ -721,7 +721,7 @@ int gmx_do_dssp(int argc, char *argv[])
         }
         gmx_rmpbc(gpbc, natoms, box, x);
         tapein = gmx_ffopen(pdbfile, "w");
-        write_pdbfile_indexed(tapein, nullptr, atoms, x, ePBC, box, ' ', -1, gnx, index, nullptr, TRUE, FALSE);
+        write_pdbfile_indexed(tapein, nullptr, atoms, x, ePBC, box, ' ', -1, gnx, index, nullptr, FALSE);
         gmx_ffclose(tapein);
         /* strip_dssp returns the number of lines found in the dssp file, i.e.
          * the number of residues plus the separator lines */
index 2057c011ce9aced6dda03fef3a2e1260e5cfc9ef..d18edaa98c9320e2383af9d729cdfd6f8d2a9b61 100644 (file)
@@ -1885,7 +1885,7 @@ int gmx_trjconv(int argc, char *argv[])
                                             model_nr++;
                                         }
                                         write_pdbfile(out, title.c_str(), &useatoms, frout.x,
-                                                      frout.ePBC, frout.box, ' ', model_nr, gc, TRUE);
+                                                      frout.ePBC, frout.box, ' ', model_nr, gc);
                                         break;
                                     case efG96:
                                         const char *outputTitle = "";
index 59715d09f9b19694a926ca4a05cf64180937901b..6d0c73bab9f4443cacf8c28e8149b489f6acde2c 100644 (file)
@@ -1266,7 +1266,7 @@ int gmx_editconf(int argc, char *argv[])
         if (outftp == efPDB)
         {
             out = gmx_ffopen(outfile, "w");
-            write_pdbfile_indexed(out, *top_tmp.name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, TRUE, FALSE);
+            write_pdbfile_indexed(out, *top_tmp.name, &atoms, x, ePBC, box, ' ', 1, isize, index, conect, FALSE);
             gmx_ffclose(out);
         }
         else
@@ -1322,7 +1322,7 @@ int gmx_editconf(int argc, char *argv[])
                 index[i] = i;
             }
             write_pdbfile_indexed(out, *top_tmp.name, &atoms, x, ePBC, box, ' ', -1, atoms.nr, index, conect,
-                                  TRUE, outftp == efPQR);
+                                  outftp == efPQR);
             sfree(index);
             if (bLegend)
             {