Improve PDB I/O (keep occupancy, b-factor, element)
authorErik Lindahl <erik@kth.se>
Sat, 21 Jun 2014 21:57:34 +0000 (23:57 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 27 Jun 2014 12:12:50 +0000 (14:12 +0200)
Retain the extra PDB info such as occupancy, b-factors,
and element name columns with trjconv when providing a
PDB structure file. If you use a TPR file, we still
write the (given atomnumbers being present).
However, we deliberately don't store the pdbinfo in the
TPR file since that is a specification for a particular
experiment, rather than something that's valid in the
simulation.

Fixes #917, #1307.

Change-Id: I80cb8f5a250ba094fe81f32c58b4eb0298164053

src/gromacs/fileio/tpxio.c
src/gromacs/gmxana/gmx_trjconv.c

index f8db1c8a892dc35c2eb361d992f7e9c746c78355..62e777ff55f11c5a897c3159ea2d301b4a8005a5 100644 (file)
@@ -2381,10 +2381,52 @@ static void do_blocka(t_fileio *fio, t_blocka *block, gmx_bool bRead,
     gmx_fio_ndo_int(fio, block->a, block->nra);
 }
 
+/* This is a primitive routine to make it possible to translate atomic numbers
+ * to element names when reading TPR files, without making the Gromacs library
+ * directory a dependency on mdrun (which is the case if we need elements.dat).
+ */
+static const char *
+atomicnumber_to_element(int atomicnumber)
+{
+    const char * p;
+
+    /* This does not have to be complete, so we only include elements likely
+     * to occur in PDB files.
+     */
+    switch (atomicnumber)
+    {
+        case 1:  p = "H";  break;
+        case 5:  p = "B";  break;
+        case 6:  p = "C";  break;
+        case 7:  p = "N";  break;
+        case 8:  p = "O";  break;
+        case 9:  p = "F";  break;
+        case 11: p = "Na"; break;
+        case 12: p = "Mg"; break;
+        case 15: p = "P";  break;
+        case 16: p = "S";  break;
+        case 17: p = "Cl"; break;
+        case 18: p = "Ar"; break;
+        case 19: p = "K";  break;
+        case 20: p = "Ca"; break;
+        case 25: p = "Mn"; break;
+        case 26: p = "Fe"; break;
+        case 28: p = "Ni"; break;
+        case 29: p = "Cu"; break;
+        case 30: p = "Zn"; break;
+        case 35: p = "Br"; break;
+        case 47: p = "Ag"; break;
+        default: p = "";   break;
+    }
+    return p;
+}
+
+
 static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
                     int file_version, gmx_groups_t *groups, int atnr)
 {
-    int i, myngrp;
+    int    i, myngrp;
+    char * p_elem;
 
     gmx_fio_do_real(fio, atom->m);
     gmx_fio_do_real(fio, atom->q);
@@ -2397,6 +2439,15 @@ static void do_atom(t_fileio *fio, t_atom *atom, int ngrp, gmx_bool bRead,
     if (file_version >= 52)
     {
         gmx_fio_do_int(fio, atom->atomnumber);
+        if (bRead)
+        {
+            /* Set element string from atomic number if present.
+             * This routine returns an empty string if the name is not found.
+             */
+            strncpy(atom->elem, atomicnumber_to_element(atom->atomnumber), 4);
+            /* avoid warnings about potentially unterminated string */
+            atom->elem[3] = '\0';
+        }
     }
     else if (bRead)
     {
@@ -2982,8 +3033,8 @@ static void set_disres_npair(gmx_mtop_t *mtop)
 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
                     int file_version)
 {
-    int      mt, mb, i;
-    t_blocka dumb;
+    int            mt, mb, i;
+    t_blocka       dumb;
 
     if (bRead)
     {
index 1ea7e0b7fe0345ff8ddcc2a2aecee3e63d4cec6b..516bd4050f766b8ab4c9baf2985f6bcefe06b2c4 100644 (file)
@@ -1249,14 +1249,20 @@ int gmx_trjconv(int argc, char *argv[])
         /* Make atoms struct for output in GRO or PDB files */
         if ((ftp == efGRO) || ((ftp == efG96) && bTPS) || (ftp == efPDB))
         {
-            /* get memory for stuff to go in .pdb file */
-            init_t_atoms(&useatoms, atoms->nr, FALSE);
+            /* get memory for stuff to go in .pdb file, and initialize
+             * the pdbinfo structure part if the input has it.
+             */
+            init_t_atoms(&useatoms, atoms->nr, (atoms->pdbinfo != NULL));
             sfree(useatoms.resinfo);
             useatoms.resinfo = atoms->resinfo;
             for (i = 0; (i < nout); i++)
             {
                 useatoms.atomname[i] = atoms->atomname[index[i]];
                 useatoms.atom[i]     = atoms->atom[index[i]];
+                if (atoms->pdbinfo != NULL)
+                {
+                    useatoms.pdbinfo[i]  = atoms->pdbinfo[index[i]];
+                }
                 useatoms.nres        = max(useatoms.nres, useatoms.atom[i].resind+1);
             }
             useatoms.nr = nout;