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,2020, 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"
69 typedef struct gmx_conect_t
73 gmx_conection_t* conect;
76 static const char* pdbtp[epdbNR] = { "ATOM ", "HETATM", "ANISOU", "CRYST1", "COMPND", "MODEL",
77 "ENDMDL", "TER", "HEADER", "TITLE", "REMARK", "CONECT" };
79 #define REMARK_SIM_BOX "REMARK THIS IS A SIMULATION BOX"
81 static void xlate_atomname_pdb2gmx(char* name)
86 length = std::strlen(name);
87 if (length > 3 && std::isdigit(name[0]))
90 for (i = 1; i < length; i++)
92 name[i - 1] = name[i];
94 name[length - 1] = temp;
98 // Deliberately taking a copy of name to return it later
99 static std::string xlate_atomname_gmx2pdb(std::string name)
101 size_t length = name.size();
102 if (length > 3 && std::isdigit(name[length - 1]))
104 char temp = name[length - 1];
105 for (size_t i = length - 1; i > 0; --i)
107 name[i] = name[i - 1];
115 void gmx_write_pdb_box(FILE* out, int ePBC, const matrix box)
117 real alpha, beta, gamma;
121 ePBC = guess_ePBC(box);
124 if (ePBC == epbcNONE)
129 if (norm2(box[YY]) * norm2(box[ZZ]) != 0)
131 alpha = RAD2DEG * gmx_angle(box[YY], box[ZZ]);
137 if (norm2(box[XX]) * norm2(box[ZZ]) != 0)
139 beta = RAD2DEG * gmx_angle(box[XX], box[ZZ]);
145 if (norm2(box[XX]) * norm2(box[YY]) != 0)
147 gamma = RAD2DEG * gmx_angle(box[XX], box[YY]);
153 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
154 if (ePBC != epbcSCREW)
156 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n", 10 * norm(box[XX]),
157 10 * norm(box[YY]), 10 * norm(box[ZZ]), alpha, beta, gamma, "P 1", 1);
161 /* Double the a-vector length and write the correct space group */
162 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n", 20 * norm(box[XX]),
163 10 * norm(box[YY]), 10 * norm(box[ZZ]), alpha, beta, gamma, "P 21 1 1", 1);
167 static void read_cryst1(char* line, int* ePBC, matrix box)
170 char sa[12], sb[12], sc[12], sg[SG_SIZE + 1], ident;
171 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
172 int syma, symb, symc;
175 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
178 if (strlen(line) >= 55)
180 strncpy(sg, line + 55, SG_SIZE);
186 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
187 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
189 fc = strtod(sc, nullptr) * 0.1;
190 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
192 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
194 ePBC_file = epbcSCREW;
204 fa = strtod(sa, nullptr) * 0.1;
205 fb = strtod(sb, nullptr) * 0.1;
206 fc = strtod(sc, nullptr) * 0.1;
207 if (ePBC_file == epbcSCREW)
213 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
217 cosa = std::cos(alpha * DEG2RAD);
225 cosb = std::cos(beta * DEG2RAD);
233 cosg = std::cos(gamma * DEG2RAD);
234 sing = std::sin(gamma * DEG2RAD);
241 box[YY][XX] = fb * cosg;
242 box[YY][YY] = fb * sing;
243 box[ZZ][XX] = fc * cosb;
244 box[ZZ][YY] = fc * (cosa - cosb * cosg) / sing;
245 box[ZZ][ZZ] = std::sqrt(fc * fc - box[ZZ][XX] * box[ZZ][XX] - box[ZZ][YY] * box[ZZ][YY]);
255 static int gmx_fprintf_pqr_atomline(FILE* fp,
256 enum PDB_record record,
258 const char* atom_name,
259 const char* res_name,
268 GMX_RELEASE_ASSERT(record == epdbATOM || record == epdbHETATM,
269 "Can only print PQR atom lines as ATOM or HETATM records");
271 /* Check atom name */
272 GMX_RELEASE_ASSERT(atom_name != nullptr, "Need atom information to print pqr");
274 /* Check residue name */
275 GMX_RELEASE_ASSERT(res_name != nullptr, "Need residue information to print pqr");
277 /* Truncate integers so they fit */
278 atom_seq_number = atom_seq_number % 100000;
279 res_seq_number = res_seq_number % 10000;
281 int n = fprintf(fp, "%-6s%5d %-4.4s%4.4s%c%4d %8.3f %8.3f %8.3f %6.2f %6.2f\n", pdbtp[record],
282 atom_seq_number, atom_name, res_name, chain_id, res_seq_number, x, y, z,
283 occupancy, b_factor);
288 void write_pdbfile_indexed(FILE* out,
290 const t_atoms* atoms,
301 gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
302 enum PDB_record type;
308 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : gmx::bromacs().c_str());
309 if (box && ((norm2(box[XX]) != 0.0F) || (norm2(box[YY]) != 0.0F) || (norm2(box[ZZ]) != 0.0F)))
311 gmx_write_pdb_box(out, ePBC, box);
313 if (atoms->havePdbInfo)
315 /* Check whether any occupancies are set, in that case leave it as is,
316 * otherwise set them all to one
319 for (int ii = 0; (ii < nindex) && bOccup; ii++)
322 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
330 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
333 for (int ii = 0; ii < nindex; ii++)
336 int resind = atoms->atom[i].resind;
337 std::string resnm = *atoms->resinfo[resind].name;
338 std::string nm = *atoms->atomname[i];
340 /* rename HG12 to 2HG1, etc. */
341 nm = xlate_atomname_gmx2pdb(nm);
342 int resnr = atoms->resinfo[resind].nr;
343 unsigned char resic = atoms->resinfo[resind].ic;
351 ch = atoms->resinfo[resind].chainid;
360 resnr = resnr % 10000;
363 if (atoms->pdbinfo != nullptr)
365 pdbinfo = atoms->pdbinfo[i];
369 gmx_pdbinfo_init_default(&pdbinfo);
371 type = static_cast<enum PDB_record>(pdbinfo.type);
372 altloc = pdbinfo.altloc;
373 if (!isalnum(altloc))
377 occup = bOccup ? 1.0 : pdbinfo.occup;
381 gmx_fprintf_pdb_atomline(out, type, i + 1, nm.c_str(), altloc, resnm.c_str(), ch, resnr,
382 resic, 10 * x[i][XX], 10 * x[i][YY], 10 * x[i][ZZ], occup,
383 bfac, atoms->atom[i].elem);
385 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
387 fprintf(out, "ANISOU%5d %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n", (i + 1) % 100000,
388 nm.c_str(), resnm.c_str(), ch, resnr, (resic == '\0') ? ' ' : resic,
389 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1], atoms->pdbinfo[i].uij[2],
390 atoms->pdbinfo[i].uij[3], atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
395 gmx_fprintf_pqr_atomline(out, type, i + 1, nm.c_str(), resnm.c_str(), ch, resnr,
396 10 * x[i][XX], 10 * x[i][YY], 10 * x[i][ZZ], occup, bfac);
400 fprintf(out, "TER\n");
401 fprintf(out, "ENDMDL\n");
405 /* Write conect records */
406 for (int i = 0; (i < gc->nconect); i++)
408 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai + 1, gc->conect[i].aj + 1);
413 void write_pdbfile(FILE* out,
415 const t_atoms* atoms,
425 snew(index, atoms->nr);
426 for (i = 0; i < atoms->nr; i++)
430 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr, atoms->nr, index,
435 static int line2type(const char* line)
440 for (k = 0; (k < 6); k++)
446 for (k = 0; (k < epdbNR); k++)
448 if (std::strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
457 static void read_anisou(char line[], int natom, t_atoms* atoms)
461 char anr[12], anm[12];
465 for (k = 0; (k < 5); k++, j++)
471 for (k = 0; (k < 4); k++, j++)
478 /* Strip off spaces */
481 /* Search backwards for number and name only */
482 atomnr = std::strtol(anr, nullptr, 10);
483 for (i = natom - 1; (i >= 0); i--)
485 if ((std::strcmp(anm, *(atoms->atomname[i])) == 0) && (atomnr == atoms->pdbinfo[i].atomnr))
492 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n", anm, atomnr);
496 if (sscanf(line + 29, "%d%d%d%d%d%d", &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
497 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
498 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
501 atoms->pdbinfo[i].bAnisotropic = TRUE;
505 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
506 atoms->pdbinfo[i].bAnisotropic = FALSE;
511 void get_pdb_atomnumber(const t_atoms* atoms, AtomProperties* aps)
513 int i, atomnumber, len;
515 char anm[6], anm_copy[6];
521 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
523 for (i = 0; (i < atoms->nr); i++)
525 std::strcpy(anm, atoms->pdbinfo[i].atomnm);
526 std::strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
527 bool atomNumberSet = false;
529 if ((anm[0] != ' ') && ((len <= 2) || !std::isdigit(anm[2])))
532 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
534 atomnumber = gmx::roundToInt(eval);
535 atomNumberSet = true;
540 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
542 atomnumber = gmx::roundToInt(eval);
543 atomNumberSet = true;
550 while ((k < std::strlen(anm)) && (std::isspace(anm[k]) || std::isdigit(anm[k])))
554 anm_copy[0] = anm[k];
556 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
558 atomnumber = gmx::roundToInt(eval);
559 atomNumberSet = true;
565 atoms->atom[i].atomnumber = atomnumber;
566 buf = aps->elementFromAtomNumber(atomnumber);
569 fprintf(debug, "Atomnumber for atom '%s' is %d\n", anm, atomnumber);
573 std::strncpy(atoms->atom[i].elem, buf.c_str(), 4);
578 read_atom(t_symtab* symtab, const char line[], int type, int natom, t_atoms* atoms, rvec x[], int chainnum, gmx_bool bChange)
583 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
584 char xc[12], yc[12], zc[12], occup[12], bfac[12];
587 int resnr, atomnumber;
589 if (natom >= atoms->nr)
591 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)", natom + 1, atoms->nr);
596 for (k = 0; (k < 5); k++, j++)
603 for (k = 0; (k < 4); k++, j++)
608 std::strcpy(anm_copy, anm);
614 for (k = 0; (k < 4); k++, j++)
624 for (k = 0; (k < 4); k++, j++)
630 resnr = std::strtol(rnr, nullptr, 10);
634 /* X,Y,Z Coordinate */
635 for (k = 0; (k < 8); k++, j++)
640 for (k = 0; (k < 8); k++, j++)
645 for (k = 0; (k < 8); k++, j++)
652 for (k = 0; (k < 6); k++, j++)
659 for (k = 0; (k < 7); k++, j++)
669 for (k = 0; (k < 2); k++, j++)
678 atomn = &(atoms->atom[natom]);
679 if ((natom == 0) || atoms->resinfo[atoms->atom[natom - 1].resind].nr != resnr
680 || atoms->resinfo[atoms->atom[natom - 1].resind].ic != resic
681 || (strcmp(*atoms->resinfo[atoms->atom[natom - 1].resind].name, resnm) != 0))
689 atomn->resind = atoms->atom[natom - 1].resind + 1;
691 atoms->nres = atomn->resind + 1;
692 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
696 atomn->resind = atoms->atom[natom - 1].resind;
700 xlate_atomname_pdb2gmx(anm);
702 atoms->atomname[natom] = put_symtab(symtab, anm);
705 atomn->atomnumber = atomnumber;
706 strncpy(atomn->elem, elem, 4);
708 x[natom][XX] = strtod(xc, nullptr) * 0.1;
709 x[natom][YY] = strtod(yc, nullptr) * 0.1;
710 x[natom][ZZ] = strtod(zc, nullptr) * 0.1;
713 atoms->pdbinfo[natom].type = type;
714 atoms->pdbinfo[natom].atomnr = strtol(anr, nullptr, 10);
715 atoms->pdbinfo[natom].altloc = altloc;
716 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
717 atoms->pdbinfo[natom].bfac = strtod(bfac, nullptr);
718 atoms->pdbinfo[natom].occup = strtod(occup, nullptr);
725 gmx_bool is_hydrogen(const char* nm)
729 std::strcpy(buf, nm);
736 else if ((std::isdigit(buf[0])) && (buf[1] == 'H'))
743 gmx_bool is_dummymass(const char* nm)
747 std::strcpy(buf, nm);
750 return (buf[0] == 'M') && (std::isdigit(buf[strlen(buf) - 1]) != 0);
753 static void gmx_conect_addline(gmx_conect_t* con, char* line)
757 std::string form2 = "%*s";
758 std::string format = form2 + "%d";
759 if (sscanf(line, format.c_str(), &ai) == 1)
764 format = form2 + "%d";
765 n = sscanf(line, format.c_str(), &aj);
768 gmx_conect_add(con, ai - 1, aj - 1); /* to prevent duplicated records */
774 void gmx_conect_dump(FILE* fp, gmx_conect conect)
776 gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
779 for (i = 0; (i < gc->nconect); i++)
781 fprintf(fp, "%6s%5d%5d\n", "CONECT", gc->conect[i].ai + 1, gc->conect[i].aj + 1);
785 gmx_conect gmx_conect_init()
794 void gmx_conect_done(gmx_conect conect)
796 gmx_conect_t* gc = conect;
801 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
803 gmx_conect_t* gc = conect;
809 for (i = 0; (i < gc->nconect); i++)
811 if (((gc->conect[i].ai == ai) && (gc->conect[i].aj == aj))
812 || ((gc->conect[i].aj == ai) && (gc->conect[i].ai == aj)))
820 void gmx_conect_add(gmx_conect conect, int ai, int aj)
822 gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
827 if (!gmx_conect_exist(conect, ai, aj))
829 srenew(gc->conect, ++gc->nconect);
830 gc->conect[gc->nconect - 1].ai = ai;
831 gc->conect[gc->nconect - 1].aj = aj;
835 int read_pdbfile(FILE* in,
846 gmx_conect_t* gc = conect;
848 gmx_bool bConnWarn = FALSE;
849 char line[STRLEN + 1];
853 gmx_bool bStop = FALSE;
857 /* Only assume pbc when there is a CRYST1 entry */
865 atoms->haveMass = FALSE;
866 atoms->haveCharge = FALSE;
867 atoms->haveType = FALSE;
868 atoms->haveBState = FALSE;
869 atoms->havePdbInfo = (atoms->pdbinfo != nullptr);
875 while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
877 line_type = line2type(line);
883 natom = read_atom(symtab, line, line_type, natom, atoms, x, chainnum, bChange);
887 if (atoms->havePdbInfo)
889 read_anisou(line, natom, atoms);
893 case epdbCRYST1: read_cryst1(line, ePBC, box); break;
897 if (std::strlen(line) > 6)
900 /* skip HEADER or TITLE and spaces */
909 /* truncate after title */
910 d = std::strstr(c, " ");
915 if (std::strlen(c) > 0)
917 std::strcpy(title, c);
923 if ((!std::strstr(line, ": ")) || (std::strstr(line + 6, "MOLECULE:")))
925 if (!(c = std::strstr(line + 6, "MOLECULE:")))
929 /* skip 'MOLECULE:' and spaces */
938 /* truncate after title */
942 while ((d[-1] == ';') && d > c)
952 std::strcat(title, "; ");
953 std::strcat(title, c);
957 std::strcpy(title, c);
964 case epdbTER: chainnum++; break;
969 sscanf(line, "%*s%d", model_nr);
973 case epdbENDMDL: bStop = TRUE; break;
977 gmx_conect_addline(gc, line);
981 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
993 void get_pdb_coordnum(FILE* in, int* natoms)
998 while (fgets2(line, STRLEN, in))
1000 if (std::strncmp(line, "ENDMDL", 6) == 0)
1004 if ((std::strncmp(line, "ATOM ", 6) == 0) || (std::strncmp(line, "HETATM", 6) == 0))
1011 void gmx_pdb_read_conf(const char* infile, t_symtab* symtab, char** name, t_atoms* atoms, rvec x[], int* ePBC, matrix box)
1013 FILE* in = gmx_fio_fopen(infile, "r");
1015 read_pdbfile(in, title, nullptr, atoms, symtab, x, ePBC, box, TRUE, nullptr);
1016 if (name != nullptr)
1018 *name = gmx_strdup(title);
1023 gmx_conect gmx_conect_generate(const t_topology* top)
1028 /* Fill the conect records */
1029 gc = gmx_conect_init();
1031 for (f = 0; (f < F_NRE); f++)
1035 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms + 1)
1037 gmx_conect_add(gc, top->idef.il[f].iatoms[i + 1], top->idef.il[f].iatoms[i + 2]);
1044 int gmx_fprintf_pdb_atomline(FILE* fp,
1045 enum PDB_record record,
1046 int atom_seq_number,
1047 const char* atom_name,
1048 char alternate_location,
1049 const char* res_name,
1052 char res_insertion_code,
1058 const char* element)
1060 char tmp_atomname[6], tmp_resname[6];
1061 gmx_bool start_name_in_col13;
1064 if (record != epdbATOM && record != epdbHETATM)
1066 gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
1069 /* Format atom name */
1070 if (atom_name != nullptr)
1072 /* If the atom name is an element name with two chars, it should start already in column 13.
1073 * Otherwise it should start in column 14, unless the name length is 4 chars.
1075 if ((element != nullptr) && (std::strlen(element) >= 2)
1076 && (gmx_strncasecmp(atom_name, element, 2) == 0))
1078 start_name_in_col13 = TRUE;
1082 start_name_in_col13 = (std::strlen(atom_name) >= 4);
1084 snprintf(tmp_atomname, sizeof(tmp_atomname), start_name_in_col13 ? "" : " ");
1085 std::strncat(tmp_atomname, atom_name, 4);
1086 tmp_atomname[5] = '\0';
1090 tmp_atomname[0] = '\0';
1093 /* Format residue name */
1094 std::strncpy(tmp_resname, (res_name != nullptr) ? res_name : "", 4);
1095 /* Make sure the string is terminated if strlen was > 4 */
1096 tmp_resname[4] = '\0';
1097 /* String is properly terminated, so now we can use strcat. By adding a
1098 * space we can write it right-justified, and if the original name was
1099 * three characters or less there will be a space added on the right side.
1101 std::strcat(tmp_resname, " ");
1103 /* Truncate integers so they fit */
1104 atom_seq_number = atom_seq_number % 100000;
1105 res_seq_number = res_seq_number % 10000;
1107 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",
1108 pdbtp[record], atom_seq_number, tmp_atomname, alternate_location, tmp_resname,
1109 chain_id, res_seq_number, res_insertion_code, x, y, z, occupancy, b_factor,
1110 (element != nullptr) ? element : "");