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,2015,2016,2017,2018,2019, 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.
48 #include "gromacs/fileio/gmxfio.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/topology/atomprop.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/topology/residuetypes.h"
55 #include "gromacs/topology/symtab.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/coolstuff.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/snprintf.h"
68 typedef struct gmx_conect_t {
71 gmx_conection_t *conect;
74 static const char *pdbtp[epdbNR] = {
75 "ATOM ", "HETATM", "ANISOU", "CRYST1",
76 "COMPND", "MODEL", "ENDMDL", "TER", "HEADER", "TITLE", "REMARK",
80 #define REMARK_SIM_BOX "REMARK THIS IS A SIMULATION BOX"
82 static void xlate_atomname_pdb2gmx(char *name)
87 length = std::strlen(name);
88 if (length > 3 && std::isdigit(name[0]))
91 for (i = 1; i < length; i++)
95 name[length-1] = temp;
99 // Deliberately taking a copy of name to return it later
100 static std::string xlate_atomname_gmx2pdb(std::string name)
102 size_t length = name.size();
103 if (length > 3 && std::isdigit(name[length-1]))
105 char temp = name[length-1];
106 for (size_t i = length-1; i > 0; --i)
116 void gmx_write_pdb_box(FILE *out, int ePBC, const matrix box)
118 real alpha, beta, gamma;
122 ePBC = guess_ePBC(box);
125 if (ePBC == epbcNONE)
130 if (norm2(box[YY])*norm2(box[ZZ]) != 0)
132 alpha = RAD2DEG*gmx_angle(box[YY], box[ZZ]);
138 if (norm2(box[XX])*norm2(box[ZZ]) != 0)
140 beta = RAD2DEG*gmx_angle(box[XX], box[ZZ]);
146 if (norm2(box[XX])*norm2(box[YY]) != 0)
148 gamma = RAD2DEG*gmx_angle(box[XX], box[YY]);
154 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
155 if (ePBC != epbcSCREW)
157 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
158 10*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
159 alpha, beta, gamma, "P 1", 1);
163 /* Double the a-vector length and write the correct space group */
164 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
165 20*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
166 alpha, beta, gamma, "P 21 1 1", 1);
171 static void read_cryst1(char *line, int *ePBC, matrix box)
174 char sa[12], sb[12], sc[12], sg[SG_SIZE+1], ident;
175 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
176 int syma, symb, symc;
179 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
182 if (strlen(line) >= 55)
184 strncpy(sg, line+55, SG_SIZE);
190 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
191 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
193 fc = strtod(sc, nullptr)*0.1;
194 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
196 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
198 ePBC_file = epbcSCREW;
208 fa = strtod(sa, nullptr)*0.1;
209 fb = strtod(sb, nullptr)*0.1;
210 fc = strtod(sc, nullptr)*0.1;
211 if (ePBC_file == epbcSCREW)
217 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
221 cosa = std::cos(alpha*DEG2RAD);
229 cosb = std::cos(beta*DEG2RAD);
237 cosg = std::cos(gamma*DEG2RAD);
238 sing = std::sin(gamma*DEG2RAD);
245 box[YY][XX] = fb*cosg;
246 box[YY][YY] = fb*sing;
247 box[ZZ][XX] = fc*cosb;
248 box[ZZ][YY] = fc*(cosa - cosb*cosg)/sing;
249 box[ZZ][ZZ] = std::sqrt(fc*fc
250 - box[ZZ][XX]*box[ZZ][XX] - box[ZZ][YY]*box[ZZ][YY]);
261 gmx_fprintf_pqr_atomline(FILE * fp,
262 enum PDB_record record,
264 const char * atom_name,
265 const char * res_name,
274 GMX_RELEASE_ASSERT(record == epdbATOM || record == epdbHETATM,
275 "Can only print PQR atom lines as ATOM or HETATM records");
277 /* Check atom name */
278 GMX_RELEASE_ASSERT(atom_name != nullptr,
279 "Need atom information to print pqr");
281 /* Check residue name */
282 GMX_RELEASE_ASSERT(res_name != nullptr,
283 "Need residue information to print pqr");
285 /* Truncate integers so they fit */
286 atom_seq_number = atom_seq_number % 100000;
287 res_seq_number = res_seq_number % 10000;
290 "%s %d %s %s %c %d %8.3f %8.3f %8.3f %6.2f %6.2f\n",
304 void write_pdbfile_indexed(FILE *out, const char *title,
305 const t_atoms *atoms, const rvec x[],
306 int ePBC, const matrix box, char chainid,
307 int model_nr, int nindex, const int index[],
311 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
312 enum PDB_record type;
318 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : gmx::bromacs().c_str());
319 if (box && ( (norm2(box[XX]) != 0.0f) || (norm2(box[YY]) != 0.0f) || (norm2(box[ZZ]) != 0.0f) ) )
321 gmx_write_pdb_box(out, ePBC, box);
323 if (atoms->havePdbInfo)
325 /* Check whether any occupancies are set, in that case leave it as is,
326 * otherwise set them all to one
329 for (int ii = 0; (ii < nindex) && bOccup; ii++)
332 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
340 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
343 for (int ii = 0; ii < nindex; ii++)
346 int resind = atoms->atom[i].resind;
347 std::string resnm = *atoms->resinfo[resind].name;
348 std::string nm = *atoms->atomname[i];
350 /* rename HG12 to 2HG1, etc. */
351 nm = xlate_atomname_gmx2pdb(nm);
352 int resnr = atoms->resinfo[resind].nr;
353 unsigned char resic = atoms->resinfo[resind].ic;
361 ch = atoms->resinfo[resind].chainid;
370 resnr = resnr % 10000;
373 if (atoms->pdbinfo != nullptr)
375 pdbinfo = atoms->pdbinfo[i];
379 gmx_pdbinfo_init_default(&pdbinfo);
381 type = static_cast<enum PDB_record>(pdbinfo.type);
382 altloc = pdbinfo.altloc;
383 if (!isalnum(altloc))
387 occup = bOccup ? 1.0 : pdbinfo.occup;
391 gmx_fprintf_pdb_atomline(out,
400 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
403 atoms->atom[i].elem);
405 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
407 fprintf(out, "ANISOU%5d %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n",
408 (i+1)%100000, nm.c_str(), resnm.c_str(), ch, resnr,
409 (resic == '\0') ? ' ' : resic,
410 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
411 atoms->pdbinfo[i].uij[2], atoms->pdbinfo[i].uij[3],
412 atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
417 gmx_fprintf_pqr_atomline(out,
424 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
430 fprintf(out, "TER\n");
431 fprintf(out, "ENDMDL\n");
435 /* Write conect records */
436 for (int i = 0; (i < gc->nconect); i++)
438 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai+1, gc->conect[i].aj+1);
443 void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms, const rvec x[],
444 int ePBC, const matrix box, char chainid, int model_nr, gmx_conect conect)
448 snew(index, atoms->nr);
449 for (i = 0; i < atoms->nr; i++)
453 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
454 atoms->nr, index, conect, false);
458 static int line2type(const char *line)
463 for (k = 0; (k < 6); k++)
469 for (k = 0; (k < epdbNR); k++)
471 if (std::strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
480 static void read_anisou(char line[], int natom, t_atoms *atoms)
484 char anr[12], anm[12];
488 for (k = 0; (k < 5); k++, j++)
494 for (k = 0; (k < 4); k++, j++)
501 /* Strip off spaces */
504 /* Search backwards for number and name only */
505 atomnr = std::strtol(anr, nullptr, 10);
506 for (i = natom-1; (i >= 0); i--)
508 if ((std::strcmp(anm, *(atoms->atomname[i])) == 0) &&
509 (atomnr == atoms->pdbinfo[i].atomnr))
516 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n",
521 if (sscanf(line+29, "%d%d%d%d%d%d",
522 &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
523 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
524 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
527 atoms->pdbinfo[i].bAnisotropic = TRUE;
531 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
532 atoms->pdbinfo[i].bAnisotropic = FALSE;
537 void get_pdb_atomnumber(const t_atoms *atoms, AtomProperties *aps)
539 int i, atomnumber, len;
541 char anm[6], anm_copy[6];
547 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
549 for (i = 0; (i < atoms->nr); i++)
551 std::strcpy(anm, atoms->pdbinfo[i].atomnm);
552 std::strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
553 bool atomNumberSet = false;
555 if ((anm[0] != ' ') && ((len <= 2) || !std::isdigit(anm[2])))
558 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
560 atomnumber = gmx::roundToInt(eval);
561 atomNumberSet = true;
566 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
568 atomnumber = gmx::roundToInt(eval);
569 atomNumberSet = true;
576 while ((k < std::strlen(anm)) && (std::isspace(anm[k]) || std::isdigit(anm[k])))
580 anm_copy[0] = anm[k];
582 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
584 atomnumber = gmx::roundToInt(eval);
585 atomNumberSet = true;
591 atoms->atom[i].atomnumber = atomnumber;
592 buf = aps->elementFromAtomNumber(atomnumber);
595 fprintf(debug, "Atomnumber for atom '%s' is %d\n",
600 std::strncpy(atoms->atom[i].elem, buf.c_str(), 4);
604 static int read_atom(t_symtab *symtab,
605 const char line[], int type, int natom,
606 t_atoms *atoms, rvec x[], int chainnum, gmx_bool bChange)
611 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
612 char xc[12], yc[12], zc[12], occup[12], bfac[12];
615 int resnr, atomnumber;
617 if (natom >= atoms->nr)
619 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)",
625 for (k = 0; (k < 5); k++, j++)
632 for (k = 0; (k < 4); k++, j++)
637 std::strcpy(anm_copy, anm);
643 for (k = 0; (k < 4); k++, j++)
653 for (k = 0; (k < 4); k++, j++)
659 resnr = std::strtol(rnr, nullptr, 10);
663 /* X,Y,Z Coordinate */
664 for (k = 0; (k < 8); k++, j++)
669 for (k = 0; (k < 8); k++, j++)
674 for (k = 0; (k < 8); k++, j++)
681 for (k = 0; (k < 6); k++, j++)
688 for (k = 0; (k < 7); k++, j++)
698 for (k = 0; (k < 2); k++, j++)
707 atomn = &(atoms->atom[natom]);
709 atoms->resinfo[atoms->atom[natom-1].resind].nr != resnr ||
710 atoms->resinfo[atoms->atom[natom-1].resind].ic != resic ||
711 (strcmp(*atoms->resinfo[atoms->atom[natom-1].resind].name, resnm) != 0))
719 atomn->resind = atoms->atom[natom-1].resind + 1;
721 atoms->nres = atomn->resind + 1;
722 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
726 atomn->resind = atoms->atom[natom-1].resind;
730 xlate_atomname_pdb2gmx(anm);
732 atoms->atomname[natom] = put_symtab(symtab, anm);
735 atomn->atomnumber = atomnumber;
736 strncpy(atomn->elem, elem, 4);
738 x[natom][XX] = strtod(xc, nullptr)*0.1;
739 x[natom][YY] = strtod(yc, nullptr)*0.1;
740 x[natom][ZZ] = strtod(zc, nullptr)*0.1;
743 atoms->pdbinfo[natom].type = type;
744 atoms->pdbinfo[natom].atomnr = strtol(anr, nullptr, 10);
745 atoms->pdbinfo[natom].altloc = altloc;
746 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
747 atoms->pdbinfo[natom].bfac = strtod(bfac, nullptr);
748 atoms->pdbinfo[natom].occup = strtod(occup, nullptr);
755 gmx_bool is_hydrogen(const char *nm)
759 std::strcpy(buf, nm);
766 else if ((std::isdigit(buf[0])) && (buf[1] == 'H'))
773 gmx_bool is_dummymass(const char *nm)
777 std::strcpy(buf, nm);
780 return (buf[0] == 'M') && (std::isdigit(buf[strlen(buf)-1]) != 0);
783 static void gmx_conect_addline(gmx_conect_t *con, char *line)
787 std::string form2 = "%%*s";
788 std::string format = form2 + "%%d";
789 if (sscanf(line, format.c_str(), &ai) == 1)
794 format = form2 + "%%d";
795 n = sscanf(line, format.c_str(), &aj);
798 srenew(con->conect, ++con->nconect);
799 con->conect[con->nconect-1].ai = ai-1;
800 con->conect[con->nconect-1].aj = aj-1;
807 void gmx_conect_dump(FILE *fp, gmx_conect conect)
809 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
812 for (i = 0; (i < gc->nconect); i++)
814 fprintf(fp, "%6s%5d%5d\n", "CONECT",
815 gc->conect[i].ai+1, gc->conect[i].aj+1);
819 gmx_conect gmx_conect_init()
828 void gmx_conect_done(gmx_conect conect)
830 gmx_conect_t *gc = conect;
835 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
837 gmx_conect_t *gc = conect;
843 for (i = 0; (i < gc->nconect); i++)
845 if (((gc->conect[i].ai == ai) &&
846 (gc->conect[i].aj == aj)) ||
847 ((gc->conect[i].aj == ai) &&
848 (gc->conect[i].ai == aj)))
856 void gmx_conect_add(gmx_conect conect, int ai, int aj)
858 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
863 if (!gmx_conect_exist(conect, ai, aj))
865 srenew(gc->conect, ++gc->nconect);
866 gc->conect[gc->nconect-1].ai = ai;
867 gc->conect[gc->nconect-1].aj = aj;
871 int read_pdbfile(FILE *in, char *title, int *model_nr,
872 t_atoms *atoms, t_symtab *symtab, rvec x[], int *ePBC,
873 matrix box, gmx_bool bChange, gmx_conect conect)
875 gmx_conect_t *gc = conect;
877 gmx_bool bConnWarn = FALSE;
882 gmx_bool bStop = FALSE;
886 /* Only assume pbc when there is a CRYST1 entry */
894 atoms->haveMass = FALSE;
895 atoms->haveCharge = FALSE;
896 atoms->haveType = FALSE;
897 atoms->haveBState = FALSE;
898 atoms->havePdbInfo = (atoms->pdbinfo != nullptr);
904 while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
906 line_type = line2type(line);
912 natom = read_atom(symtab, line, line_type, natom, atoms, x, chainnum, bChange);
916 if (atoms->havePdbInfo)
918 read_anisou(line, natom, atoms);
923 read_cryst1(line, ePBC, box);
928 if (std::strlen(line) > 6)
931 /* skip HEADER or TITLE and spaces */
940 /* truncate after title */
941 d = std::strstr(c, " ");
946 if (std::strlen(c) > 0)
948 std::strcpy(title, c);
954 if ((!std::strstr(line, ": ")) || (std::strstr(line+6, "MOLECULE:")))
956 if (!(c = std::strstr(line+6, "MOLECULE:")) )
960 /* skip 'MOLECULE:' and spaces */
969 /* truncate after title */
973 while ( (d[-1] == ';') && d > c)
983 std::strcat(title, "; ");
984 std::strcat(title, c);
988 std::strcpy(title, c);
1002 sscanf(line, "%*s%d", model_nr);
1012 gmx_conect_addline(gc, line);
1014 else if (!bConnWarn)
1016 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
1029 void get_pdb_coordnum(FILE *in, int *natoms)
1034 while (fgets2(line, STRLEN, in))
1036 if (std::strncmp(line, "ENDMDL", 6) == 0)
1040 if ((std::strncmp(line, "ATOM ", 6) == 0) || (std::strncmp(line, "HETATM", 6) == 0))
1047 void gmx_pdb_read_conf(const char *infile,
1048 t_symtab *symtab, char **name, t_atoms *atoms,
1049 rvec x[], int *ePBC, matrix box)
1051 FILE *in = gmx_fio_fopen(infile, "r");
1053 read_pdbfile(in, title, nullptr, atoms, symtab, x, ePBC, box, TRUE, nullptr);
1054 if (name != nullptr)
1056 *name = gmx_strdup(title);
1061 gmx_conect gmx_conect_generate(const t_topology *top)
1066 /* Fill the conect records */
1067 gc = gmx_conect_init();
1069 for (f = 0; (f < F_NRE); f++)
1073 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms+1)
1075 gmx_conect_add(gc, top->idef.il[f].iatoms[i+1],
1076 top->idef.il[f].iatoms[i+2]);
1084 gmx_fprintf_pdb_atomline(FILE * fp,
1085 enum PDB_record record,
1086 int atom_seq_number,
1087 const char * atom_name,
1088 char alternate_location,
1089 const char * res_name,
1092 char res_insertion_code,
1098 const char * element)
1100 char tmp_atomname[6], tmp_resname[6];
1101 gmx_bool start_name_in_col13;
1104 if (record != epdbATOM && record != epdbHETATM)
1106 gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
1109 /* Format atom name */
1110 if (atom_name != nullptr)
1112 /* If the atom name is an element name with two chars, it should start already in column 13.
1113 * Otherwise it should start in column 14, unless the name length is 4 chars.
1115 if ( (element != nullptr) && (std::strlen(element) >= 2) && (gmx_strncasecmp(atom_name, element, 2) == 0) )
1117 start_name_in_col13 = TRUE;
1121 start_name_in_col13 = (std::strlen(atom_name) >= 4);
1123 snprintf(tmp_atomname, sizeof(tmp_atomname), start_name_in_col13 ? "" : " ");
1124 std::strncat(tmp_atomname, atom_name, 4);
1125 tmp_atomname[5] = '\0';
1129 tmp_atomname[0] = '\0';
1132 /* Format residue name */
1133 std::strncpy(tmp_resname, (res_name != nullptr) ? res_name : "", 4);
1134 /* Make sure the string is terminated if strlen was > 4 */
1135 tmp_resname[4] = '\0';
1136 /* String is properly terminated, so now we can use strcat. By adding a
1137 * space we can write it right-justified, and if the original name was
1138 * three characters or less there will be a space added on the right side.
1140 std::strcat(tmp_resname, " ");
1142 /* Truncate integers so they fit */
1143 atom_seq_number = atom_seq_number % 100000;
1144 res_seq_number = res_seq_number % 10000;
1147 "%-6s%5d %-4.4s%c%4.4s%c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f %2s\n",
1159 (element != nullptr) ? element : "");