2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "gromacs/utility/futil.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
64 typedef struct gmx_conect_t {
67 gmx_conection_t *conect;
70 static const char *pdbtp[epdbNR] = {
71 "ATOM ", "HETATM", "ANISOU", "CRYST1",
72 "COMPND", "MODEL", "ENDMDL", "TER", "HEADER", "TITLE", "REMARK",
77 /* this is not very good,
78 but these are only used in gmx_trjconv and gmx_editconv */
79 static gmx_bool bWideFormat = FALSE;
80 #define REMARK_SIM_BOX "REMARK THIS IS A SIMULATION BOX"
82 void set_pdb_wide_format(gmx_bool bSet)
87 static void xlate_atomname_pdb2gmx(char *name)
92 length = strlen(name);
93 if (length > 3 && isdigit(name[0]))
96 for (i = 1; i < length; i++)
100 name[length-1] = temp;
104 static void xlate_atomname_gmx2pdb(char *name)
109 length = strlen(name);
110 if (length > 3 && isdigit(name[length-1]))
112 temp = name[length-1];
113 for (i = length-1; i > 0; --i)
122 void gmx_write_pdb_box(FILE *out, int ePBC, matrix box)
124 real alpha, beta, gamma;
128 ePBC = guess_ePBC(box);
131 if (ePBC == epbcNONE)
136 if (norm2(box[YY])*norm2(box[ZZ]) != 0)
138 alpha = RAD2DEG*acos(cos_angle_no_table(box[YY], box[ZZ]));
144 if (norm2(box[XX])*norm2(box[ZZ]) != 0)
146 beta = RAD2DEG*acos(cos_angle_no_table(box[XX], box[ZZ]));
152 if (norm2(box[XX])*norm2(box[YY]) != 0)
154 gamma = RAD2DEG*acos(cos_angle_no_table(box[XX], box[YY]));
160 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
161 if (ePBC != epbcSCREW)
163 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
164 10*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
165 alpha, beta, gamma, "P 1", 1);
169 /* Double the a-vector length and write the correct space group */
170 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
171 20*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
172 alpha, beta, gamma, "P 21 1 1", 1);
177 static void read_cryst1(char *line, int *ePBC, matrix box)
180 char sa[12], sb[12], sc[12], sg[SG_SIZE+1], ident;
181 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
182 int syma, symb, symc;
185 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
188 if (strlen(line) >= 55)
190 strncpy(sg, line+55, SG_SIZE);
196 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
197 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
199 fc = strtod(sc, NULL)*0.1;
200 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
202 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
204 ePBC_file = epbcSCREW;
214 fa = strtod(sa, NULL)*0.1;
215 fb = strtod(sb, NULL)*0.1;
216 fc = strtod(sc, NULL)*0.1;
217 if (ePBC_file == epbcSCREW)
223 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
227 cosa = cos(alpha*DEG2RAD);
235 cosb = cos(beta*DEG2RAD);
243 cosg = cos(gamma*DEG2RAD);
244 sing = sin(gamma*DEG2RAD);
251 box[YY][XX] = fb*cosg;
252 box[YY][YY] = fb*sing;
253 box[ZZ][XX] = fc*cosb;
254 box[ZZ][YY] = fc*(cosa - cosb*cosg)/sing;
255 box[ZZ][ZZ] = sqrt(fc*fc
256 - box[ZZ][XX]*box[ZZ][XX] - box[ZZ][YY]*box[ZZ][YY]);
266 void write_pdbfile_indexed(FILE *out, const char *title,
267 t_atoms *atoms, rvec x[],
268 int ePBC, matrix box, char chainid,
269 int model_nr, atom_id nindex, const atom_id index[],
270 gmx_conect conect, gmx_bool bTerSepChains)
272 gmx_conect_t *gc = (gmx_conect_t *)conect;
273 char resnm[6], nm[6], pdbform[128], pukestring[100];
275 int resind, resnr, type;
276 unsigned char resic, ch;
280 int chainnum, lastchainnum;
281 int lastresind, lastchainresind;
282 gmx_residuetype_t rt;
283 const char *p_restype;
284 const char *p_lastrestype;
286 gmx_residuetype_init(&rt);
288 bromacs(pukestring, 99);
289 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : pukestring);
292 fprintf(out, "REMARK This file does not adhere to the PDB standard\n");
293 fprintf(out, "REMARK As a result of, some programs may not like it\n");
295 if (box && ( norm2(box[XX]) || norm2(box[YY]) || norm2(box[ZZ]) ) )
297 gmx_write_pdb_box(out, ePBC, box);
301 /* Check whether any occupancies are set, in that case leave it as is,
302 * otherwise set them all to one
305 for (ii = 0; (ii < nindex) && bOccup; ii++)
308 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
316 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
318 lastchainresind = -1;
323 for (ii = 0; ii < nindex; ii++)
327 resind = atoms->atom[i].resind;
328 chainnum = atoms->resinfo[resind].chainnum;
329 p_lastrestype = p_restype;
330 gmx_residuetype_get_type(rt, *atoms->resinfo[resind].name, &p_restype);
332 /* Add a TER record if we changed chain, and if either the previous or this chain is protein/DNA/RNA. */
333 if (bTerSepChains && ii > 0 && chainnum != lastchainnum)
335 /* Only add TER if the previous chain contained protein/DNA/RNA. */
336 if (gmx_residuetype_is_protein(rt, p_lastrestype) || gmx_residuetype_is_dna(rt, p_lastrestype) || gmx_residuetype_is_rna(rt, p_lastrestype))
338 fprintf(out, "TER\n");
340 lastchainnum = chainnum;
341 lastchainresind = lastresind;
344 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
345 strncpy(nm, *atoms->atomname[i], sizeof(nm)-1);
346 /* rename HG12 to 2HG1, etc. */
347 xlate_atomname_gmx2pdb(nm);
348 resnr = atoms->resinfo[resind].nr;
349 resic = atoms->resinfo[resind].ic;
356 ch = atoms->resinfo[resind].chainid;
365 resnr = resnr % 10000;
369 type = atoms->pdbinfo[i].type;
370 occup = bOccup ? 1.0 : atoms->pdbinfo[i].occup;
371 bfac = atoms->pdbinfo[i].bfac;
382 "%-6s%5u %-4.4s %3.3s %c%4d%c %10.5f%10.5f%10.5f%8.4f%8.4f %2s\n");
386 /* Check whether atomname is an element name */
387 if ((strlen(nm) < 4) && (gmx_strcasecmp(nm, atoms->atom[i].elem) != 0))
389 strcpy(pdbform, get_pdbformat());
393 strcpy(pdbform, get_pdbformat4());
397 if (nlongname < maxwln)
399 fprintf(stderr, "WARNING: Writing out atom name (%s) longer than 4 characters to .pdb file\n", nm);
401 else if (nlongname == maxwln)
403 fprintf(stderr, "WARNING: More than %d long atom names, will not write more warnings\n", maxwln);
408 strcat(pdbform, "%6.2f%6.2f %2s\n");
410 fprintf(out, pdbform, pdbtp[type], (i+1)%100000, nm, resnm, ch, resnr,
411 (resic == '\0') ? ' ' : resic,
412 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], occup, bfac, atoms->atom[i].elem);
413 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
415 fprintf(out, "ANISOU%5u %-4.4s%3.3s %c%4d%c %7d%7d%7d%7d%7d%7d\n",
416 (i+1)%100000, nm, resnm, ch, resnr,
417 (resic == '\0') ? ' ' : resic,
418 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
419 atoms->pdbinfo[i].uij[2], atoms->pdbinfo[i].uij[3],
420 atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
424 fprintf(out, "TER\n");
425 fprintf(out, "ENDMDL\n");
429 /* Write conect records */
430 for (i = 0; (i < gc->nconect); i++)
432 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai+1, gc->conect[i].aj+1);
436 gmx_residuetype_destroy(rt);
439 void write_pdbfile(FILE *out, const char *title, t_atoms *atoms, rvec x[],
440 int ePBC, matrix box, char chainid, int model_nr, gmx_conect conect, gmx_bool bTerSepChains)
444 snew(index, atoms->nr);
445 for (i = 0; i < atoms->nr; i++)
449 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
450 atoms->nr, index, conect, bTerSepChains);
454 int line2type(char *line)
459 for (k = 0; (k < 6); k++)
465 for (k = 0; (k < epdbNR); k++)
467 if (strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
476 static void read_anisou(char line[], int natom, t_atoms *atoms)
480 char anr[12], anm[12];
484 for (k = 0; (k < 5); k++, j++)
490 for (k = 0; (k < 4); k++, j++)
497 /* Strip off spaces */
500 /* Search backwards for number and name only */
501 atomnr = strtol(anr, NULL, 10);
502 for (i = natom-1; (i >= 0); i--)
504 if ((strcmp(anm, *(atoms->atomname[i])) == 0) &&
505 (atomnr == atoms->pdbinfo[i].atomnr))
512 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n",
517 if (sscanf(line+29, "%d%d%d%d%d%d",
518 &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
519 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
520 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
523 atoms->pdbinfo[i].bAnisotropic = TRUE;
527 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
528 atoms->pdbinfo[i].bAnisotropic = FALSE;
533 void get_pdb_atomnumber(t_atoms *atoms, gmx_atomprop_t aps)
535 int i, atomnumber, len;
537 char anm[6], anm_copy[6], *ptr;
543 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
545 for (i = 0; (i < atoms->nr); i++)
547 strcpy(anm, atoms->pdbinfo[i].atomnm);
548 strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
551 if ((anm[0] != ' ') && ((len <= 2) || ((len > 2) && !isdigit(anm[2]))))
554 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
556 atomnumber = gmx_nint(eval);
561 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
563 atomnumber = gmx_nint(eval);
567 if (atomnumber == NOTSET)
570 while ((k < strlen(anm)) && (isspace(anm[k]) || isdigit(anm[k])))
574 anm_copy[0] = anm[k];
576 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
578 atomnumber = gmx_nint(eval);
581 atoms->atom[i].atomnumber = atomnumber;
582 ptr = gmx_atomprop_element(aps, atomnumber);
583 strncpy(atoms->atom[i].elem, ptr == NULL ? "" : ptr, 4);
586 fprintf(debug, "Atomnumber for atom '%s' is %d\n", anm, atomnumber);
591 static int read_atom(t_symtab *symtab,
592 char line[], int type, int natom,
593 t_atoms *atoms, rvec x[], int chainnum, gmx_bool bChange)
598 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12];
599 char xc[12], yc[12], zc[12], occup[12], bfac[12];
602 int resnr, atomnumber;
604 if (natom >= atoms->nr)
606 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)",
612 for (k = 0; (k < 5); k++, j++)
619 for (k = 0; (k < 4); k++, j++)
624 strcpy(anm_copy, anm);
630 for (k = 0; (k < 4); k++, j++)
640 for (k = 0; (k < 4); k++, j++)
646 resnr = strtol(rnr, NULL, 10);
650 /* X,Y,Z Coordinate */
651 for (k = 0; (k < 8); k++, j++)
656 for (k = 0; (k < 8); k++, j++)
661 for (k = 0; (k < 8); k++, j++)
668 for (k = 0; (k < 6); k++, j++)
675 for (k = 0; (k < 7); k++, j++)
683 atomn = &(atoms->atom[natom]);
685 atoms->resinfo[atoms->atom[natom-1].resind].nr != resnr ||
686 atoms->resinfo[atoms->atom[natom-1].resind].ic != resic ||
687 (strcmp(*atoms->resinfo[atoms->atom[natom-1].resind].name, resnm) != 0))
695 atomn->resind = atoms->atom[natom-1].resind + 1;
697 atoms->nres = atomn->resind + 1;
698 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
702 atomn->resind = atoms->atom[natom-1].resind;
706 xlate_atomname_pdb2gmx(anm);
708 atoms->atomname[natom] = put_symtab(symtab, anm);
711 atomn->atomnumber = atomnumber;
712 atomn->elem[0] = '\0';
714 x[natom][XX] = strtod(xc, NULL)*0.1;
715 x[natom][YY] = strtod(yc, NULL)*0.1;
716 x[natom][ZZ] = strtod(zc, NULL)*0.1;
719 atoms->pdbinfo[natom].type = type;
720 atoms->pdbinfo[natom].atomnr = strtol(anr, NULL, 10);
721 atoms->pdbinfo[natom].altloc = altloc;
722 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
723 atoms->pdbinfo[natom].bfac = strtod(bfac, NULL);
724 atoms->pdbinfo[natom].occup = strtod(occup, NULL);
731 gmx_bool is_hydrogen(const char *nm)
742 else if ((isdigit(buf[0])) && (buf[1] == 'H'))
749 gmx_bool is_dummymass(const char *nm)
756 if ((buf[0] == 'M') && isdigit(buf[strlen(buf)-1]))
764 static void gmx_conect_addline(gmx_conect_t *con, char *line)
767 char format[32], form2[32];
769 sprintf(form2, "%%*s");
770 sprintf(format, "%s%%d", form2);
771 if (sscanf(line, format, &ai) == 1)
775 strcat(form2, "%*s");
776 sprintf(format, "%s%%d", form2);
777 n = sscanf(line, format, &aj);
780 srenew(con->conect, ++con->nconect);
781 con->conect[con->nconect-1].ai = ai-1;
782 con->conect[con->nconect-1].aj = aj-1;
789 void gmx_conect_dump(FILE *fp, gmx_conect conect)
791 gmx_conect_t *gc = (gmx_conect_t *)conect;
794 for (i = 0; (i < gc->nconect); i++)
796 fprintf(fp, "%6s%5d%5d\n", "CONECT",
797 gc->conect[i].ai+1, gc->conect[i].aj+1);
801 gmx_conect gmx_conect_init()
807 return (gmx_conect) gc;
810 void gmx_conect_done(gmx_conect conect)
812 gmx_conect_t *gc = (gmx_conect_t *)conect;
817 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
819 gmx_conect_t *gc = (gmx_conect_t *)conect;
825 for (i = 0; (i < gc->nconect); i++)
827 if (((gc->conect[i].ai == ai) &&
828 (gc->conect[i].aj == aj)) ||
829 ((gc->conect[i].aj == ai) &&
830 (gc->conect[i].ai == aj)))
838 void gmx_conect_add(gmx_conect conect, int ai, int aj)
840 gmx_conect_t *gc = (gmx_conect_t *)conect;
846 if (!gmx_conect_exist(conect, ai, aj))
848 srenew(gc->conect, ++gc->nconect);
849 gc->conect[gc->nconect-1].ai = ai;
850 gc->conect[gc->nconect-1].aj = aj;
854 int read_pdbfile(FILE *in, char *title, int *model_nr,
855 t_atoms *atoms, rvec x[], int *ePBC, matrix box, gmx_bool bChange,
858 gmx_conect_t *gc = (gmx_conect_t *)conect;
861 gmx_bool bConnWarn = FALSE;
865 int natom, chainnum, nres_ter_prev = 0;
867 gmx_bool bStop = FALSE;
871 /* Only assume pbc when there is a CRYST1 entry */
879 open_symtab(&symtab);
885 while (!bStop && (fgets2(line, STRLEN, in) != NULL))
887 line_type = line2type(line);
893 natom = read_atom(&symtab, line, line_type, natom, atoms, x, chainnum, bChange);
899 read_anisou(line, natom, atoms);
904 read_cryst1(line, ePBC, box);
909 if (strlen(line) > 6)
912 /* skip HEADER or TITLE and spaces */
921 /* truncate after title */
935 if ((!strstr(line, ": ")) || (strstr(line+6, "MOLECULE:")))
937 if (!(c = strstr(line+6, "MOLECULE:")) )
941 /* skip 'MOLECULE:' and spaces */
950 /* truncate after title */
954 while ( (d[-1] == ';') && d > c)
983 sscanf(line, "%*s%d", model_nr);
993 gmx_conect_addline(gc, line);
997 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
1007 free_symtab(&symtab);
1011 void get_pdb_coordnum(FILE *in, int *natoms)
1016 while (fgets2(line, STRLEN, in))
1018 if (strncmp(line, "ENDMDL", 6) == 0)
1022 if ((strncmp(line, "ATOM ", 6) == 0) || (strncmp(line, "HETATM", 6) == 0))
1029 void read_pdb_conf(const char *infile, char *title,
1030 t_atoms *atoms, rvec x[], int *ePBC, matrix box, gmx_bool bChange,
1035 in = gmx_fio_fopen(infile, "r");
1036 read_pdbfile(in, title, NULL, atoms, x, ePBC, box, bChange, conect);
1040 gmx_conect gmx_conect_generate(t_topology *top)
1045 /* Fill the conect records */
1046 gc = gmx_conect_init();
1048 for (f = 0; (f < F_NRE); f++)
1052 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms+1)
1054 gmx_conect_add(gc, top->idef.il[f].iatoms[i+1],
1055 top->idef.il[f].iatoms[i+2]);
1062 const char* get_pdbformat()
1064 static const char *pdbformat = "%-6s%5u %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f";
1068 const char* get_pdbformat4()
1070 static const char *pdbformat4 = "%-6s%5u %-4.4s %3.3s %c%4d%c %8.3f%8.3f%8.3f";