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, 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",
81 #define REMARK_SIM_BOX "REMARK THIS IS A SIMULATION BOX"
83 static void xlate_atomname_pdb2gmx(char *name)
88 length = std::strlen(name);
89 if (length > 3 && std::isdigit(name[0]))
92 for (i = 1; i < length; i++)
96 name[length-1] = temp;
100 // Deliberately taking a copy of name to return it later
101 static std::string xlate_atomname_gmx2pdb(std::string name)
103 size_t length = name.size();
104 if (length > 3 && std::isdigit(name[length-1]))
106 char temp = name[length-1];
107 for (size_t i = length-1; i > 0; --i)
117 void gmx_write_pdb_box(FILE *out, int ePBC, const matrix box)
119 real alpha, beta, gamma;
123 ePBC = guess_ePBC(box);
126 if (ePBC == epbcNONE)
131 if (norm2(box[YY])*norm2(box[ZZ]) != 0)
133 alpha = RAD2DEG*gmx_angle(box[YY], box[ZZ]);
139 if (norm2(box[XX])*norm2(box[ZZ]) != 0)
141 beta = RAD2DEG*gmx_angle(box[XX], box[ZZ]);
147 if (norm2(box[XX])*norm2(box[YY]) != 0)
149 gamma = RAD2DEG*gmx_angle(box[XX], box[YY]);
155 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
156 if (ePBC != epbcSCREW)
158 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
159 10*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
160 alpha, beta, gamma, "P 1", 1);
164 /* Double the a-vector length and write the correct space group */
165 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
166 20*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
167 alpha, beta, gamma, "P 21 1 1", 1);
172 static void read_cryst1(char *line, int *ePBC, matrix box)
175 char sa[12], sb[12], sc[12], sg[SG_SIZE+1], ident;
176 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
177 int syma, symb, symc;
180 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
183 if (strlen(line) >= 55)
185 strncpy(sg, line+55, SG_SIZE);
191 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
192 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
194 fc = strtod(sc, nullptr)*0.1;
195 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
197 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
199 ePBC_file = epbcSCREW;
209 fa = strtod(sa, nullptr)*0.1;
210 fb = strtod(sb, nullptr)*0.1;
211 fc = strtod(sc, nullptr)*0.1;
212 if (ePBC_file == epbcSCREW)
218 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
222 cosa = std::cos(alpha*DEG2RAD);
230 cosb = std::cos(beta*DEG2RAD);
238 cosg = std::cos(gamma*DEG2RAD);
239 sing = std::sin(gamma*DEG2RAD);
246 box[YY][XX] = fb*cosg;
247 box[YY][YY] = fb*sing;
248 box[ZZ][XX] = fc*cosb;
249 box[ZZ][YY] = fc*(cosa - cosb*cosg)/sing;
250 box[ZZ][ZZ] = std::sqrt(fc*fc
251 - box[ZZ][XX]*box[ZZ][XX] - box[ZZ][YY]*box[ZZ][YY]);
262 gmx_fprintf_pqr_atomline(FILE * fp,
263 enum PDB_record record,
265 const char * atom_name,
266 const char * res_name,
275 GMX_RELEASE_ASSERT(record == epdbATOM || record == epdbHETATM,
276 "Can only print PQR atom lines as ATOM or HETATM records");
278 /* Check atom name */
279 GMX_RELEASE_ASSERT(atom_name != nullptr,
280 "Need atom information to print pqr");
282 /* Check residue name */
283 GMX_RELEASE_ASSERT(res_name != nullptr,
284 "Need residue information to print pqr");
286 /* Truncate integers so they fit */
287 atom_seq_number = atom_seq_number % 100000;
288 res_seq_number = res_seq_number % 10000;
291 "%s %d %s %s %c %d %8.3f %8.3f %8.3f %6.2f %6.2f\n",
305 void write_pdbfile_indexed(FILE *out, const char *title,
306 const t_atoms *atoms, const rvec x[],
307 int ePBC, const matrix box, char chainid,
308 int model_nr, int nindex, const int index[],
309 gmx_conect conect, gmx_bool bTerSepChains,
312 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
315 enum PDB_record type;
316 unsigned char resic, ch;
320 int chainnum, lastchainnum;
321 gmx_residuetype_t*rt;
322 const char *p_restype;
323 const char *p_lastrestype;
325 gmx_residuetype_init(&rt);
327 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : gmx::bromacs().c_str());
328 if (box && ( (norm2(box[XX]) != 0.0f) || (norm2(box[YY]) != 0.0f) || (norm2(box[ZZ]) != 0.0f) ) )
330 gmx_write_pdb_box(out, ePBC, box);
332 if (atoms->havePdbInfo)
334 /* Check whether any occupancies are set, in that case leave it as is,
335 * otherwise set them all to one
338 for (ii = 0; (ii < nindex) && bOccup; ii++)
341 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
349 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
354 for (ii = 0; ii < nindex; ii++)
357 resind = atoms->atom[i].resind;
358 chainnum = atoms->resinfo[resind].chainnum;
359 p_lastrestype = p_restype;
360 gmx_residuetype_get_type(rt, *atoms->resinfo[resind].name, &p_restype);
362 /* Add a TER record if we changed chain, and if either the previous or this chain is protein/DNA/RNA. */
363 if (bTerSepChains && ii > 0 && chainnum != lastchainnum)
365 /* Only add TER if the previous chain contained protein/DNA/RNA. */
366 if (gmx_residuetype_is_protein(rt, p_lastrestype) || gmx_residuetype_is_dna(rt, p_lastrestype) || gmx_residuetype_is_rna(rt, p_lastrestype))
368 fprintf(out, "TER\n");
370 lastchainnum = chainnum;
373 std::string resnm = *atoms->resinfo[resind].name;
374 std::string nm = *atoms->atomname[i];
376 /* rename HG12 to 2HG1, etc. */
377 nm = xlate_atomname_gmx2pdb(nm);
378 resnr = atoms->resinfo[resind].nr;
379 resic = atoms->resinfo[resind].ic;
386 ch = atoms->resinfo[resind].chainid;
395 resnr = resnr % 10000;
398 if (atoms->pdbinfo != nullptr)
400 pdbinfo = atoms->pdbinfo[i];
404 gmx_pdbinfo_init_default(&pdbinfo);
406 type = static_cast<enum PDB_record>(pdbinfo.type);
407 altloc = pdbinfo.altloc;
408 if (!isalnum(altloc))
412 occup = bOccup ? 1.0 : pdbinfo.occup;
416 gmx_fprintf_pdb_atomline(out,
425 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
428 atoms->atom[i].elem);
430 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
432 fprintf(out, "ANISOU%5d %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n",
433 (i+1)%100000, nm.c_str(), resnm.c_str(), ch, resnr,
434 (resic == '\0') ? ' ' : resic,
435 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
436 atoms->pdbinfo[i].uij[2], atoms->pdbinfo[i].uij[3],
437 atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
442 gmx_fprintf_pqr_atomline(out,
449 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
455 fprintf(out, "TER\n");
456 fprintf(out, "ENDMDL\n");
460 /* Write conect records */
461 for (i = 0; (i < gc->nconect); i++)
463 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai+1, gc->conect[i].aj+1);
467 gmx_residuetype_destroy(rt);
470 void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms, const rvec x[],
471 int ePBC, const matrix box, char chainid, int model_nr, gmx_conect conect, gmx_bool bTerSepChains)
475 snew(index, atoms->nr);
476 for (i = 0; i < atoms->nr; i++)
480 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
481 atoms->nr, index, conect, bTerSepChains, false);
485 static int line2type(const char *line)
490 for (k = 0; (k < 6); k++)
496 for (k = 0; (k < epdbNR); k++)
498 if (std::strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
507 static void read_anisou(char line[], int natom, t_atoms *atoms)
511 char anr[12], anm[12];
515 for (k = 0; (k < 5); k++, j++)
521 for (k = 0; (k < 4); k++, j++)
528 /* Strip off spaces */
531 /* Search backwards for number and name only */
532 atomnr = std::strtol(anr, nullptr, 10);
533 for (i = natom-1; (i >= 0); i--)
535 if ((std::strcmp(anm, *(atoms->atomname[i])) == 0) &&
536 (atomnr == atoms->pdbinfo[i].atomnr))
543 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n",
548 if (sscanf(line+29, "%d%d%d%d%d%d",
549 &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
550 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
551 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
554 atoms->pdbinfo[i].bAnisotropic = TRUE;
558 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
559 atoms->pdbinfo[i].bAnisotropic = FALSE;
564 void get_pdb_atomnumber(const t_atoms *atoms, gmx_atomprop_t aps)
566 int i, atomnumber, len;
568 char anm[6], anm_copy[6], *ptr;
574 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
576 for (i = 0; (i < atoms->nr); i++)
578 std::strcpy(anm, atoms->pdbinfo[i].atomnm);
579 std::strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
580 bool atomNumberSet = false;
582 if ((anm[0] != ' ') && ((len <= 2) || !std::isdigit(anm[2])))
585 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
587 atomnumber = gmx::roundToInt(eval);
588 atomNumberSet = true;
593 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
595 atomnumber = gmx::roundToInt(eval);
596 atomNumberSet = true;
603 while ((k < std::strlen(anm)) && (std::isspace(anm[k]) || std::isdigit(anm[k])))
607 anm_copy[0] = anm[k];
609 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
611 atomnumber = gmx::roundToInt(eval);
612 atomNumberSet = true;
617 atoms->atom[i].atomnumber = atomnumber;
618 ptr = gmx_atomprop_element(aps, atomnumber);
621 fprintf(debug, "Atomnumber for atom '%s' is %d\n",
629 std::strncpy(atoms->atom[i].elem, ptr == nullptr ? "" : ptr, 4);
633 static int read_atom(t_symtab *symtab,
634 const char line[], int type, int natom,
635 t_atoms *atoms, rvec x[], int chainnum, gmx_bool bChange)
640 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
641 char xc[12], yc[12], zc[12], occup[12], bfac[12];
644 int resnr, atomnumber;
646 if (natom >= atoms->nr)
648 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)",
654 for (k = 0; (k < 5); k++, j++)
661 for (k = 0; (k < 4); k++, j++)
666 std::strcpy(anm_copy, anm);
672 for (k = 0; (k < 4); k++, j++)
682 for (k = 0; (k < 4); k++, j++)
688 resnr = std::strtol(rnr, nullptr, 10);
692 /* X,Y,Z Coordinate */
693 for (k = 0; (k < 8); k++, j++)
698 for (k = 0; (k < 8); k++, j++)
703 for (k = 0; (k < 8); k++, j++)
710 for (k = 0; (k < 6); k++, j++)
717 for (k = 0; (k < 7); k++, j++)
727 for (k = 0; (k < 2); k++, j++)
736 atomn = &(atoms->atom[natom]);
738 atoms->resinfo[atoms->atom[natom-1].resind].nr != resnr ||
739 atoms->resinfo[atoms->atom[natom-1].resind].ic != resic ||
740 (strcmp(*atoms->resinfo[atoms->atom[natom-1].resind].name, resnm) != 0))
748 atomn->resind = atoms->atom[natom-1].resind + 1;
750 atoms->nres = atomn->resind + 1;
751 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
755 atomn->resind = atoms->atom[natom-1].resind;
759 xlate_atomname_pdb2gmx(anm);
761 atoms->atomname[natom] = put_symtab(symtab, anm);
764 atomn->atomnumber = atomnumber;
765 strncpy(atomn->elem, elem, 4);
767 x[natom][XX] = strtod(xc, nullptr)*0.1;
768 x[natom][YY] = strtod(yc, nullptr)*0.1;
769 x[natom][ZZ] = strtod(zc, nullptr)*0.1;
772 atoms->pdbinfo[natom].type = type;
773 atoms->pdbinfo[natom].atomnr = strtol(anr, nullptr, 10);
774 atoms->pdbinfo[natom].altloc = altloc;
775 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
776 atoms->pdbinfo[natom].bfac = strtod(bfac, nullptr);
777 atoms->pdbinfo[natom].occup = strtod(occup, nullptr);
784 gmx_bool is_hydrogen(const char *nm)
788 std::strcpy(buf, nm);
795 else if ((std::isdigit(buf[0])) && (buf[1] == 'H'))
802 gmx_bool is_dummymass(const char *nm)
806 std::strcpy(buf, nm);
809 return (buf[0] == 'M') && (std::isdigit(buf[strlen(buf)-1]) != 0);
812 static void gmx_conect_addline(gmx_conect_t *con, char *line)
816 std::string form2 = "%%*s";
817 std::string format = form2 + "%%d";
818 if (sscanf(line, format.c_str(), &ai) == 1)
823 format = form2 + "%%d";
824 n = sscanf(line, format.c_str(), &aj);
827 srenew(con->conect, ++con->nconect);
828 con->conect[con->nconect-1].ai = ai-1;
829 con->conect[con->nconect-1].aj = aj-1;
836 void gmx_conect_dump(FILE *fp, gmx_conect conect)
838 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
841 for (i = 0; (i < gc->nconect); i++)
843 fprintf(fp, "%6s%5d%5d\n", "CONECT",
844 gc->conect[i].ai+1, gc->conect[i].aj+1);
848 gmx_conect gmx_conect_init()
857 void gmx_conect_done(gmx_conect conect)
859 gmx_conect_t *gc = conect;
864 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
866 gmx_conect_t *gc = conect;
872 for (i = 0; (i < gc->nconect); i++)
874 if (((gc->conect[i].ai == ai) &&
875 (gc->conect[i].aj == aj)) ||
876 ((gc->conect[i].aj == ai) &&
877 (gc->conect[i].ai == aj)))
885 void gmx_conect_add(gmx_conect conect, int ai, int aj)
887 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
892 if (!gmx_conect_exist(conect, ai, aj))
894 srenew(gc->conect, ++gc->nconect);
895 gc->conect[gc->nconect-1].ai = ai;
896 gc->conect[gc->nconect-1].aj = aj;
900 int read_pdbfile(FILE *in, char *title, int *model_nr,
901 t_atoms *atoms, t_symtab *symtab, rvec x[], int *ePBC,
902 matrix box, gmx_bool bChange, gmx_conect conect)
904 gmx_conect_t *gc = conect;
906 gmx_bool bConnWarn = FALSE;
911 gmx_bool bStop = FALSE;
915 /* Only assume pbc when there is a CRYST1 entry */
923 atoms->haveMass = FALSE;
924 atoms->haveCharge = FALSE;
925 atoms->haveType = FALSE;
926 atoms->haveBState = FALSE;
927 atoms->havePdbInfo = (atoms->pdbinfo != nullptr);
933 while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
935 line_type = line2type(line);
941 natom = read_atom(symtab, line, line_type, natom, atoms, x, chainnum, bChange);
945 if (atoms->havePdbInfo)
947 read_anisou(line, natom, atoms);
952 read_cryst1(line, ePBC, box);
957 if (std::strlen(line) > 6)
960 /* skip HEADER or TITLE and spaces */
969 /* truncate after title */
970 d = std::strstr(c, " ");
975 if (std::strlen(c) > 0)
977 std::strcpy(title, c);
983 if ((!std::strstr(line, ": ")) || (std::strstr(line+6, "MOLECULE:")))
985 if (!(c = std::strstr(line+6, "MOLECULE:")) )
989 /* skip 'MOLECULE:' and spaces */
998 /* truncate after title */
1002 while ( (d[-1] == ';') && d > c)
1012 std::strcat(title, "; ");
1013 std::strcat(title, c);
1017 std::strcpy(title, c);
1031 sscanf(line, "%*s%d", model_nr);
1041 gmx_conect_addline(gc, line);
1043 else if (!bConnWarn)
1045 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
1058 void get_pdb_coordnum(FILE *in, int *natoms)
1063 while (fgets2(line, STRLEN, in))
1065 if (std::strncmp(line, "ENDMDL", 6) == 0)
1069 if ((std::strncmp(line, "ATOM ", 6) == 0) || (std::strncmp(line, "HETATM", 6) == 0))
1076 void gmx_pdb_read_conf(const char *infile,
1077 t_symtab *symtab, char **name, t_atoms *atoms,
1078 rvec x[], int *ePBC, matrix box)
1080 FILE *in = gmx_fio_fopen(infile, "r");
1082 read_pdbfile(in, title, nullptr, atoms, symtab, x, ePBC, box, TRUE, nullptr);
1083 if (name != nullptr)
1085 *name = gmx_strdup(title);
1090 gmx_conect gmx_conect_generate(const t_topology *top)
1095 /* Fill the conect records */
1096 gc = gmx_conect_init();
1098 for (f = 0; (f < F_NRE); f++)
1102 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms+1)
1104 gmx_conect_add(gc, top->idef.il[f].iatoms[i+1],
1105 top->idef.il[f].iatoms[i+2]);
1113 gmx_fprintf_pdb_atomline(FILE * fp,
1114 enum PDB_record record,
1115 int atom_seq_number,
1116 const char * atom_name,
1117 char alternate_location,
1118 const char * res_name,
1121 char res_insertion_code,
1127 const char * element)
1129 char tmp_atomname[6], tmp_resname[6];
1130 gmx_bool start_name_in_col13;
1133 if (record != epdbATOM && record != epdbHETATM)
1135 gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
1138 /* Format atom name */
1139 if (atom_name != nullptr)
1141 /* If the atom name is an element name with two chars, it should start already in column 13.
1142 * Otherwise it should start in column 14, unless the name length is 4 chars.
1144 if ( (element != nullptr) && (std::strlen(element) >= 2) && (gmx_strncasecmp(atom_name, element, 2) == 0) )
1146 start_name_in_col13 = TRUE;
1150 start_name_in_col13 = (std::strlen(atom_name) >= 4);
1152 snprintf(tmp_atomname, sizeof(tmp_atomname), start_name_in_col13 ? "" : " ");
1153 std::strncat(tmp_atomname, atom_name, 4);
1154 tmp_atomname[5] = '\0';
1158 tmp_atomname[0] = '\0';
1161 /* Format residue name */
1162 std::strncpy(tmp_resname, (res_name != nullptr) ? res_name : "", 4);
1163 /* Make sure the string is terminated if strlen was > 4 */
1164 tmp_resname[4] = '\0';
1165 /* String is properly terminated, so now we can use strcat. By adding a
1166 * space we can write it right-justified, and if the original name was
1167 * three characters or less there will be a space added on the right side.
1169 std::strcat(tmp_resname, " ");
1171 /* Truncate integers so they fit */
1172 atom_seq_number = atom_seq_number % 100000;
1173 res_seq_number = res_seq_number % 10000;
1176 "%-6s%5d %-4.4s%c%4.4s%c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f %2s\n",
1188 (element != nullptr) ? element : "");