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 static void xlate_atomname_gmx2pdb(char *name)
105 length = std::strlen(name);
106 if (length > 3 && std::isdigit(name[length-1]))
108 temp = name[length-1];
109 for (i = length-1; i > 0; --i)
118 void gmx_write_pdb_box(FILE *out, int ePBC, const matrix box)
120 real alpha, beta, gamma;
124 ePBC = guess_ePBC(box);
127 if (ePBC == epbcNONE)
132 if (norm2(box[YY])*norm2(box[ZZ]) != 0)
134 alpha = RAD2DEG*gmx_angle(box[YY], box[ZZ]);
140 if (norm2(box[XX])*norm2(box[ZZ]) != 0)
142 beta = RAD2DEG*gmx_angle(box[XX], box[ZZ]);
148 if (norm2(box[XX])*norm2(box[YY]) != 0)
150 gamma = RAD2DEG*gmx_angle(box[XX], box[YY]);
156 fprintf(out, "REMARK THIS IS A SIMULATION BOX\n");
157 if (ePBC != epbcSCREW)
159 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
160 10*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
161 alpha, beta, gamma, "P 1", 1);
165 /* Double the a-vector length and write the correct space group */
166 fprintf(out, "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
167 20*norm(box[XX]), 10*norm(box[YY]), 10*norm(box[ZZ]),
168 alpha, beta, gamma, "P 21 1 1", 1);
173 static void read_cryst1(char *line, int *ePBC, matrix box)
176 char sa[12], sb[12], sc[12], sg[SG_SIZE+1], ident;
177 double fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
178 int syma, symb, symc;
181 sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
184 if (strlen(line) >= 55)
186 strncpy(sg, line+55, SG_SIZE);
192 sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
193 if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
195 fc = strtod(sc, nullptr)*0.1;
196 ePBC_file = (fc > 0 ? epbcXYZ : epbcXY);
198 if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
200 ePBC_file = epbcSCREW;
210 fa = strtod(sa, nullptr)*0.1;
211 fb = strtod(sb, nullptr)*0.1;
212 fc = strtod(sc, nullptr)*0.1;
213 if (ePBC_file == epbcSCREW)
219 if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
223 cosa = std::cos(alpha*DEG2RAD);
231 cosb = std::cos(beta*DEG2RAD);
239 cosg = std::cos(gamma*DEG2RAD);
240 sing = std::sin(gamma*DEG2RAD);
247 box[YY][XX] = fb*cosg;
248 box[YY][YY] = fb*sing;
249 box[ZZ][XX] = fc*cosb;
250 box[ZZ][YY] = fc*(cosa - cosb*cosg)/sing;
251 box[ZZ][ZZ] = std::sqrt(fc*fc
252 - box[ZZ][XX]*box[ZZ][XX] - box[ZZ][YY]*box[ZZ][YY]);
263 gmx_fprintf_pqr_atomline(FILE * fp,
264 enum PDB_record record,
266 const char * atom_name,
267 const char * res_name,
276 GMX_RELEASE_ASSERT(record == epdbATOM || record == epdbHETATM,
277 "Can only print PQR atom lines as ATOM or HETATM records");
279 /* Check atom name */
280 GMX_RELEASE_ASSERT(atom_name != nullptr,
281 "Need atom information to print pqr");
283 /* Check residue name */
284 GMX_RELEASE_ASSERT(res_name != nullptr,
285 "Need residue information to print pqr");
287 /* Truncate integers so they fit */
288 atom_seq_number = atom_seq_number % 100000;
289 res_seq_number = res_seq_number % 10000;
292 "%s %d %s %s %c %d %8.3f %8.3f %8.3f %6.2f %6.2f\n",
306 void write_pdbfile_indexed(FILE *out, const char *title,
307 const t_atoms *atoms, const rvec x[],
308 int ePBC, const matrix box, char chainid,
309 int model_nr, int nindex, const int index[],
310 gmx_conect conect, gmx_bool bTerSepChains,
313 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
314 char resnm[6], nm[6];
317 enum PDB_record type;
318 unsigned char resic, ch;
322 int chainnum, lastchainnum;
323 gmx_residuetype_t*rt;
324 const char *p_restype;
325 const char *p_lastrestype;
327 gmx_residuetype_init(&rt);
329 fprintf(out, "TITLE %s\n", (title && title[0]) ? title : gmx::bromacs().c_str());
330 if (box && ( (norm2(box[XX]) != 0.0f) || (norm2(box[YY]) != 0.0f) || (norm2(box[ZZ]) != 0.0f) ) )
332 gmx_write_pdb_box(out, ePBC, box);
334 if (atoms->havePdbInfo)
336 /* Check whether any occupancies are set, in that case leave it as is,
337 * otherwise set them all to one
340 for (ii = 0; (ii < nindex) && bOccup; ii++)
343 bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
351 fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
356 for (ii = 0; ii < nindex; ii++)
359 resind = atoms->atom[i].resind;
360 chainnum = atoms->resinfo[resind].chainnum;
361 p_lastrestype = p_restype;
362 gmx_residuetype_get_type(rt, *atoms->resinfo[resind].name, &p_restype);
364 /* Add a TER record if we changed chain, and if either the previous or this chain is protein/DNA/RNA. */
365 if (bTerSepChains && ii > 0 && chainnum != lastchainnum)
367 /* Only add TER if the previous chain contained protein/DNA/RNA. */
368 if (gmx_residuetype_is_protein(rt, p_lastrestype) || gmx_residuetype_is_dna(rt, p_lastrestype) || gmx_residuetype_is_rna(rt, p_lastrestype))
370 fprintf(out, "TER\n");
372 lastchainnum = chainnum;
375 strncpy(resnm, *atoms->resinfo[resind].name, sizeof(resnm)-1);
376 resnm[sizeof(resnm)-1] = 0;
377 strncpy(nm, *atoms->atomname[i], sizeof(nm)-1);
378 nm[sizeof(nm)-1] = 0;
380 /* rename HG12 to 2HG1, etc. */
381 xlate_atomname_gmx2pdb(nm);
382 resnr = atoms->resinfo[resind].nr;
383 resic = atoms->resinfo[resind].ic;
390 ch = atoms->resinfo[resind].chainid;
399 resnr = resnr % 10000;
402 if (atoms->pdbinfo != nullptr)
404 pdbinfo = atoms->pdbinfo[i];
408 gmx_pdbinfo_init_default(&pdbinfo);
410 type = static_cast<enum PDB_record>(pdbinfo.type);
411 altloc = pdbinfo.altloc;
412 if (!isalnum(altloc))
416 occup = bOccup ? 1.0 : pdbinfo.occup;
420 gmx_fprintf_pdb_atomline(out,
429 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
432 atoms->atom[i].elem);
434 if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
436 fprintf(out, "ANISOU%5d %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n",
437 (i+1)%100000, nm, resnm, ch, resnr,
438 (resic == '\0') ? ' ' : resic,
439 atoms->pdbinfo[i].uij[0], atoms->pdbinfo[i].uij[1],
440 atoms->pdbinfo[i].uij[2], atoms->pdbinfo[i].uij[3],
441 atoms->pdbinfo[i].uij[4], atoms->pdbinfo[i].uij[5]);
446 gmx_fprintf_pqr_atomline(out,
453 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ],
459 fprintf(out, "TER\n");
460 fprintf(out, "ENDMDL\n");
464 /* Write conect records */
465 for (i = 0; (i < gc->nconect); i++)
467 fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai+1, gc->conect[i].aj+1);
471 gmx_residuetype_destroy(rt);
474 void write_pdbfile(FILE *out, const char *title, const t_atoms *atoms, const rvec x[],
475 int ePBC, const matrix box, char chainid, int model_nr, gmx_conect conect, gmx_bool bTerSepChains)
479 snew(index, atoms->nr);
480 for (i = 0; i < atoms->nr; i++)
484 write_pdbfile_indexed(out, title, atoms, x, ePBC, box, chainid, model_nr,
485 atoms->nr, index, conect, bTerSepChains, false);
489 static int line2type(const char *line)
494 for (k = 0; (k < 6); k++)
500 for (k = 0; (k < epdbNR); k++)
502 if (std::strncmp(type, pdbtp[k], strlen(pdbtp[k])) == 0)
511 static void read_anisou(char line[], int natom, t_atoms *atoms)
515 char anr[12], anm[12];
519 for (k = 0; (k < 5); k++, j++)
525 for (k = 0; (k < 4); k++, j++)
532 /* Strip off spaces */
535 /* Search backwards for number and name only */
536 atomnr = std::strtol(anr, nullptr, 10);
537 for (i = natom-1; (i >= 0); i--)
539 if ((std::strcmp(anm, *(atoms->atomname[i])) == 0) &&
540 (atomnr == atoms->pdbinfo[i].atomnr))
547 fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n",
552 if (sscanf(line+29, "%d%d%d%d%d%d",
553 &atoms->pdbinfo[i].uij[U11], &atoms->pdbinfo[i].uij[U22],
554 &atoms->pdbinfo[i].uij[U33], &atoms->pdbinfo[i].uij[U12],
555 &atoms->pdbinfo[i].uij[U13], &atoms->pdbinfo[i].uij[U23])
558 atoms->pdbinfo[i].bAnisotropic = TRUE;
562 fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
563 atoms->pdbinfo[i].bAnisotropic = FALSE;
568 void get_pdb_atomnumber(const t_atoms *atoms, gmx_atomprop_t aps)
570 int i, atomnumber, len;
572 char anm[6], anm_copy[6], *ptr;
578 gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
580 for (i = 0; (i < atoms->nr); i++)
582 std::strcpy(anm, atoms->pdbinfo[i].atomnm);
583 std::strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
584 bool atomNumberSet = false;
586 if ((anm[0] != ' ') && ((len <= 2) || !std::isdigit(anm[2])))
589 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
591 atomnumber = std::round(eval);
592 atomNumberSet = true;
597 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
599 atomnumber = std::round(eval);
600 atomNumberSet = true;
607 while ((k < std::strlen(anm)) && (std::isspace(anm[k]) || std::isdigit(anm[k])))
611 anm_copy[0] = anm[k];
613 if (gmx_atomprop_query(aps, epropElement, "???", anm_copy, &eval))
615 atomnumber = std::round(eval);
616 atomNumberSet = true;
621 atoms->atom[i].atomnumber = atomnumber;
622 ptr = gmx_atomprop_element(aps, atomnumber);
625 fprintf(debug, "Atomnumber for atom '%s' is %d\n",
633 std::strncpy(atoms->atom[i].elem, ptr == nullptr ? "" : ptr, 4);
637 static int read_atom(t_symtab *symtab,
638 const char line[], int type, int natom,
639 t_atoms *atoms, rvec x[], int chainnum, gmx_bool bChange)
644 char anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
645 char xc[12], yc[12], zc[12], occup[12], bfac[12];
648 int resnr, atomnumber;
650 if (natom >= atoms->nr)
652 gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)",
658 for (k = 0; (k < 5); k++, j++)
665 for (k = 0; (k < 4); k++, j++)
670 std::strcpy(anm_copy, anm);
676 for (k = 0; (k < 4); k++, j++)
686 for (k = 0; (k < 4); k++, j++)
692 resnr = std::strtol(rnr, nullptr, 10);
696 /* X,Y,Z Coordinate */
697 for (k = 0; (k < 8); k++, j++)
702 for (k = 0; (k < 8); k++, j++)
707 for (k = 0; (k < 8); k++, j++)
714 for (k = 0; (k < 6); k++, j++)
721 for (k = 0; (k < 7); k++, j++)
731 for (k = 0; (k < 2); k++, j++)
740 atomn = &(atoms->atom[natom]);
742 atoms->resinfo[atoms->atom[natom-1].resind].nr != resnr ||
743 atoms->resinfo[atoms->atom[natom-1].resind].ic != resic ||
744 (strcmp(*atoms->resinfo[atoms->atom[natom-1].resind].name, resnm) != 0))
752 atomn->resind = atoms->atom[natom-1].resind + 1;
754 atoms->nres = atomn->resind + 1;
755 t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
759 atomn->resind = atoms->atom[natom-1].resind;
763 xlate_atomname_pdb2gmx(anm);
765 atoms->atomname[natom] = put_symtab(symtab, anm);
768 atomn->atomnumber = atomnumber;
769 strncpy(atomn->elem, elem, 4);
771 x[natom][XX] = strtod(xc, nullptr)*0.1;
772 x[natom][YY] = strtod(yc, nullptr)*0.1;
773 x[natom][ZZ] = strtod(zc, nullptr)*0.1;
776 atoms->pdbinfo[natom].type = type;
777 atoms->pdbinfo[natom].atomnr = strtol(anr, nullptr, 10);
778 atoms->pdbinfo[natom].altloc = altloc;
779 strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
780 atoms->pdbinfo[natom].bfac = strtod(bfac, nullptr);
781 atoms->pdbinfo[natom].occup = strtod(occup, nullptr);
788 gmx_bool is_hydrogen(const char *nm)
792 std::strcpy(buf, nm);
799 else if ((std::isdigit(buf[0])) && (buf[1] == 'H'))
806 gmx_bool is_dummymass(const char *nm)
810 std::strcpy(buf, nm);
813 return (buf[0] == 'M') && (std::isdigit(buf[strlen(buf)-1]) != 0);
816 static void gmx_conect_addline(gmx_conect_t *con, char *line)
819 char format[32], form2[32];
821 sprintf(form2, "%%*s");
822 sprintf(format, "%s%%d", form2);
823 if (sscanf(line, format, &ai) == 1)
827 std::strcat(form2, "%*s");
828 sprintf(format, "%s%%d", form2);
829 n = sscanf(line, format, &aj);
832 srenew(con->conect, ++con->nconect);
833 con->conect[con->nconect-1].ai = ai-1;
834 con->conect[con->nconect-1].aj = aj-1;
841 void gmx_conect_dump(FILE *fp, gmx_conect conect)
843 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
846 for (i = 0; (i < gc->nconect); i++)
848 fprintf(fp, "%6s%5d%5d\n", "CONECT",
849 gc->conect[i].ai+1, gc->conect[i].aj+1);
853 gmx_conect gmx_conect_init()
862 void gmx_conect_done(gmx_conect conect)
864 gmx_conect_t *gc = conect;
869 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
871 gmx_conect_t *gc = conect;
877 for (i = 0; (i < gc->nconect); i++)
879 if (((gc->conect[i].ai == ai) &&
880 (gc->conect[i].aj == aj)) ||
881 ((gc->conect[i].aj == ai) &&
882 (gc->conect[i].ai == aj)))
890 void gmx_conect_add(gmx_conect conect, int ai, int aj)
892 gmx_conect_t *gc = static_cast<gmx_conect_t *>(conect);
897 if (!gmx_conect_exist(conect, ai, aj))
899 srenew(gc->conect, ++gc->nconect);
900 gc->conect[gc->nconect-1].ai = ai;
901 gc->conect[gc->nconect-1].aj = aj;
905 int read_pdbfile(FILE *in, char *title, int *model_nr,
906 t_atoms *atoms, t_symtab *symtab, rvec x[], int *ePBC,
907 matrix box, gmx_bool bChange, gmx_conect conect)
909 gmx_conect_t *gc = conect;
911 gmx_bool bConnWarn = FALSE;
916 gmx_bool bStop = FALSE;
920 /* Only assume pbc when there is a CRYST1 entry */
928 atoms->haveMass = FALSE;
929 atoms->haveCharge = FALSE;
930 atoms->haveType = FALSE;
931 atoms->haveBState = FALSE;
932 atoms->havePdbInfo = (atoms->pdbinfo != nullptr);
938 while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
940 line_type = line2type(line);
946 natom = read_atom(symtab, line, line_type, natom, atoms, x, chainnum, bChange);
950 if (atoms->havePdbInfo)
952 read_anisou(line, natom, atoms);
957 read_cryst1(line, ePBC, box);
962 if (std::strlen(line) > 6)
965 /* skip HEADER or TITLE and spaces */
974 /* truncate after title */
975 d = std::strstr(c, " ");
980 if (std::strlen(c) > 0)
982 std::strcpy(title, c);
988 if ((!std::strstr(line, ": ")) || (std::strstr(line+6, "MOLECULE:")))
990 if (!(c = std::strstr(line+6, "MOLECULE:")) )
994 /* skip 'MOLECULE:' and spaces */
1003 /* truncate after title */
1007 while ( (d[-1] == ';') && d > c)
1017 std::strcat(title, "; ");
1018 std::strcat(title, c);
1022 std::strcpy(title, c);
1036 sscanf(line, "%*s%d", model_nr);
1046 gmx_conect_addline(gc, line);
1048 else if (!bConnWarn)
1050 fprintf(stderr, "WARNING: all CONECT records are ignored\n");
1063 void get_pdb_coordnum(FILE *in, int *natoms)
1068 while (fgets2(line, STRLEN, in))
1070 if (std::strncmp(line, "ENDMDL", 6) == 0)
1074 if ((std::strncmp(line, "ATOM ", 6) == 0) || (std::strncmp(line, "HETATM", 6) == 0))
1081 void gmx_pdb_read_conf(const char *infile,
1082 t_symtab *symtab, char **name, t_atoms *atoms,
1083 rvec x[], int *ePBC, matrix box)
1085 FILE *in = gmx_fio_fopen(infile, "r");
1087 read_pdbfile(in, title, nullptr, atoms, symtab, x, ePBC, box, TRUE, nullptr);
1088 if (name != nullptr)
1090 *name = gmx_strdup(title);
1095 gmx_conect gmx_conect_generate(const t_topology *top)
1100 /* Fill the conect records */
1101 gc = gmx_conect_init();
1103 for (f = 0; (f < F_NRE); f++)
1107 for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms+1)
1109 gmx_conect_add(gc, top->idef.il[f].iatoms[i+1],
1110 top->idef.il[f].iatoms[i+2]);
1118 gmx_fprintf_pdb_atomline(FILE * fp,
1119 enum PDB_record record,
1120 int atom_seq_number,
1121 const char * atom_name,
1122 char alternate_location,
1123 const char * res_name,
1126 char res_insertion_code,
1132 const char * element)
1134 char tmp_atomname[6], tmp_resname[6];
1135 gmx_bool start_name_in_col13;
1138 if (record != epdbATOM && record != epdbHETATM)
1140 gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
1143 /* Format atom name */
1144 if (atom_name != nullptr)
1146 /* If the atom name is an element name with two chars, it should start already in column 13.
1147 * Otherwise it should start in column 14, unless the name length is 4 chars.
1149 if ( (element != nullptr) && (std::strlen(element) >= 2) && (gmx_strncasecmp(atom_name, element, 2) == 0) )
1151 start_name_in_col13 = TRUE;
1155 start_name_in_col13 = (std::strlen(atom_name) >= 4);
1157 snprintf(tmp_atomname, sizeof(tmp_atomname), start_name_in_col13 ? "" : " ");
1158 std::strncat(tmp_atomname, atom_name, 4);
1159 tmp_atomname[5] = '\0';
1163 tmp_atomname[0] = '\0';
1166 /* Format residue name */
1167 std::strncpy(tmp_resname, (res_name != nullptr) ? res_name : "", 4);
1168 /* Make sure the string is terminated if strlen was > 4 */
1169 tmp_resname[4] = '\0';
1170 /* String is properly terminated, so now we can use strcat. By adding a
1171 * space we can write it right-justified, and if the original name was
1172 * three characters or less there will be a space added on the right side.
1174 std::strcat(tmp_resname, " ");
1176 /* Truncate integers so they fit */
1177 atom_seq_number = atom_seq_number % 100000;
1178 res_seq_number = res_seq_number % 10000;
1181 "%-6s%5d %-4.4s%c%4.4s%c%4d%c %8.3f%8.3f%8.3f%6.2f%6.2f %2s\n",
1193 (element != nullptr) ? element : "");