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, 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.
62 typedef struct gmx_conect_t {
65 gmx_conection_t *conect;
68 static const char *pdbtp[epdbNR] = {
69 "ATOM ", "HETATM", "ANISOU", "CRYST1",
70 "COMPND", "MODEL", "ENDMDL", "TER", "HEADER", "TITLE", "REMARK",
75 /* this is not very good,
76 but these are only used in gmx_trjconv and gmx_editconv */
77 static gmx_bool bWideFormat = FALSE;
78 #define REMARK_SIM_BOX "REMARK THIS IS A SIMULATION BOX"
80 void set_pdb_wide_format(gmx_bool bSet)
85 static void xlate_atomname_pdb2gmx(char *name)
90 length = strlen(name);
91 if (length > 3 && isdigit(name[0]))
94 for (i = 1; i < length; i++)
98 name[length-1] = temp;
102 static void xlate_atomname_gmx2pdb(char *name)
107 length = strlen(name);
108 if (length > 3 && isdigit(name[length-1]))
110 temp = name[length-1];
111 for (i = length-1; i > 0; --i)
120 void gmx_write_pdb_box(FILE *out, int ePBC, matrix box)
122 real alpha, beta, gamma;
126 ePBC = guess_ePBC(box);
129 if (ePBC == epbcNONE)
134 if (norm2(box[YY])*norm2(box[ZZ]) != 0)
136 alpha = RAD2DEG*acos(cos_angle_no_table(box[YY], box[ZZ]));
142 if (norm2(box[XX])*norm2(box[ZZ]) != 0)
144 beta = RAD2DEG*acos(cos_angle_no_table(box[XX], box[ZZ]));
150 if (norm2(box[XX])*norm2(box[YY]) != 0)
152 gamma = RAD2DEG*acos(cos_angle_no_table(box[XX], box[YY]));
158 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
159 if (ePBC != epbcSCREW)
161 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
162 10*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
163 alpha, beta, gamma, "P 1", 1);
167 /* Double the a-vector length and write the correct space group */
168 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
169 20*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
170 alpha, beta, gamma, "P 21 1 1", 1);
175 static void read_cryst1(char *line, int *ePBC, matrix box)
178 char sa[12], sb[12], sc[12], sg[SG_SIZE+1], ident;
179 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
180 int syma, symb, symc;
183 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
186 if (strlen(line) >= 55)
188 strncpy(sg, line+55, SG_SIZE);
194 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
195 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
197 fc = strtod(sc, NULL)*0.1;
198 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
200 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
202 ePBC_file = epbcSCREW;
212 fa = strtod(sa, NULL)*0.1;
213 fb = strtod(sb, NULL)*0.1;
214 fc = strtod(sc, NULL)*0.1;
215 if (ePBC_file == epbcSCREW)
221 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
225 cosa = cos(alpha*DEG2RAD);
233 cosb = cos(beta*DEG2RAD);
241 cosg = cos(gamma*DEG2RAD);
242 sing = sin(gamma*DEG2RAD);
249 box[YY][XX] = fb*cosg;
250 box[YY][YY] = fb*sing;
251 box[ZZ][XX] = fc*cosb;
252 box[ZZ][YY] = fc*(cosa - cosb*cosg)/sing;
253 box[ZZ][ZZ] = sqrt(fc*fc
254 - box[ZZ][XX]*box[ZZ][XX] - box[ZZ][YY]*box[ZZ][YY]);
264 void write_pdbfile_indexed(FILE *out, const char *title,
265 t_atoms *atoms, rvec x[],
266 int ePBC, matrix box, char chainid,
267 int model_nr, atom_id nindex, const atom_id index[],
268 gmx_conect conect, gmx_bool bTerSepChains)
270 gmx_conect_t *gc = (gmx_conect_t *)conect;
271 char resnm[6], nm[6], pdbform[128], pukestring[100];
273 int resind, resnr, type;
274 unsigned char resic, ch;
278 int chainnum, lastchainnum;
279 int lastresind, lastchainresind;
280 gmx_residuetype_t rt;
281 const char *p_restype;
282 const char *p_lastrestype;
284 gmx_residuetype_init(&rt);
286 bromacs(pukestring, 99);
287 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : pukestring);
290 fprintf(out, "REMARK This file does not adhere to the PDB standard\n");
291 fprintf(out, "REMARK As a result of, some programs may not like it\n");
293 if (box && ( norm2(box[XX]) || norm2(box[YY]) || norm2(box[ZZ]) ) )
295 gmx_write_pdb_box(out, ePBC, box);
299 /* Check whether any occupancies are set, in that case leave it as is,
300 * otherwise set them all to one
303 for (ii = 0; (ii < nindex) && bOccup; ii++)
306 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
314 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
316 lastchainresind = -1;
321 for (ii = 0; ii < nindex; ii++)
325 resind = atoms->atom[i].resind;
326 chainnum = atoms->resinfo[resind].chainnum;
327 p_lastrestype = p_restype;
328 gmx_residuetype_get_type(rt, *atoms->resinfo[resind].name, &p_restype);
330 /* Add a TER record if we changed chain, and if either the previous or this chain is protein/DNA/RNA. */
331 if (bTerSepChains && ii > 0 && chainnum != lastchainnum)
333 /* Only add TER if the previous chain contained protein/DNA/RNA. */
334 if (gmx_residuetype_is_protein(rt, p_lastrestype) || gmx_residuetype_is_dna(rt, p_lastrestype) || gmx_residuetype_is_rna(rt, p_lastrestype))
336 fprintf(out, "TER\n");
338 lastchainnum = chainnum;
339 lastchainresind = lastresind;
342 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
343 strncpy(nm, *atoms->atomname[i], sizeof(nm)-1);
344 /* rename HG12 to 2HG1, etc. */
345 xlate_atomname_gmx2pdb(nm);
346 resnr = atoms->resinfo[resind].nr;
347 resic = atoms->resinfo[resind].ic;
354 ch = atoms->resinfo[resind].chainid;
363 resnr = resnr % 10000;
367 type = atoms->pdbinfo[i].type;
368 occup = bOccup ? 1.0 : atoms->pdbinfo[i].occup;
369 bfac = atoms->pdbinfo[i].bfac;
380 "%-6s%5u %-4.4s %3.3s %c%4d%c %10.5f%10.5f%10.5f%8.4f%8.4f %2s\n");
384 /* Check whether atomname is an element name */
385 if ((strlen(nm) < 4) && (gmx_strcasecmp(nm, atoms->atom[i].elem) != 0))
387 strcpy(pdbform, get_pdbformat());
391 strcpy(pdbform, get_pdbformat4());
395 if (nlongname < maxwln)
397 fprintf(stderr, "WARNING: Writing out atom name (%s) longer than 4 characters to .pdb file\n", nm);
399 else if (nlongname == maxwln)
401 fprintf(stderr, "WARNING: More than %d long atom names, will not write more warnings\n", maxwln);
406 strcat(pdbform, "%6.2f%6.2f %2s\n");
408 fprintf(out, pdbform, pdbtp[type], (i+1)%100000, nm, resnm, ch, resnr,
409 (resic == '\0') ? ' ' : resic,
410 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], occup, bfac, atoms->atom[i].elem);
411 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
413 fprintf(out, "ANISOU%5u %-4.4s%3.3s %c%4d%c %7d%7d%7d%7d%7d%7d\n",
414 (i+1)%100000, nm, resnm, ch, resnr,
415 (resic == '\0') ? ' ' : resic,
416 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
417 atoms->pdbinfo[i].uij[2], atoms->pdbinfo[i].uij[3],
418 atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
422 fprintf(out, "TER\n");
423 fprintf(out, "ENDMDL\n");
427 /* Write conect records */
428 for (i = 0; (i < gc->nconect); i++)
430 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai+1, gc->conect[i].aj+1);
434 gmx_residuetype_destroy(rt);
437 void write_pdbfile(FILE *out, const char *title, t_atoms *atoms, rvec x[],
438 int ePBC, matrix box, char chainid, int model_nr, gmx_conect conect, gmx_bool bTerSepChains)
442 snew(index, atoms->nr);
443 for (i = 0; i < atoms->nr; i++)
447 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
448 atoms->nr, index, conect, bTerSepChains);
452 int line2type(char *line)
457 for (k = 0; (k < 6); k++)
463 for (k = 0; (k < epdbNR); k++)
465 if (strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
474 static void read_anisou(char line[], int natom, t_atoms *atoms)
478 char anr[12], anm[12];
482 for (k = 0; (k < 5); k++, j++)
488 for (k = 0; (k < 4); k++, j++)
495 /* Strip off spaces */
498 /* Search backwards for number and name only */
499 atomnr = strtol(anr, NULL, 10);
500 for (i = natom-1; (i >= 0); i--)
502 if ((strcmp(anm, *(atoms->atomname[i])) == 0) &&
503 (atomnr == atoms->pdbinfo[i].atomnr))
510 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n",
515 if (sscanf(line+29, "%d%d%d%d%d%d",
516 &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
517 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
518 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
521 atoms->pdbinfo[i].bAnisotropic = TRUE;
525 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
526 atoms->pdbinfo[i].bAnisotropic = FALSE;
531 void get_pdb_atomnumber(t_atoms *atoms, gmx_atomprop_t aps)
533 int i, atomnumber, len;
535 char anm[6], anm_copy[6], *ptr;
541 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
543 for (i = 0; (i < atoms->nr); i++)
545 strcpy(anm, atoms->pdbinfo[i].atomnm);
546 strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
549 if ((anm[0] != ' ') && ((len <= 2) || ((len > 2) && !isdigit(anm[2]))))
552 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
554 atomnumber = gmx_nint(eval);
559 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
561 atomnumber = gmx_nint(eval);
565 if (atomnumber == NOTSET)
568 while ((k < strlen(anm)) && (isspace(anm[k]) || isdigit(anm[k])))
572 anm_copy[0] = anm[k];
574 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
576 atomnumber = gmx_nint(eval);
579 atoms->atom[i].atomnumber = atomnumber;
580 ptr = gmx_atomprop_element(aps, atomnumber);
581 strncpy(atoms->atom[i].elem, ptr == NULL ? "" : ptr, 4);
584 fprintf(debug, "Atomnumber for atom '%s' is %d\n", anm, atomnumber);
589 static int read_atom(t_symtab *symtab,
590 char line[], int type, int natom,
591 t_atoms *atoms, rvec x[], int chainnum, gmx_bool bChange)
596 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12];
597 char xc[12], yc[12], zc[12], occup[12], bfac[12];
600 int resnr, atomnumber;
602 if (natom >= atoms->nr)
604 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)",
610 for (k = 0; (k < 5); k++, j++)
617 for (k = 0; (k < 4); k++, j++)
622 strcpy(anm_copy, anm);
628 for (k = 0; (k < 4); k++, j++)
638 for (k = 0; (k < 4); k++, j++)
644 resnr = strtol(rnr, NULL, 10);
648 /* X,Y,Z Coordinate */
649 for (k = 0; (k < 8); k++, j++)
654 for (k = 0; (k < 8); k++, j++)
659 for (k = 0; (k < 8); k++, j++)
666 for (k = 0; (k < 6); k++, j++)
673 for (k = 0; (k < 7); k++, j++)
681 atomn = &(atoms->atom[natom]);
683 atoms->resinfo[atoms->atom[natom-1].resind].nr != resnr ||
684 atoms->resinfo[atoms->atom[natom-1].resind].ic != resic ||
685 (strcmp(*atoms->resinfo[atoms->atom[natom-1].resind].name, resnm) != 0))
693 atomn->resind = atoms->atom[natom-1].resind + 1;
695 atoms->nres = atomn->resind + 1;
696 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
700 atomn->resind = atoms->atom[natom-1].resind;
704 xlate_atomname_pdb2gmx(anm);
706 atoms->atomname[natom] = put_symtab(symtab, anm);
709 atomn->atomnumber = atomnumber;
710 atomn->elem[0] = '\0';
712 x[natom][XX] = strtod(xc, NULL)*0.1;
713 x[natom][YY] = strtod(yc, NULL)*0.1;
714 x[natom][ZZ] = strtod(zc, NULL)*0.1;
717 atoms->pdbinfo[natom].type = type;
718 atoms->pdbinfo[natom].atomnr = strtol(anr, NULL, 10);
719 atoms->pdbinfo[natom].altloc = altloc;
720 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
721 atoms->pdbinfo[natom].bfac = strtod(bfac, NULL);
722 atoms->pdbinfo[natom].occup = strtod(occup, NULL);
729 gmx_bool is_hydrogen(const char *nm)
740 else if ((isdigit(buf[0])) && (buf[1] == 'H'))
747 gmx_bool is_dummymass(const char *nm)
754 if ((buf[0] == 'M') && isdigit(buf[strlen(buf)-1]))
762 static void gmx_conect_addline(gmx_conect_t *con, char *line)
765 char format[32], form2[32];
767 sprintf(form2, "%%*s");
768 sprintf(format, "%s%%d", form2);
769 if (sscanf(line, format, &ai) == 1)
773 strcat(form2, "%*s");
774 sprintf(format, "%s%%d", form2);
775 n = sscanf(line, format, &aj);
778 srenew(con->conect, ++con->nconect);
779 con->conect[con->nconect-1].ai = ai-1;
780 con->conect[con->nconect-1].aj = aj-1;
787 void gmx_conect_dump(FILE *fp, gmx_conect conect)
789 gmx_conect_t *gc = (gmx_conect_t *)conect;
792 for (i = 0; (i < gc->nconect); i++)
794 fprintf(fp, "%6s%5d%5d\n", "CONECT",
795 gc->conect[i].ai+1, gc->conect[i].aj+1);
799 gmx_conect gmx_conect_init()
805 return (gmx_conect) gc;
808 void gmx_conect_done(gmx_conect conect)
810 gmx_conect_t *gc = (gmx_conect_t *)conect;
815 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
817 gmx_conect_t *gc = (gmx_conect_t *)conect;
823 for (i = 0; (i < gc->nconect); i++)
825 if (((gc->conect[i].ai == ai) &&
826 (gc->conect[i].aj == aj)) ||
827 ((gc->conect[i].aj == ai) &&
828 (gc->conect[i].ai == aj)))
836 void gmx_conect_add(gmx_conect conect, int ai, int aj)
838 gmx_conect_t *gc = (gmx_conect_t *)conect;
844 if (!gmx_conect_exist(conect, ai, aj))
846 srenew(gc->conect, ++gc->nconect);
847 gc->conect[gc->nconect-1].ai = ai;
848 gc->conect[gc->nconect-1].aj = aj;
852 int read_pdbfile(FILE *in, char *title, int *model_nr,
853 t_atoms *atoms, rvec x[], int *ePBC, matrix box, gmx_bool bChange,
856 gmx_conect_t *gc = (gmx_conect_t *)conect;
859 gmx_bool bConnWarn = FALSE;
863 int natom, chainnum, nres_ter_prev = 0;
865 gmx_bool bStop = FALSE;
869 /* Only assume pbc when there is a CRYST1 entry */
877 open_symtab(&symtab);
883 while (!bStop && (fgets2(line, STRLEN, in) != NULL))
885 line_type = line2type(line);
891 natom = read_atom(&symtab, line, line_type, natom, atoms, x, chainnum, bChange);
897 read_anisou(line, natom, atoms);
902 read_cryst1(line, ePBC, box);
907 if (strlen(line) > 6)
910 /* skip HEADER or TITLE and spaces */
919 /* truncate after title */
933 if ((!strstr(line, ": ")) || (strstr(line+6, "MOLECULE:")))
935 if (!(c = strstr(line+6, "MOLECULE:")) )
939 /* skip 'MOLECULE:' and spaces */
948 /* truncate after title */
952 while ( (d[-1] == ';') && d > c)
981 sscanf(line, "%*s%d", model_nr);
991 gmx_conect_addline(gc, line);
995 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
1005 free_symtab(&symtab);
1009 void get_pdb_coordnum(FILE *in, int *natoms)
1014 while (fgets2(line, STRLEN, in))
1016 if (strncmp(line, "ENDMDL", 6) == 0)
1020 if ((strncmp(line, "ATOM ", 6) == 0) || (strncmp(line, "HETATM", 6) == 0))
1027 void read_pdb_conf(const char *infile, char *title,
1028 t_atoms *atoms, rvec x[], int *ePBC, matrix box, gmx_bool bChange,
1033 in = gmx_fio_fopen(infile, "r");
1034 read_pdbfile(in, title, NULL, atoms, x, ePBC, box, bChange, conect);
1038 gmx_conect gmx_conect_generate(t_topology *top)
1043 /* Fill the conect records */
1044 gc = gmx_conect_init();
1046 for (f = 0; (f < F_NRE); f++)
1050 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms+1)
1052 gmx_conect_add(gc, top->idef.il[f].iatoms[i+1],
1053 top->idef.il[f].iatoms[i+2]);
1060 const char* get_pdbformat()
1062 static const char *pdbformat = "%-6s%5u %-4.4s%3.3s %c%4d%c %8.3f%8.3f%8.3f";
1066 const char* get_pdbformat4()
1068 static const char *pdbformat4 = "%-6s%5u %-4.4s %3.3s %c%4d%c %8.3f%8.3f%8.3f";