Enable 4-letter resname in PDB output, keeps more pdbinfo.
authorErik Lindahl <erik@kth.se>
Sat, 21 Jun 2014 21:30:50 +0000 (23:30 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 23 Jun 2014 15:49:12 +0000 (17:49 +0200)
This still fully adheres to the PDB standard since column 21
is not used by the standard. All common programs (PyMol, VMD, etc)
understand the 4-letter format, and programs that only read three
letters will still read the same filename as they used to. In
particular, this conserves most residue names during pdb<->gro
format conversions. We have also killed the non-standard
wide pdb format to avoid writing broken PDB files.

Fixes #725. Refs #917.

Change-Id: I9b6b8f2e191acdfb65ca2b5d96f39249cd71ea98

src/gromacs/fileio/pdbio.c
src/gromacs/fileio/pdbio.h
src/gromacs/gmxana/gmx_anaeig.c
src/gromacs/gmxana/gmx_chi.c
src/gromacs/gmxana/gmx_editconf.c
src/gromacs/mdlib/constr.c
src/gromacs/mdlib/domdec.c
src/gromacs/mdlib/pme.c
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesMaxPDBOutput.xml
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesPDBOutputWithPDBInput.xml
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesSelectedPDBOutput.xml

index 16fa4a23dba0b04afd903b9d02e6520116a5ea18..439b80a59e5bad4b3557ce447e0e0a3ff78f5bc4 100644 (file)
@@ -75,16 +75,8 @@ static const char *pdbtp[epdbNR] = {
 };
 
 
-/* this is not very good,
-   but these are only used in gmx_trjconv and gmx_editconv */
-static gmx_bool bWideFormat = FALSE;
 #define REMARK_SIM_BOX "REMARK    THIS IS A SIMULATION BOX"
 
-void set_pdb_wide_format(gmx_bool bSet)
-{
-    bWideFormat = bSet;
-}
-
 static void xlate_atomname_pdb2gmx(char *name)
 {
     int  i, length;
@@ -271,10 +263,12 @@ void write_pdbfile_indexed(FILE *out, const char *title,
                            gmx_conect conect, gmx_bool bTerSepChains)
 {
     gmx_conect_t     *gc = (gmx_conect_t *)conect;
-    char              resnm[6], nm[6], pdbform[128], pukestring[100];
+    char              resnm[6], nm[6], pukestring[100];
     atom_id           i, ii;
-    int               resind, resnr, type;
+    int               resind, resnr;
+    enum PDB_record   type;
     unsigned char     resic, ch;
+    char              altloc;
     real              occup, bfac;
     gmx_bool          bOccup;
     int               nlongname = 0;
@@ -288,11 +282,6 @@ void write_pdbfile_indexed(FILE *out, const char *title,
 
     bromacs(pukestring, 99);
     fprintf(out, "TITLE     %s\n", (title && title[0]) ? title : pukestring);
-    if (bWideFormat)
-    {
-        fprintf(out, "REMARK    This file does not adhere to the PDB standard\n");
-        fprintf(out, "REMARK    As a result of, some programs may not like it\n");
-    }
     if (box && ( norm2(box[XX]) || norm2(box[YY]) || norm2(box[ZZ]) ) )
     {
         gmx_write_pdb_box(out, ePBC, box);
@@ -343,7 +332,10 @@ void write_pdbfile_indexed(FILE *out, const char *title,
         }
 
         strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
+        resnm[sizeof(resnm)-1] = 0;
         strncpy(nm, *atoms->atomname[i], sizeof(nm)-1);
+        nm[sizeof(nm)-1] = 0;
+
         /* rename HG12 to 2HG1, etc. */
         xlate_atomname_gmx2pdb(nm);
         resnr = atoms->resinfo[resind].nr;
@@ -367,53 +359,40 @@ void write_pdbfile_indexed(FILE *out, const char *title,
         }
         if (atoms->pdbinfo)
         {
-            type  = atoms->pdbinfo[i].type;
+            type   = (enum PDB_record)(atoms->pdbinfo[i].type);
+            altloc = atoms->pdbinfo[i].altloc;
+            if (!isalnum(altloc))
+            {
+                altloc = ' ';
+            }
             occup = bOccup ? 1.0 : atoms->pdbinfo[i].occup;
             bfac  = atoms->pdbinfo[i].bfac;
         }
         else
         {
-            type  = 0;
-            occup = 1.0;
-            bfac  = 0.0;
-        }
-        if (bWideFormat)
-        {
-            strcpy(pdbform,
-                   "%-6s%5u %-4.4s %3.3s %c%4d%c   %10.5f%10.5f%10.5f%8.4f%8.4f    %2s\n");
-        }
-        else
-        {
-            /* Check whether atomname is an element name */
-            if ((strlen(nm) < 4) && (gmx_strcasecmp(nm, atoms->atom[i].elem) != 0))
-            {
-                strcpy(pdbform, get_pdbformat());
-            }
-            else
-            {
-                strcpy(pdbform, get_pdbformat4());
-                if (strlen(nm) > 4)
-                {
-                    int maxwln = 20;
-                    if (nlongname < maxwln)
-                    {
-                        fprintf(stderr, "WARNING: Writing out atom name (%s) longer than 4 characters to .pdb file\n", nm);
-                    }
-                    else if (nlongname == maxwln)
-                    {
-                        fprintf(stderr, "WARNING: More than %d long atom names, will not write more warnings\n", maxwln);
-                    }
-                    nlongname++;
-                }
-            }
-            strcat(pdbform, "%6.2f%6.2f          %2s\n");
+            type   = epdbATOM;
+            occup  = 1.0;
+            bfac   = 0.0;
+            altloc = ' ';
         }
-        fprintf(out, pdbform, pdbtp[type], (i+1)%100000, nm, resnm, ch, resnr,
-                (resic == '\0') ? ' ' : resic,
-                10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], occup, bfac, atoms->atom[i].elem);
+
+        gmx_fprintf_pdb_atomline(out,
+                                 type,
+                                 i+1,
+                                 nm,
+                                 altloc,
+                                 resnm,
+                                 ch,
+                                 resnr,
+                                 resic,
+                                 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
+                                 occup,
+                                 bfac,
+                                 atoms->atom[i].elem);
+
         if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
         {
-            fprintf(out, "ANISOU%5u  %-4.4s%3.3s %c%4d%c %7d%7d%7d%7d%7d%7d\n",
+            fprintf(out, "ANISOU%5u  %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n",
                     (i+1)%100000, nm, resnm, ch, resnr,
                     (resic == '\0') ? ' ' : resic,
                     atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
@@ -596,7 +575,7 @@ static int read_atom(t_symtab *symtab,
     t_atom       *atomn;
     int           j, k;
     char          nc = '\0';
-    char          anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12];
+    char          anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
     char          xc[12], yc[12], zc[12], occup[12], bfac[12];
     unsigned char resic;
     char          chainid;
@@ -679,6 +658,17 @@ static int read_atom(t_symtab *symtab,
     }
     bfac[k] = nc;
 
+    /* 10 blanks */
+    j += 10;
+
+    /* Element name */
+    for (k = 0; (k < 2); k++, j++)
+    {
+        elem[k] = line[j];
+    }
+    elem[k] = nc;
+    trim(elem);
+
     if (atoms->atom)
     {
         atomn = &(atoms->atom[natom]);
@@ -710,7 +700,7 @@ static int read_atom(t_symtab *symtab,
         atomn->m               = 0.0;
         atomn->q               = 0.0;
         atomn->atomnumber      = atomnumber;
-        atomn->elem[0]         = '\0';
+        strncpy(atomn->elem, elem, 4);
     }
     x[natom][XX] = strtod(xc, NULL)*0.1;
     x[natom][YY] = strtod(yc, NULL)*0.1;
@@ -1060,14 +1050,83 @@ gmx_conect gmx_conect_generate(t_topology *top)
     return gc;
 }
 
-const char* get_pdbformat()
+int
+gmx_fprintf_pdb_atomline(FILE *            fp,
+                         enum PDB_record   record,
+                         int               atom_seq_number,
+                         const char *      atom_name,
+                         char              alternate_location,
+                         const char *      res_name,
+                         char              chain_id,
+                         int               res_seq_number,
+                         char              res_insertion_code,
+                         real              x,
+                         real              y,
+                         real              z,
+                         real              occupancy,
+                         real              b_factor,
+                         const char *      element)
 {
-    static const char *pdbformat = "%-6s%5u  %-4.4s%3.3s %c%4d%c   %8.3f%8.3f%8.3f";
-    return pdbformat;
-}
+    char     tmp_atomname[6], tmp_resname[6];
+    gmx_bool start_name_in_col13;
+    int      n;
 
-const char* get_pdbformat4()
-{
-    static const char *pdbformat4 = "%-6s%5u %-4.4s %3.3s %c%4d%c   %8.3f%8.3f%8.3f";
-    return pdbformat4;
+    if (record != epdbATOM && record != epdbHETATM)
+    {
+        gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
+    }
+
+    /* Format atom name */
+    if (atom_name != NULL)
+    {
+        /* If the atom name is an element name with two chars, it should start already in column 13.
+         * Otherwise it should start in column 14, unless the name length is 4 chars.
+         */
+        if ( (element != NULL) && (strlen(element) >= 2) && (gmx_strncasecmp(atom_name, element, 2) == 0) )
+        {
+            start_name_in_col13 = TRUE;
+        }
+        else
+        {
+            start_name_in_col13 = (strlen(atom_name) >= 4);
+        }
+        sprintf(tmp_atomname, start_name_in_col13 ? "" : " ");
+        strncat(tmp_atomname, atom_name, 4);
+        tmp_atomname[5] = '\0';
+    }
+    else
+    {
+        tmp_atomname[0] = '\0';
+    }
+
+    /* Format residue name */
+    strncpy(tmp_resname, (res_name != NULL) ? res_name : "", 4);
+    /* Make sure the string is terminated if strlen was > 4 */
+    tmp_resname[4] = '\0';
+    /* String is properly terminated, so now we can use strcat. By adding a
+     * space we can write it right-justified, and if the original name was
+     * three characters or less there will be a space added on the right side.
+     */
+    strcat(tmp_resname, " ");
+
+    /* Truncate integers so they fit */
+    atom_seq_number = atom_seq_number % 100000;
+    res_seq_number  = res_seq_number % 10000;
+
+    n = fprintf(fp,
+                "%-6s%5d %-4.4s%c%4.4s%c%4d%c   %8.3f%8.3f%8.3f%6.2f%6.2f          %2s\n",
+                pdbtp[record],
+                atom_seq_number,
+                tmp_atomname,
+                alternate_location,
+                tmp_resname,
+                chain_id,
+                res_seq_number,
+                res_insertion_code,
+                x, y, z,
+                occupancy,
+                b_factor,
+                (element != NULL) ? element : "");
+
+    return n;
 }
index 18745d2c87b8971da088822ea34f2bb2ae94dbfa..ab53be32925cc2863acbc36776a5e02db3e4213f 100644 (file)
@@ -2,8 +2,8 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
- * Copyright (c) 2001-2004, The GROMACS development team,
- * Copyright (c) 2013, by the GROMACS development team, led by
+ * Copyright (c) 2001-2004, The GROMACS development team.
+ * Copyright (c) 2013,2014, 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.
@@ -49,27 +49,41 @@ extern "C" {
 
 typedef struct gmx_conect_t *gmx_conect;
 
-/* THE pdb format (for ATOM/HETATOM lines) */
-const char* get_pdbformat(void);
-const char* get_pdbformat4(void);
-
 /* Enumerated type for pdb records. The other entries are ignored
  * when reading a pdb file
  */
-enum {
+enum PDB_record {
     epdbATOM,   epdbHETATM, epdbANISOU, epdbCRYST1, epdbCOMPND,
     epdbMODEL,  epdbENDMDL, epdbTER,    epdbHEADER, epdbTITLE, epdbREMARK,
     epdbCONECT, epdbNR
 };
 
+/* Write a PDB line with an ATOM or HETATM record directly to a file pointer.
+ *
+ * Returns the number of characters printed.
+ */
+int
+gmx_fprintf_pdb_atomline(FILE *            fp,
+                         enum PDB_record   record,
+                         int               atom_seq_number,
+                         const char *      atom_name,
+                         char              alternate_location,
+                         const char *      res_name,
+                         char              chain_id,
+                         int               res_seq_number,
+                         char              res_insertion_code,
+                         real              x,
+                         real              y,
+                         real              z,
+                         real              occupancy,
+                         real              b_factor,
+                         const char *      element);
+
 /* Enumerated value for indexing an uij entry (anisotropic temperature factors) */
 enum {
     U11, U22, U33, U12, U13, U23
 };
 
-void set_pdb_wide_format(gmx_bool bSet);
-/* If bSet, use wider format for occupancy and bfactor */
-
 void pdb_use_ter(gmx_bool bSet);
 /* set read_pdbatoms to read upto 'TER' or 'ENDMDL' (default, bSet=FALSE).
    This function is fundamentally broken as far as thread-safety is concerned.*/
index 053e9c585ffde54f6ed2589769737327db30d6b2..21d4f836ceb3fa2930b5113986c01454450ffcd7 100644 (file)
@@ -685,7 +685,7 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
         rvec       *x;
         real       *b = NULL;
         matrix      box;
-        char       *resnm, *atnm, pdbform[STRLEN];
+        char       *resnm, *atnm;
         gmx_bool    bPDB, b4D;
         FILE       *out;
 
@@ -746,9 +746,6 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
         }
         if ( ( b4D || bSplit ) && bPDB)
         {
-            strcpy(pdbform, get_pdbformat());
-            strcat(pdbform, "%8.4f%8.4f\n");
-
             out = gmx_ffopen(threedplotfile, "w");
             fprintf(out, "HEADER    %s\n", str);
             if (b4D)
@@ -763,8 +760,8 @@ static void project(const char *trajfile, t_topology *top, int ePBC, matrix topb
                     fprintf(out, "TER\n");
                     j = 0;
                 }
-                fprintf(out, pdbform, "ATOM", i+1, "C", "PRJ", ' ', j+1,
-                        10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 10*b[i]);
+                gmx_fprintf_pdb_atomline(out, epdbATOM, i+1, "C", ' ', "PRJ", ' ', j+1, ' ',
+                                         10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 10*b[i], "");
                 if (j > 0)
                 {
                     fprintf(out, "CONECT%5d%5d\n", i, i+1);
index be8f57e5d75bfeaee8a27283c6d7adbe94304f1a..90362c04d07e4e940727ffd0d90bfa4dff2fe589 100644 (file)
@@ -1057,7 +1057,6 @@ static void order_params(FILE *log,
 {
     FILE *fp;
     int   nh[edMax];
-    char  buf[STRLEN];
     int   i, Dih, Xi;
     real  S2Max, S2Min;
 
@@ -1170,11 +1169,10 @@ static void order_params(FILE *log,
         x0 *= 10.0; /* nm -> angstrom */
         y0 *= 10.0; /* nm -> angstrom */
         z0 *= 10.0; /* nm -> angstrom */
-        sprintf(buf, "%s%%6.f%%6.2f\n", get_pdbformat());
         for (i = 0; (i < 10); i++)
         {
-            fprintf(fp, buf, "ATOM  ", atoms->nr+1+i, "CA", "LEG", ' ',
-                    atoms->nres+1, ' ', x0, y0, z0+(1.2*i), 0.0, -0.1*i);
+            gmx_fprintf_pdb_atomline(fp, epdbATOM, atoms->nr+1+i, "CA", ' ', "LEG", ' ', atoms->nres+1, ' ',
+                                     x0, y0, z0+(1.2*i), 0.0, -0.1*i, "");
         }
         gmx_ffclose(fp);
     }
index 28adc9e5a2f4a6f1c706b6b74f8255f0087f7ef9..cd55910022a5f869a8c1c40bd0f6dbb390f3ef78 100644 (file)
@@ -428,10 +428,8 @@ void visualize_box(FILE *out, int a0, int r0, matrix box, rvec gridsize)
 
         for (i = 0; i < nat; i++)
         {
-            fprintf(out, get_pdbformat(), "ATOM", a0 + i, "C", "BOX", 'K' + i
-                    / NCUCVERT, r0 + i, 10 * vert[i][XX], 10 * vert[i][YY], 10
-                    * vert[i][ZZ]);
-            fprintf(out, "\n");
+            gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i / NCUCVERT, r0 + i, ' ',
+                                     10*vert[i][XX], 10*vert[i][YY], 10*vert[i][ZZ], 1.0, 0.0, "");
         }
 
         edge = compact_unitcell_edges();
@@ -455,10 +453,8 @@ void visualize_box(FILE *out, int a0, int r0, matrix box, rvec gridsize)
             {
                 for (x = 0; x <= 1; x++)
                 {
-                    fprintf(out, get_pdbformat(), "ATOM", a0 + i, "C", "BOX", 'K' + i
-                            / 8, r0 + i, x * 10 * box[XX][XX],
-                            y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ]);
-                    fprintf(out, "\n");
+                    gmx_fprintf_pdb_atomline(out, epdbATOM, a0 + i, "C", ' ', "BOX", 'K' + i/8, r0+i, ' ',
+                                             x * 10 * box[XX][XX], y * 10 * box[YY][YY], z * 10 * box[ZZ][ZZ], 1.0, 0.0, "");
                     i++;
                 }
             }
@@ -1295,7 +1291,6 @@ int gmx_editconf(int argc, char *argv[])
             out = gmx_ffopen(outfile, "w");
             if (bMead)
             {
-                set_pdb_wide_format(TRUE);
                 fprintf(out, "REMARK    "
                         "The B-factors in this file hold atomic radii\n");
                 fprintf(out, "REMARK    "
index 1f5addc99b990da065e764d3c2095f16f9f982eb..c12d013040d651233ca9cc61427ef0752d6684f5 100644 (file)
@@ -218,7 +218,7 @@ static void write_constr_pdb(const char *fn, const char *title,
                              int start, int homenr, t_commrec *cr,
                              rvec x[], matrix box)
 {
-    char          fname[STRLEN], format[STRLEN];
+    char          fname[STRLEN];
     FILE         *out;
     int           dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
     gmx_domdec_t *dd;
@@ -241,7 +241,6 @@ static void write_constr_pdb(const char *fn, const char *title,
     {
         sprintf(fname, "%s.pdb", fn);
     }
-    sprintf(format, "%s\n", get_pdbformat());
 
     out = gmx_fio_fopen(fname, "w");
 
@@ -262,9 +261,8 @@ static void write_constr_pdb(const char *fn, const char *title,
             ii = i;
         }
         gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
-        fprintf(out, format, "ATOM", (ii+1)%100000,
-                anm, resnm, ' ', resnr%10000, ' ',
-                10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
+        gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
+                                 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
     }
     fprintf(out, "TER\n");
 
index 102847d111d01f35d24d7467c965deef99ca4f8a..ea2e7a6d42ea38854f796c72d98d8d594cbd5c3b 100644 (file)
@@ -1925,7 +1925,7 @@ static void write_dd_grid_pdb(const char *fn, gmx_int64_t step,
                               gmx_domdec_t *dd, matrix box, gmx_ddbox_t *ddbox)
 {
     rvec   grid_s[2], *grid_r = NULL, cx, r;
-    char   fname[STRLEN], format[STRLEN], buf[22];
+    char   fname[STRLEN], buf[22];
     FILE  *out;
     int    a, i, d, z, y, x;
     matrix tric;
@@ -1965,7 +1965,6 @@ static void write_dd_grid_pdb(const char *fn, gmx_int64_t step,
             }
         }
         sprintf(fname, "%s_%s.pdb", fn, gmx_step_str(step, buf));
-        sprintf(format, "%s%s\n", get_pdbformat(), "%6.2f%6.2f");
         out = gmx_fio_fopen(fname, "w");
         gmx_write_pdb_box(out, dd->bScrewPBC ? epbcSCREW : epbcXYZ, box);
         a = 1;
@@ -1986,8 +1985,8 @@ static void write_dd_grid_pdb(const char *fn, gmx_int64_t step,
                         cx[YY] = grid_r[i*2+y][YY];
                         cx[ZZ] = grid_r[i*2+z][ZZ];
                         mvmul(tric, cx, r);
-                        fprintf(out, format, "ATOM", a++, "CA", "GLY", ' ', 1+i,
-                                ' ', 10*r[XX], 10*r[YY], 10*r[ZZ], 1.0, vol);
+                        gmx_fprintf_pdb_atomline(out, epdbATOM, a++, "CA", ' ', "GLY", ' ', i+1, ' ',
+                                                 10*r[XX], 10*r[YY], 10*r[ZZ], 1.0, vol, "");
                     }
                 }
             }
@@ -2014,7 +2013,7 @@ void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
                   gmx_mtop_t *mtop, t_commrec *cr,
                   int natoms, rvec x[], matrix box)
 {
-    char          fname[STRLEN], format[STRLEN], format4[STRLEN], buf[22];
+    char          fname[STRLEN], buf[22];
     FILE         *out;
     int           i, ii, resnr, c;
     char         *atomname, *resname;
@@ -2029,9 +2028,6 @@ void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
 
     sprintf(fname, "%s_%s_n%d.pdb", fn, gmx_step_str(step, buf), cr->sim_nodeid);
 
-    sprintf(format, "%s%s\n", get_pdbformat(), "%6.2f%6.2f");
-    sprintf(format4, "%s%s\n", get_pdbformat4(), "%6.2f%6.2f");
-
     out = gmx_fio_fopen(fname, "w");
 
     fprintf(out, "TITLE     %s\n", title);
@@ -2057,10 +2053,8 @@ void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
         {
             b = dd->comm->zones.n + 1;
         }
-        fprintf(out, strlen(atomname) < 4 ? format : format4,
-                "ATOM", (ii+1)%100000,
-                atomname, resname, ' ', resnr%10000, ' ',
-                10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, b);
+        gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, atomname, ' ', resname, ' ', resnr, ' ',
+                                 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, b, "");
     }
     fprintf(out, "TER\n");
 
index af339e7ec8b03e009e648c96fd11d9cdd8b55e42..90a43f97ef1232991218c5c3ede2ca9800abab33 100644 (file)
@@ -1097,13 +1097,12 @@ static int copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid,
     {
 #ifdef DEBUG_PME
         FILE *fp, *fp2;
-        char  fn[STRLEN], format[STRLEN];
+        char  fn[STRLEN];
         real  val;
         sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
         fp = gmx_ffopen(fn, "w");
         sprintf(fn, "pmegrid%d.txt", pme->nodeid);
         fp2 = gmx_ffopen(fn, "w");
-        sprintf(format, "%s%s\n", pdbformat, "%6.2f%6.2f");
 #endif
 
         for (ix = 0; ix < local_fft_ndata[XX]; ix++)
@@ -1119,8 +1118,8 @@ static int copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid,
                     val = 100*pmegrid[pmeidx];
                     if (pmegrid[pmeidx] != 0)
                     {
-                        fprintf(fp, format, "ATOM", pmeidx, "CA", "GLY", ' ', pmeidx, ' ',
-                                5.0*ix, 5.0*iy, 5.0*iz, 1.0, val);
+                        gmx_fprintf_pdb_atomline(fp, epdbATOM, pmeidx, "CA", ' ', "GLY", ' ', pmeidx, ' ',
+                                                 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val, "");
                     }
                     if (pmegrid[pmeidx] != 0)
                     {
index c6d5223262c9688104f352126b9af71dc6e10e7b..0e930e2f22e4f902d4e14a92a199567d007dd54d 100644 (file)
@@ -147,9 +147,9 @@ ATOM      6  S2   RB A   2      20.000  20.000   0.000  1.00  0.50
 ATOM      7  CB   RA B   1      20.000  30.000   0.000  1.50  0.00            
 ATOM      8  S1   RA B   1      20.000  40.000   0.000  1.50  0.30            
 ATOM      9  S2   RA B   1      30.000  10.000   0.000  2.00  0.30            
-ATOM     13  CB   RD B   2      40.000  10.000   0.000  1.00  0.00            
-ATOM     14  S1   RD B   2      40.000  20.000   0.000  0.50  0.00            
-ATOM     15  S2   RD B   2      40.000  30.000   0.000  0.00  0.00            
+ATOM     13  CB A RD B   2      40.000  10.000   0.000  1.00  0.00            
+ATOM     14  S1 A RD B   2      40.000  20.000   0.000  0.50  0.00            
+ATOM     15  S2 A RD B   2      40.000  30.000   0.000  0.00  0.00            
 TER
 ENDMDL
 ]]></String>
index b859a5346f6b351cbfc360ea71c31ce74569384e..90921b3a97f8a817492e5d4cd2a9c1ac8337f768 100644 (file)
@@ -114,9 +114,9 @@ ATOM      9  S2   RA B   1      30.000  10.000   0.000  1.00  0.30
 ATOM     10  CB   RC B   1A     30.000  20.000   0.000  0.00  0.00            
 ATOM     11  S1   RC B   1A     30.000  30.000   0.000  0.00  0.00            
 ATOM     12  S2   RC B   1A     30.000  40.000   0.000  0.00  0.30            
-ATOM     13  CB   RD B   2      40.000  10.000   0.000  1.00  0.00            
-ATOM     14  S1   RD B   2      40.000  20.000   0.000  0.50  0.00            
-ATOM     15  S2   RD B   2      40.000  30.000   0.000  0.00  0.00            
+ATOM     13  CB A RD B   2      40.000  10.000   0.000  1.00  0.00            
+ATOM     14  S1 A RD B   2      40.000  20.000   0.000  0.50  0.00            
+ATOM     15  S2 A RD B   2      40.000  30.000   0.000  0.00  0.00            
 TER
 ENDMDL
 ]]></String>
index b92ae8e201cc65ddb8b9f18de467bd90058913e3..103ca163b31ce2db04aca283a31ebed49e3fcbf6 100644 (file)
@@ -147,8 +147,8 @@ ATOM      6  S2   RB A   2      20.000  20.000   0.000  1.00  0.50
 ATOM      7  CB   RA B   1      20.000  30.000   0.000  1.50  0.00            
 ATOM      8  S1   RA B   1      20.000  40.000   0.000  1.50  0.30            
 ATOM      9  S2   RA B   1      30.000  10.000   0.000  2.00  0.30            
-ATOM     13  CB   RD B   2      40.000  10.000   0.000  1.00  0.00            
-ATOM     14  S1   RD B   2      40.000  20.000   0.000  0.50  0.00            
+ATOM     13  CB A RD B   2      40.000  10.000   0.000  1.00  0.00            
+ATOM     14  S1 A RD B   2      40.000  20.000   0.000  0.50  0.00            
 TER
 ENDMDL
 ]]></String>