07de14f92c63b6bab2e895be6d0cf2a44ee11d7c
[alexxy/gromacs.git] / src / gromacs / fileio / pdbio.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "gromacs/fileio/pdbio.h"
41
42 #include <algorithm>
43 #include <cctype>
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
47
48 #include <string>
49
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/topology/atomprop.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/topology/symtab.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/coolstuff.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/enumerationhelpers.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/snprintf.h"
66 #include "gromacs/utility/stringtoenumvalueconverter.h"
67
68 typedef struct
69 {
70     int ai, aj;
71 } gmx_conection_t;
72
73 typedef struct gmx_conect_t
74 {
75     int              nconect;
76     gmx_bool         bSorted;
77     gmx_conection_t* conect;
78 } gmx_conect_t;
79
80 const char* enumValueToString(PdbRecordType enumValue)
81 {
82     static constexpr gmx::EnumerationArray<PdbRecordType, const char*> pdbRecordTypeName = {
83         "ATOM  ", "HETATM", "ANISOU", "CRYST1", "COMPND", "MODEL",
84         "ENDMDL", "TER",    "HEADER", "TITLE",  "REMARK", "CONECT"
85     };
86     return pdbRecordTypeName[enumValue];
87 }
88
89 void gmx_write_pdb_box(FILE* out, PbcType pbcType, const matrix box)
90 {
91     real alpha, beta, gamma;
92
93     if (pbcType == PbcType::Unset)
94     {
95         pbcType = guessPbcType(box);
96     }
97
98     if (pbcType == PbcType::No)
99     {
100         return;
101     }
102
103     if (norm2(box[YY]) * norm2(box[ZZ]) != 0)
104     {
105         alpha = gmx::c_rad2Deg * gmx_angle(box[YY], box[ZZ]);
106     }
107     else
108     {
109         alpha = 90;
110     }
111     if (norm2(box[XX]) * norm2(box[ZZ]) != 0)
112     {
113         beta = gmx::c_rad2Deg * gmx_angle(box[XX], box[ZZ]);
114     }
115     else
116     {
117         beta = 90;
118     }
119     if (norm2(box[XX]) * norm2(box[YY]) != 0)
120     {
121         gamma = gmx::c_rad2Deg * gmx_angle(box[XX], box[YY]);
122     }
123     else
124     {
125         gamma = 90;
126     }
127     fprintf(out, "REMARK    THIS IS A SIMULATION BOX\n");
128     if (pbcType != PbcType::Screw)
129     {
130         fprintf(out,
131                 "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
132                 10 * norm(box[XX]),
133                 10 * norm(box[YY]),
134                 10 * norm(box[ZZ]),
135                 alpha,
136                 beta,
137                 gamma,
138                 "P 1",
139                 1);
140     }
141     else
142     {
143         /* Double the a-vector length and write the correct space group */
144         fprintf(out,
145                 "CRYST1%9.3f%9.3f%9.3f%7.2f%7.2f%7.2f %-11s%4d\n",
146                 20 * norm(box[XX]),
147                 10 * norm(box[YY]),
148                 10 * norm(box[ZZ]),
149                 alpha,
150                 beta,
151                 gamma,
152                 "P 21 1 1",
153                 1);
154     }
155 }
156
157 static void read_cryst1(char* line, PbcType* pbcType, matrix box)
158 {
159 #define SG_SIZE 11
160     char    sa[12], sb[12], sc[12], sg[SG_SIZE + 1], ident;
161     double  fa, fb, fc, alpha, beta, gamma, cosa, cosb, cosg, sing;
162     int     syma, symb, symc;
163     PbcType pbcTypeFile;
164
165     sscanf(line, "%*s%s%s%s%lf%lf%lf", sa, sb, sc, &alpha, &beta, &gamma);
166
167     pbcTypeFile = PbcType::Unset;
168     if (strlen(line) >= 55)
169     {
170         strncpy(sg, line + 55, SG_SIZE);
171         sg[SG_SIZE] = '\0';
172         ident       = ' ';
173         syma        = 0;
174         symb        = 0;
175         symc        = 0;
176         sscanf(sg, "%c %d %d %d", &ident, &syma, &symb, &symc);
177         if (ident == 'P' && syma == 1 && symb <= 1 && symc <= 1)
178         {
179             fc          = strtod(sc, nullptr) * 0.1;
180             pbcTypeFile = (fc > 0 ? PbcType::Xyz : PbcType::XY);
181         }
182         if (ident == 'P' && syma == 21 && symb == 1 && symc == 1)
183         {
184             pbcTypeFile = PbcType::Screw;
185         }
186     }
187     if (pbcType)
188     {
189         *pbcType = pbcTypeFile;
190     }
191
192     if (box)
193     {
194         fa = strtod(sa, nullptr) * 0.1;
195         fb = strtod(sb, nullptr) * 0.1;
196         fc = strtod(sc, nullptr) * 0.1;
197         if (pbcTypeFile == PbcType::Screw)
198         {
199             fa *= 0.5;
200         }
201         clear_mat(box);
202         box[XX][XX] = fa;
203         if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0))
204         {
205             if (alpha != 90.0)
206             {
207                 cosa = std::cos(alpha * gmx::c_deg2Rad);
208             }
209             else
210             {
211                 cosa = 0;
212             }
213             if (beta != 90.0)
214             {
215                 cosb = std::cos(beta * gmx::c_deg2Rad);
216             }
217             else
218             {
219                 cosb = 0;
220             }
221             if (gamma != 90.0)
222             {
223                 cosg = std::cos(gamma * gmx::c_deg2Rad);
224                 sing = std::sin(gamma * gmx::c_deg2Rad);
225             }
226             else
227             {
228                 cosg = 0;
229                 sing = 1;
230             }
231             box[YY][XX] = fb * cosg;
232             box[YY][YY] = fb * sing;
233             box[ZZ][XX] = fc * cosb;
234             box[ZZ][YY] = fc * (cosa - cosb * cosg) / sing;
235             box[ZZ][ZZ] = std::sqrt(fc * fc - box[ZZ][XX] * box[ZZ][XX] - box[ZZ][YY] * box[ZZ][YY]);
236         }
237         else
238         {
239             box[YY][YY] = fb;
240             box[ZZ][ZZ] = fc;
241         }
242     }
243 }
244
245 static int gmx_fprintf_pqr_atomline(FILE*         fp,
246                                     PdbRecordType record,
247                                     int           atom_seq_number,
248                                     const char*   atom_name,
249                                     const char*   res_name,
250                                     char          chain_id,
251                                     int           res_seq_number,
252                                     real          x,
253                                     real          y,
254                                     real          z,
255                                     real          occupancy,
256                                     real          b_factor)
257 {
258     GMX_RELEASE_ASSERT(record == PdbRecordType::Atom || record == PdbRecordType::Hetatm,
259                        "Can only print PQR atom lines as ATOM or HETATM records");
260
261     /* Check atom name */
262     GMX_RELEASE_ASSERT(atom_name != nullptr, "Need atom information to print pqr");
263
264     /* Check residue name */
265     GMX_RELEASE_ASSERT(res_name != nullptr, "Need residue information to print pqr");
266
267     /* Truncate integers so they fit */
268     atom_seq_number = atom_seq_number % 100000;
269     res_seq_number  = res_seq_number % 10000;
270
271     int n = fprintf(fp,
272                     "%-6s%5d %-4.4s%4.4s%c%4d %8.3f %8.3f %8.3f %6.2f %6.2f\n",
273                     enumValueToString(record),
274                     atom_seq_number,
275                     atom_name,
276                     res_name,
277                     chain_id,
278                     res_seq_number,
279                     x,
280                     y,
281                     z,
282                     occupancy,
283                     b_factor);
284
285     return n;
286 }
287
288 void write_pdbfile_indexed(FILE*          out,
289                            const char*    title,
290                            const t_atoms* atoms,
291                            const rvec     x[],
292                            PbcType        pbcType,
293                            const matrix   box,
294                            char           chainid,
295                            int            model_nr,
296                            int            nindex,
297                            const int      index[],
298                            gmx_conect     conect,
299                            bool           usePqrFormat)
300 {
301     gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
302     PdbRecordType type;
303     char          altloc;
304     real          occup, bfac;
305     gmx_bool      bOccup;
306
307
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)))
310     {
311         gmx_write_pdb_box(out, pbcType, box);
312     }
313     if (atoms->havePdbInfo)
314     {
315         /* Check whether any occupancies are set, in that case leave it as is,
316          * otherwise set them all to one
317          */
318         bOccup = TRUE;
319         for (int ii = 0; (ii < nindex) && bOccup; ii++)
320         {
321             int i  = index[ii];
322             bOccup = bOccup && (atoms->pdbinfo[i].occup == 0.0);
323         }
324     }
325     else
326     {
327         bOccup = FALSE;
328     }
329
330     fprintf(out, "MODEL %8d\n", model_nr > 0 ? model_nr : 1);
331
332     for (int ii = 0; ii < nindex; ii++)
333     {
334         int         i      = index[ii];
335         int         resind = atoms->atom[i].resind;
336         std::string resnm  = *atoms->resinfo[resind].name;
337         std::string nm     = *atoms->atomname[i];
338
339         int           resnr = atoms->resinfo[resind].nr;
340         unsigned char resic = atoms->resinfo[resind].ic;
341         unsigned char ch;
342         if (chainid != ' ')
343         {
344             ch = chainid;
345         }
346         else
347         {
348             ch = atoms->resinfo[resind].chainid;
349
350             if (ch == 0)
351             {
352                 ch = ' ';
353             }
354         }
355         if (resnr >= 10000)
356         {
357             resnr = resnr % 10000;
358         }
359         t_pdbinfo pdbinfo;
360         if (atoms->pdbinfo != nullptr)
361         {
362             pdbinfo = atoms->pdbinfo[i];
363         }
364         else
365         {
366             gmx_pdbinfo_init_default(&pdbinfo);
367         }
368         type   = pdbinfo.type;
369         altloc = pdbinfo.altloc;
370         if (!isalnum(altloc))
371         {
372             altloc = ' ';
373         }
374         occup = bOccup ? 1.0 : pdbinfo.occup;
375         bfac  = pdbinfo.bfac;
376         if (!usePqrFormat)
377         {
378             gmx_fprintf_pdb_atomline(out,
379                                      type,
380                                      i + 1,
381                                      nm.c_str(),
382                                      altloc,
383                                      resnm.c_str(),
384                                      ch,
385                                      resnr,
386                                      resic,
387                                      10 * x[i][XX],
388                                      10 * x[i][YY],
389                                      10 * x[i][ZZ],
390                                      occup,
391                                      bfac,
392                                      atoms->atom[i].elem);
393
394             if (atoms->pdbinfo && atoms->pdbinfo[i].bAnisotropic)
395             {
396                 fprintf(out,
397                         "ANISOU%5d  %-4.4s%4.4s%c%4d%c %7d%7d%7d%7d%7d%7d\n",
398                         (i + 1) % 100000,
399                         nm.c_str(),
400                         resnm.c_str(),
401                         ch,
402                         resnr,
403                         (resic == '\0') ? ' ' : resic,
404                         atoms->pdbinfo[i].uij[0],
405                         atoms->pdbinfo[i].uij[1],
406                         atoms->pdbinfo[i].uij[2],
407                         atoms->pdbinfo[i].uij[3],
408                         atoms->pdbinfo[i].uij[4],
409                         atoms->pdbinfo[i].uij[5]);
410             }
411         }
412         else
413         {
414             gmx_fprintf_pqr_atomline(out,
415                                      type,
416                                      i + 1,
417                                      nm.c_str(),
418                                      resnm.c_str(),
419                                      ch,
420                                      resnr,
421                                      10 * x[i][XX],
422                                      10 * x[i][YY],
423                                      10 * x[i][ZZ],
424                                      occup,
425                                      bfac);
426         }
427     }
428
429     fprintf(out, "TER\n");
430     fprintf(out, "ENDMDL\n");
431
432     if (nullptr != gc)
433     {
434         /* Write conect records */
435         for (int i = 0; (i < gc->nconect); i++)
436         {
437             fprintf(out, "CONECT%5d%5d\n", gc->conect[i].ai + 1, gc->conect[i].aj + 1);
438         }
439     }
440 }
441
442 void write_pdbfile(FILE*          out,
443                    const char*    title,
444                    const t_atoms* atoms,
445                    const rvec     x[],
446                    PbcType        pbcType,
447                    const matrix   box,
448                    char           chainid,
449                    int            model_nr,
450                    gmx_conect     conect)
451 {
452     int i, *index;
453
454     snew(index, atoms->nr);
455     for (i = 0; i < atoms->nr; i++)
456     {
457         index[i] = i;
458     }
459     write_pdbfile_indexed(
460             out, title, atoms, x, pbcType, box, chainid, model_nr, atoms->nr, index, conect, false);
461     sfree(index);
462 }
463
464 static void read_anisou(char line[], int natom, t_atoms* atoms)
465 {
466     int  i, j, k, atomnr;
467     char nc = '\0';
468     char anr[12], anm[12];
469
470     /* Skip over type */
471     j = 6;
472     for (k = 0; (k < 5); k++, j++)
473     {
474         anr[k] = line[j];
475     }
476     anr[k] = nc;
477     j++;
478     for (k = 0; (k < 4); k++, j++)
479     {
480         anm[k] = line[j];
481     }
482     anm[k] = nc;
483     j++;
484
485     /* Strip off spaces */
486     trim(anm);
487
488     /* Search backwards for number and name only */
489     atomnr = std::strtol(anr, nullptr, 10);
490     for (i = natom - 1; (i >= 0); i--)
491     {
492         if ((std::strcmp(anm, *(atoms->atomname[i])) == 0) && (atomnr == atoms->pdbinfo[i].atomnr))
493         {
494             break;
495         }
496     }
497     if (i < 0)
498     {
499         fprintf(stderr, "Skipping ANISOU record (atom %s %d not found)\n", anm, atomnr);
500     }
501     else
502     {
503         if (sscanf(line + 29,
504                    "%d%d%d%d%d%d",
505                    &atoms->pdbinfo[i].uij[U11],
506                    &atoms->pdbinfo[i].uij[U22],
507                    &atoms->pdbinfo[i].uij[U33],
508                    &atoms->pdbinfo[i].uij[U12],
509                    &atoms->pdbinfo[i].uij[U13],
510                    &atoms->pdbinfo[i].uij[U23])
511             == 6)
512         {
513             atoms->pdbinfo[i].bAnisotropic = TRUE;
514         }
515         else
516         {
517             fprintf(stderr, "Invalid ANISOU record for atom %d\n", i);
518             atoms->pdbinfo[i].bAnisotropic = FALSE;
519         }
520     }
521 }
522
523 void get_pdb_atomnumber(const t_atoms* atoms, AtomProperties* aps)
524 {
525     int    i, atomnumber, len;
526     size_t k;
527     char   anm[6], anm_copy[6];
528     char   nc = '\0';
529     real   eval;
530
531     if (!atoms->pdbinfo)
532     {
533         gmx_incons("Trying to deduce atomnumbers when no pdb information is present");
534     }
535     for (i = 0; (i < atoms->nr); i++)
536     {
537         std::strcpy(anm, atoms->pdbinfo[i].atomnm);
538         std::strcpy(anm_copy, atoms->pdbinfo[i].atomnm);
539         bool atomNumberSet = false;
540         len                = strlen(anm);
541         if ((anm[0] != ' ') && ((len <= 2) || !std::isdigit(anm[2])))
542         {
543             anm_copy[2] = nc;
544             if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
545             {
546                 atomnumber    = gmx::roundToInt(eval);
547                 atomNumberSet = true;
548             }
549             else
550             {
551                 anm_copy[1] = nc;
552                 if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
553                 {
554                     atomnumber    = gmx::roundToInt(eval);
555                     atomNumberSet = true;
556                 }
557             }
558         }
559         if (!atomNumberSet)
560         {
561             k = 0;
562             while ((k < std::strlen(anm)) && (std::isspace(anm[k]) || std::isdigit(anm[k])))
563             {
564                 k++;
565             }
566             anm_copy[0] = anm[k];
567             anm_copy[1] = nc;
568             if (aps->setAtomProperty(epropElement, "???", anm_copy, &eval))
569             {
570                 atomnumber    = gmx::roundToInt(eval);
571                 atomNumberSet = true;
572             }
573         }
574         static constexpr size_t sc_maxElementNameLength = 3;
575         static_assert(sizeof(atoms->atom[i].elem) >= sc_maxElementNameLength + 1);
576         std::string element;
577         if (atomNumberSet)
578         {
579             atoms->atom[i].atomnumber = atomnumber;
580             element                   = aps->elementFromAtomNumber(atomnumber);
581             if (debug)
582             {
583                 fprintf(debug, "Atomnumber for atom '%s' is %d\n", anm, atomnumber);
584             }
585         }
586         element.resize(sc_maxElementNameLength);
587         std::strcpy(atoms->atom[i].elem, element.c_str());
588     }
589 }
590
591 static int
592 read_atom(t_symtab* symtab, const char line[], PdbRecordType type, int natom, t_atoms* atoms, rvec x[], int chainnum)
593 {
594     t_atom*       atomn;
595     int           j, k;
596     char          nc = '\0';
597     char          anr[12], anm[12], anm_copy[12], altloc, resnm[12], rnr[12], elem[3];
598     char          xc[12], yc[12], zc[12], occup[12], bfac[12];
599     unsigned char resic;
600     char          chainid;
601     int           resnr, atomnumber;
602
603     if (natom >= atoms->nr)
604     {
605         gmx_fatal(FARGS, "\nFound more atoms (%d) in pdb file than expected (%d)", natom + 1, atoms->nr);
606     }
607
608     /* Skip over type */
609     j = 6;
610     for (k = 0; (k < 5); k++, j++)
611     {
612         anr[k] = line[j];
613     }
614     anr[k] = nc;
615     trim(anr);
616     j++;
617     for (k = 0; (k < 4); k++, j++)
618     {
619         anm[k] = line[j];
620     }
621     anm[k] = nc;
622     std::strcpy(anm_copy, anm);
623     rtrim(anm_copy);
624     atomnumber = 0;
625     trim(anm);
626     altloc = line[j];
627     j++;
628     for (k = 0; (k < 4); k++, j++)
629     {
630         resnm[k] = line[j];
631     }
632     resnm[k] = nc;
633     trim(resnm);
634
635     chainid = line[j];
636     j++;
637
638     for (k = 0; (k < 4); k++, j++)
639     {
640         rnr[k] = line[j];
641     }
642     rnr[k] = nc;
643     trim(rnr);
644     resnr = std::strtol(rnr, nullptr, 10);
645     resic = line[j];
646     j += 4;
647
648     /* X,Y,Z Coordinate */
649     for (k = 0; (k < 8); k++, j++)
650     {
651         xc[k] = line[j];
652     }
653     xc[k] = nc;
654     for (k = 0; (k < 8); k++, j++)
655     {
656         yc[k] = line[j];
657     }
658     yc[k] = nc;
659     for (k = 0; (k < 8); k++, j++)
660     {
661         zc[k] = line[j];
662     }
663     zc[k] = nc;
664
665     /* Weight */
666     for (k = 0; (k < 6); k++, j++)
667     {
668         occup[k] = line[j];
669     }
670     occup[k] = nc;
671
672     /* B-Factor */
673     for (k = 0; (k < 7); k++, j++)
674     {
675         bfac[k] = line[j];
676     }
677     bfac[k] = nc;
678
679     /* 10 blanks */
680     j += 10;
681
682     /* Element name */
683     for (k = 0; (k < 2); k++, j++)
684     {
685         elem[k] = line[j];
686     }
687     elem[k] = nc;
688     trim(elem);
689
690     if (atoms->atom)
691     {
692         atomn = &(atoms->atom[natom]);
693         if ((natom == 0) || atoms->resinfo[atoms->atom[natom - 1].resind].nr != resnr
694             || atoms->resinfo[atoms->atom[natom - 1].resind].ic != resic
695             || (strcmp(*atoms->resinfo[atoms->atom[natom - 1].resind].name, resnm) != 0))
696         {
697             if (natom == 0)
698             {
699                 atomn->resind = 0;
700             }
701             else
702             {
703                 atomn->resind = atoms->atom[natom - 1].resind + 1;
704             }
705             atoms->nres = atomn->resind + 1;
706             t_atoms_set_resinfo(atoms, natom, symtab, resnm, resnr, resic, chainnum, chainid);
707         }
708         else
709         {
710             atomn->resind = atoms->atom[natom - 1].resind;
711         }
712         atoms->atomname[natom] = put_symtab(symtab, anm);
713         atomn->m               = 0.0;
714         atomn->q               = 0.0;
715         atomn->atomnumber      = atomnumber;
716         strncpy(atomn->elem, elem, 4);
717     }
718     x[natom][XX] = strtod(xc, nullptr) * 0.1;
719     x[natom][YY] = strtod(yc, nullptr) * 0.1;
720     x[natom][ZZ] = strtod(zc, nullptr) * 0.1;
721     if (atoms->pdbinfo)
722     {
723         atoms->pdbinfo[natom].type   = type;
724         atoms->pdbinfo[natom].atomnr = strtol(anr, nullptr, 10);
725         atoms->pdbinfo[natom].altloc = altloc;
726         strcpy(atoms->pdbinfo[natom].atomnm, anm_copy);
727         atoms->pdbinfo[natom].bfac  = strtod(bfac, nullptr);
728         atoms->pdbinfo[natom].occup = strtod(occup, nullptr);
729     }
730     natom++;
731
732     return natom;
733 }
734
735 gmx_bool is_hydrogen(const char* nm)
736 {
737     char buf[30];
738
739     std::strcpy(buf, nm);
740     trim(buf);
741
742     return ((buf[0] == 'H') || ((std::isdigit(buf[0]) != 0) && (buf[1] == 'H')));
743 }
744
745 gmx_bool is_dummymass(const char* nm)
746 {
747     char buf[30];
748
749     std::strcpy(buf, nm);
750     trim(buf);
751
752     return (buf[0] == 'M') && (std::isdigit(buf[strlen(buf) - 1]) != 0);
753 }
754
755 static void gmx_conect_addline(gmx_conect_t* con, char* line)
756 {
757     int n, ai, aj;
758
759     std::string form2  = "%*s";
760     std::string format = form2 + "%d";
761     if (sscanf(line, format.c_str(), &ai) == 1)
762     {
763         do
764         {
765             form2 += "%*s";
766             format = form2 + "%d";
767             n      = sscanf(line, format.c_str(), &aj);
768             if (n == 1)
769             {
770                 gmx_conect_add(con, ai - 1, aj - 1); /* to prevent duplicated records */
771             }
772         } while (n == 1);
773     }
774 }
775
776 void gmx_conect_dump(FILE* fp, gmx_conect conect)
777 {
778     gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
779     int           i;
780
781     for (i = 0; (i < gc->nconect); i++)
782     {
783         fprintf(fp, "%6s%5d%5d\n", "CONECT", gc->conect[i].ai + 1, gc->conect[i].aj + 1);
784     }
785 }
786
787 gmx_conect gmx_conect_init()
788 {
789     gmx_conect_t* gc;
790
791     snew(gc, 1);
792
793     return gc;
794 }
795
796 void gmx_conect_done(gmx_conect conect)
797 {
798     gmx_conect_t* gc = conect;
799
800     sfree(gc->conect);
801 }
802
803 gmx_bool gmx_conect_exist(gmx_conect conect, int ai, int aj)
804 {
805     gmx_conect_t* gc = conect;
806     int           i;
807
808     /* if (!gc->bSorted)
809        sort_conect(gc);*/
810
811     for (i = 0; (i < gc->nconect); i++)
812     {
813         if (((gc->conect[i].ai == ai) && (gc->conect[i].aj == aj))
814             || ((gc->conect[i].aj == ai) && (gc->conect[i].ai == aj)))
815         {
816             return TRUE;
817         }
818     }
819     return FALSE;
820 }
821
822 void gmx_conect_add(gmx_conect conect, int ai, int aj)
823 {
824     gmx_conect_t* gc = static_cast<gmx_conect_t*>(conect);
825
826     /* if (!gc->bSorted)
827        sort_conect(gc);*/
828
829     if (!gmx_conect_exist(conect, ai, aj))
830     {
831         srenew(gc->conect, ++gc->nconect);
832         gc->conect[gc->nconect - 1].ai = ai;
833         gc->conect[gc->nconect - 1].aj = aj;
834     }
835 }
836
837 int read_pdbfile(FILE*      in,
838                  char*      title,
839                  int*       model_nr,
840                  t_atoms*   atoms,
841                  t_symtab*  symtab,
842                  rvec       x[],
843                  PbcType*   pbcType,
844                  matrix     box,
845                  gmx_conect conect)
846 {
847     gmx_conect_t* gc = conect;
848     gmx_bool      bCOMPND;
849     gmx_bool      bConnWarn = FALSE;
850     char          line[STRLEN + 1];
851     char *        c, *d;
852     int           natom, chainnum;
853     gmx_bool      bStop = FALSE;
854
855     if (pbcType)
856     {
857         /* Only assume pbc when there is a CRYST1 entry */
858         *pbcType = PbcType::No;
859     }
860     if (box != nullptr)
861     {
862         clear_mat(box);
863     }
864
865     atoms->haveMass    = FALSE;
866     atoms->haveCharge  = FALSE;
867     atoms->haveType    = FALSE;
868     atoms->haveBState  = FALSE;
869     atoms->havePdbInfo = (atoms->pdbinfo != nullptr);
870
871     bCOMPND  = FALSE;
872     title[0] = '\0';
873     natom    = 0;
874     chainnum = 0;
875     static const gmx::StringToEnumValueConverter<PdbRecordType, enumValueToString, gmx::StringCompareType::Exact, gmx::StripStrings::Yes>
876             sc_pdbRecordTypeIdentifier;
877     while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
878     {
879         // Convert line to a null-terminated string, then take a substring of at most 6 chars
880         std::string                  lineAsString(line);
881         std::optional<PdbRecordType> line_type =
882                 sc_pdbRecordTypeIdentifier.valueFrom(lineAsString.substr(0, 6));
883         if (!line_type)
884         {
885             // Skip this line because it does not start with a
886             // recognized leading string describing a PDB record type.
887             continue;
888         }
889
890         switch (line_type.value())
891         {
892             case PdbRecordType::Atom:
893             case PdbRecordType::Hetatm:
894                 natom = read_atom(symtab, line, line_type.value(), natom, atoms, x, chainnum);
895                 break;
896
897             case PdbRecordType::Anisou:
898                 if (atoms->havePdbInfo)
899                 {
900                     read_anisou(line, natom, atoms);
901                 }
902                 break;
903
904             case PdbRecordType::Cryst1: read_cryst1(line, pbcType, box); break;
905
906             case PdbRecordType::Title:
907             case PdbRecordType::Header:
908                 if (std::strlen(line) > 6)
909                 {
910                     c = line + 6;
911                     /* skip HEADER or TITLE and spaces */
912                     while (c[0] != ' ')
913                     {
914                         c++;
915                     }
916                     while (c[0] == ' ')
917                     {
918                         c++;
919                     }
920                     /* truncate after title */
921                     d = std::strstr(c, "      ");
922                     if (d)
923                     {
924                         d[0] = '\0';
925                     }
926                     if (std::strlen(c) > 0)
927                     {
928                         std::strcpy(title, c);
929                     }
930                 }
931                 break;
932
933             case PdbRecordType::Compound:
934                 if ((!std::strstr(line, ": ")) || (std::strstr(line + 6, "MOLECULE:")))
935                 {
936                     if (!(c = std::strstr(line + 6, "MOLECULE:")))
937                     {
938                         c = line;
939                     }
940                     /* skip 'MOLECULE:' and spaces */
941                     while (c[0] != ' ')
942                     {
943                         c++;
944                     }
945                     while (c[0] == ' ')
946                     {
947                         c++;
948                     }
949                     /* truncate after title */
950                     d = strstr(c, "   ");
951                     if (d)
952                     {
953                         while ((d[-1] == ';') && d > c)
954                         {
955                             d--;
956                         }
957                         d[0] = '\0';
958                     }
959                     if (strlen(c) > 0)
960                     {
961                         if (bCOMPND)
962                         {
963                             std::strcat(title, "; ");
964                             std::strcat(title, c);
965                         }
966                         else
967                         {
968                             std::strcpy(title, c);
969                         }
970                     }
971                     bCOMPND = TRUE;
972                 }
973                 break;
974
975             case PdbRecordType::Ter: chainnum++; break;
976
977             case PdbRecordType::Model:
978                 if (model_nr)
979                 {
980                     sscanf(line, "%*s%d", model_nr);
981                 }
982                 break;
983
984             case PdbRecordType::EndModel: bStop = TRUE; break;
985             case PdbRecordType::Conect:
986                 if (gc)
987                 {
988                     gmx_conect_addline(gc, line);
989                 }
990                 else if (!bConnWarn)
991                 {
992                     fprintf(stderr, "WARNING: all CONECT records are ignored\n");
993                     bConnWarn = TRUE;
994                 }
995                 break;
996
997             default: break;
998         }
999     }
1000
1001     return natom;
1002 }
1003
1004 void get_pdb_coordnum(FILE* in, int* natoms)
1005 {
1006     char line[STRLEN];
1007
1008     *natoms = 0;
1009     while (fgets2(line, STRLEN, in))
1010     {
1011         if (std::strncmp(line, "ENDMDL", 6) == 0)
1012         {
1013             break;
1014         }
1015         if ((std::strncmp(line, "ATOM  ", 6) == 0) || (std::strncmp(line, "HETATM", 6) == 0))
1016         {
1017             (*natoms)++;
1018         }
1019     }
1020 }
1021
1022 void gmx_pdb_read_conf(const char* infile, t_symtab* symtab, char** name, t_atoms* atoms, rvec x[], PbcType* pbcType, matrix box)
1023 {
1024     FILE* in = gmx_fio_fopen(infile, "r");
1025     char  title[STRLEN];
1026     read_pdbfile(in, title, nullptr, atoms, symtab, x, pbcType, box, nullptr);
1027     if (name != nullptr)
1028     {
1029         *name = gmx_strdup(title);
1030     }
1031     gmx_fio_fclose(in);
1032 }
1033
1034 gmx_conect gmx_conect_generate(const t_topology* top)
1035 {
1036     int        f, i;
1037     gmx_conect gc;
1038
1039     /* Fill the conect records */
1040     gc = gmx_conect_init();
1041
1042     for (f = 0; (f < F_NRE); f++)
1043     {
1044         if (IS_CHEMBOND(f))
1045         {
1046             for (i = 0; (i < top->idef.il[f].nr); i += interaction_function[f].nratoms + 1)
1047             {
1048                 gmx_conect_add(gc, top->idef.il[f].iatoms[i + 1], top->idef.il[f].iatoms[i + 2]);
1049             }
1050         }
1051     }
1052     return gc;
1053 }
1054
1055 int gmx_fprintf_pdb_atomline(FILE*         fp,
1056                              PdbRecordType record,
1057                              int           atom_seq_number,
1058                              const char*   atom_name,
1059                              char          alternate_location,
1060                              const char*   res_name,
1061                              char          chain_id,
1062                              int           res_seq_number,
1063                              char          res_insertion_code,
1064                              real          x,
1065                              real          y,
1066                              real          z,
1067                              real          occupancy,
1068                              real          b_factor,
1069                              const char*   element)
1070 {
1071     char     tmp_atomname[6], tmp_resname[6];
1072     gmx_bool start_name_in_col13;
1073     int      n;
1074
1075     if (record != PdbRecordType::Atom && record != PdbRecordType::Hetatm)
1076     {
1077         gmx_fatal(FARGS, "Can only print PDB atom lines as ATOM or HETATM records");
1078     }
1079
1080     /* Format atom name */
1081     if (atom_name != nullptr)
1082     {
1083         /* If the atom name is an element name with two chars, it should start already in column 13.
1084          * Otherwise it should start in column 14, unless the name length is 4 chars.
1085          */
1086         if ((element != nullptr) && (std::strlen(element) >= 2)
1087             && (gmx_strncasecmp(atom_name, element, 2) == 0))
1088         {
1089             start_name_in_col13 = TRUE;
1090         }
1091         else
1092         {
1093             start_name_in_col13 = (std::strlen(atom_name) >= 4);
1094         }
1095         snprintf(tmp_atomname, sizeof(tmp_atomname), start_name_in_col13 ? "" : " ");
1096         std::strncat(tmp_atomname, atom_name, 4);
1097         tmp_atomname[5] = '\0';
1098     }
1099     else
1100     {
1101         tmp_atomname[0] = '\0';
1102     }
1103
1104     /* Format residue name */
1105     std::strncpy(tmp_resname, (res_name != nullptr) ? res_name : "", 4);
1106     /* Make sure the string is terminated if strlen was > 4 */
1107     tmp_resname[4] = '\0';
1108     /* String is properly terminated, so now we can use strcat. By adding a
1109      * space we can write it right-justified, and if the original name was
1110      * three characters or less there will be a space added on the right side.
1111      */
1112     std::strcat(tmp_resname, " ");
1113
1114     /* Truncate integers so they fit */
1115     atom_seq_number = atom_seq_number % 100000;
1116     res_seq_number  = res_seq_number % 10000;
1117
1118     n = fprintf(fp,
1119                 "%-6s%5d %-4.4s%c%4.4s%c%4d%c   %8.3f%8.3f%8.3f%6.2f%6.2f          %2s\n",
1120                 enumValueToString(record),
1121                 atom_seq_number,
1122                 tmp_atomname,
1123                 alternate_location,
1124                 tmp_resname,
1125                 chain_id,
1126                 res_seq_number,
1127                 res_insertion_code,
1128                 x,
1129                 y,
1130                 z,
1131                 occupancy,
1132                 b_factor,
1133                 (element != nullptr) ? element : "");
1134
1135     return n;
1136 }