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