Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / pdb2top.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, 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 "pdb2top.h"
41
42 #include <cctype>
43 #include <cmath>
44 #include <cstdio>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <string>
49 #include <vector>
50
51 #include "gromacs/fileio/pdbio.h"
52 #include "gromacs/gmxpreprocess/add_par.h"
53 #include "gromacs/gmxpreprocess/fflibutil.h"
54 #include "gromacs/gmxpreprocess/gen_ad.h"
55 #include "gromacs/gmxpreprocess/gen_vsite.h"
56 #include "gromacs/gmxpreprocess/grompp_impl.h"
57 #include "gromacs/gmxpreprocess/h_db.h"
58 #include "gromacs/gmxpreprocess/notset.h"
59 #include "gromacs/gmxpreprocess/pgutil.h"
60 #include "gromacs/gmxpreprocess/specbond.h"
61 #include "gromacs/gmxpreprocess/topdirs.h"
62 #include "gromacs/gmxpreprocess/topio.h"
63 #include "gromacs/gmxpreprocess/toputil.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/topology/residuetypes.h"
67 #include "gromacs/topology/symtab.h"
68 #include "gromacs/utility/binaryinformation.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/dir_separator.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/logger.h"
75 #include "gromacs/utility/niceheader.h"
76 #include "gromacs/utility/path.h"
77 #include "gromacs/utility/programcontext.h"
78 #include "gromacs/utility/smalloc.h"
79 #include "gromacs/utility/strconvert.h"
80 #include "gromacs/utility/strdb.h"
81 #include "gromacs/utility/stringutil.h"
82 #include "gromacs/utility/textwriter.h"
83
84 #include "hackblock.h"
85 #include "resall.h"
86
87 /* this must correspond to enum in pdb2top.h */
88 const char* hh[ehisNR] = { "HISD", "HISE", "HISH", "HIS1" };
89
90 static int missing_atoms(const PreprocessResidue* rp, int resind, t_atoms* at, int i0, int i, const gmx::MDLogger& logger)
91 {
92     int nmiss = 0;
93     for (int j = 0; j < rp->natom(); j++)
94     {
95         const char* name   = *(rp->atomname[j]);
96         bool        bFound = false;
97         for (int k = i0; k < i; k++)
98         {
99             bFound = (bFound || (gmx_strcasecmp(*(at->atomname[k]), name) == 0));
100         }
101         if (!bFound)
102         {
103             nmiss++;
104             GMX_LOG(logger.warning)
105                     .asParagraph()
106                     .appendTextFormatted("atom %s is missing in residue %s %d in the pdb file",
107                                          name,
108                                          *(at->resinfo[resind].name),
109                                          at->resinfo[resind].nr);
110             if (name[0] == 'H' || name[0] == 'h')
111             {
112                 GMX_LOG(logger.warning)
113                         .asParagraph()
114                         .appendTextFormatted(
115                                 "You might need to add atom %s to the hydrogen database of "
116                                 "building block %s "
117                                 "in the file %s.hdb (see the manual)",
118                                 name,
119                                 *(at->resinfo[resind].rtp),
120                                 rp->filebase.c_str());
121             }
122         }
123     }
124
125     return nmiss;
126 }
127
128 bool is_int(double x)
129 {
130     const double tol = 1e-4;
131     int          ix;
132
133     if (x < 0)
134     {
135         x = -x;
136     }
137     ix = gmx::roundToInt(x);
138
139     return (fabs(x - ix) < tol);
140 }
141
142 static void choose_ff_impl(const char*          ffsel,
143                            char*                forcefield,
144                            int                  ff_maxlen,
145                            char*                ffdir,
146                            int                  ffdir_maxlen,
147                            const gmx::MDLogger& logger)
148 {
149     std::vector<gmx::DataFileInfo> ffdirs = fflib_enumerate_forcefields();
150     const int                      nff    = ssize(ffdirs);
151
152     /* Replace with unix path separators */
153 #if DIR_SEPARATOR != '/'
154     for (int i = 0; i < nff; ++i)
155     {
156         std::replace(ffdirs[i].dir.begin(), ffdirs[i].dir.end(), DIR_SEPARATOR, '/');
157     }
158 #endif
159
160     /* Store the force field names in ffs */
161     std::vector<std::string> ffs;
162     ffs.reserve(ffdirs.size());
163     for (int i = 0; i < nff; ++i)
164     {
165         ffs.push_back(gmx::stripSuffixIfPresent(ffdirs[i].name, fflib_forcefield_dir_ext()));
166     }
167
168     int sel;
169     if (ffsel != nullptr)
170     {
171         sel        = -1;
172         int cwdsel = -1;
173         int nfound = 0;
174         for (int i = 0; i < nff; ++i)
175         {
176             if (ffs[i] == ffsel)
177             {
178                 /* Matching ff name */
179                 sel = i;
180                 nfound++;
181
182                 if (ffdirs[i].dir == ".")
183                 {
184                     cwdsel = i;
185                 }
186             }
187         }
188
189         if (cwdsel != -1)
190         {
191             sel = cwdsel;
192         }
193
194         if (nfound > 1)
195         {
196             if (cwdsel != -1)
197             {
198                 GMX_LOG(logger.warning)
199                         .asParagraph()
200                         .appendTextFormatted(
201                                 "Force field '%s' occurs in %d places. pdb2gmx is using the one in "
202                                 "the current directory. Use interactive selection "
203                                 "(not the -ff option) if you would prefer a different one.",
204                                 ffsel,
205                                 nfound);
206             }
207             else
208             {
209                 std::string message = gmx::formatString(
210                         "Force field '%s' occurs in %d places, but not in "
211                         "the current directory.\n"
212                         "Run without the -ff switch and select the force "
213                         "field interactively.",
214                         ffsel,
215                         nfound);
216                 GMX_THROW(gmx::InconsistentInputError(message));
217             }
218         }
219         else if (nfound == 0)
220         {
221             std::string message = gmx::formatString(
222                     "Could not find force field '%s' in current directory, "
223                     "install tree or GMXLIB path.",
224                     ffsel);
225             GMX_THROW(gmx::InconsistentInputError(message));
226         }
227     }
228     else if (nff > 1)
229     {
230         std::vector<std::string> desc;
231         desc.reserve(ffdirs.size());
232         for (int i = 0; i < nff; ++i)
233         {
234             std::string docFileName(gmx::Path::join(ffdirs[i].dir, ffdirs[i].name, fflib_forcefield_doc()));
235             // TODO: Just try to open the file with a method that does not
236             // throw/bail out with a fatal error instead of multiple checks.
237             if (gmx::File::exists(docFileName, gmx::File::returnFalseOnError))
238             {
239                 // TODO: Use a C++ API without such an intermediate/fixed-length buffer.
240                 char buf[STRLEN];
241                 /* We don't use fflib_open, because we don't want printf's */
242                 FILE* fp = gmx_ffopen(docFileName, "r");
243                 get_a_line(fp, buf, STRLEN);
244                 gmx_ffclose(fp);
245                 desc.emplace_back(buf);
246             }
247             else
248             {
249                 desc.push_back(ffs[i]);
250             }
251         }
252         /* Order force fields from the same dir alphabetically
253          * and put deprecated force fields at the end.
254          */
255         for (int i = 0; i < nff; ++i)
256         {
257             for (int j = i + 1; j < nff; ++j)
258             {
259                 if (ffdirs[i].dir == ffdirs[j].dir
260                     && ((desc[i][0] == '[' && desc[j][0] != '[')
261                         || ((desc[i][0] == '[' || desc[j][0] != '[')
262                             && gmx_strcasecmp(desc[i].c_str(), desc[j].c_str()) > 0)))
263                 {
264                     std::swap(ffdirs[i].name, ffdirs[j].name);
265                     std::swap(ffs[i], ffs[j]);
266                     std::swap(desc[i], desc[j]);
267                 }
268             }
269         }
270
271         GMX_LOG(logger.info).asParagraph().appendTextFormatted("Select the Force Field:");
272         for (int i = 0; i < nff; ++i)
273         {
274             if (i == 0 || ffdirs[i - 1].dir != ffdirs[i].dir)
275             {
276                 if (ffdirs[i].dir == ".")
277                 {
278                     GMX_LOG(logger.info)
279                             .asParagraph()
280                             .appendTextFormatted("From current directory:");
281                 }
282                 else
283                 {
284                     GMX_LOG(logger.info)
285                             .asParagraph()
286                             .appendTextFormatted("From '%s':", ffdirs[i].dir.c_str());
287                 }
288             }
289             GMX_LOG(logger.info).asParagraph().appendTextFormatted("%2d: %s", i + 1, desc[i].c_str());
290         }
291
292         sel = -1;
293         // TODO: Add a C++ API for this.
294         char  buf[STRLEN];
295         char* pret;
296         do
297         {
298             pret = fgets(buf, STRLEN, stdin);
299
300             if (pret != nullptr)
301             {
302                 sel = strtol(buf, nullptr, 10);
303                 sel--;
304             }
305         } while (pret == nullptr || (sel < 0) || (sel >= nff));
306
307         /* Check for a current limitation of the fflib code.
308          * It will always read from the first ff directory in the list.
309          * This check assumes that the order of ffs matches the order
310          * in which fflib_open searches ff library files.
311          */
312         for (int i = 0; i < sel; i++)
313         {
314             if (ffs[i] == ffs[sel])
315             {
316                 std::string message = gmx::formatString(
317                         "Can only select the first of multiple force "
318                         "field entries with directory name '%s%s' in "
319                         "the list. If you want to use the next entry, "
320                         "run pdb2gmx in a different directory, set GMXLIB "
321                         "to point to the desired force field first, and/or "
322                         "rename or move the force field directory present "
323                         "in the current working directory.",
324                         ffs[sel].c_str(),
325                         fflib_forcefield_dir_ext());
326                 GMX_THROW(gmx::NotImplementedError(message));
327             }
328         }
329     }
330     else
331     {
332         sel = 0;
333     }
334
335     if (ffs[sel].length() >= static_cast<size_t>(ff_maxlen))
336     {
337         std::string message = gmx::formatString("Length of force field name (%d) >= maxlen (%d)",
338                                                 static_cast<int>(ffs[sel].length()),
339                                                 ff_maxlen);
340         GMX_THROW(gmx::InvalidInputError(message));
341     }
342     strcpy(forcefield, ffs[sel].c_str());
343
344     std::string ffpath;
345     if (ffdirs[sel].bFromDefaultDir)
346     {
347         ffpath = ffdirs[sel].name;
348     }
349     else
350     {
351         ffpath = gmx::Path::join(ffdirs[sel].dir, ffdirs[sel].name);
352     }
353     if (ffpath.length() >= static_cast<size_t>(ffdir_maxlen))
354     {
355         std::string message = gmx::formatString("Length of force field dir (%d) >= maxlen (%d)",
356                                                 static_cast<int>(ffpath.length()),
357                                                 ffdir_maxlen);
358         GMX_THROW(gmx::InvalidInputError(message));
359     }
360     strcpy(ffdir, ffpath.c_str());
361 }
362
363 void choose_ff(const char* ffsel, char* forcefield, int ff_maxlen, char* ffdir, int ffdir_maxlen, const gmx::MDLogger& logger)
364 {
365     try
366     {
367         choose_ff_impl(ffsel, forcefield, ff_maxlen, ffdir, ffdir_maxlen, logger);
368     }
369     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
370 }
371
372 void choose_watermodel(const char* wmsel, const char* ffdir, char** watermodel, const gmx::MDLogger& logger)
373 {
374     const char* fn_watermodels = "watermodels.dat";
375     FILE*       fp;
376     char        buf[STRLEN];
377     int         nwm, sel, i;
378     char**      model;
379     char*       pret;
380
381     if (strcmp(wmsel, "none") == 0)
382     {
383         *watermodel = nullptr;
384
385         return;
386     }
387     else if (strcmp(wmsel, "select") != 0)
388     {
389         *watermodel = gmx_strdup(wmsel);
390
391         return;
392     }
393
394     std::string filename = gmx::Path::join(ffdir, fn_watermodels);
395     if (!fflib_fexist(filename))
396     {
397         GMX_LOG(logger.warning)
398                 .asParagraph()
399                 .appendTextFormatted("No file '%s' found, will not include a water model", fn_watermodels);
400         *watermodel = nullptr;
401
402         return;
403     }
404
405     fp = fflib_open(filename);
406     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Select the Water Model:");
407     nwm   = 0;
408     model = nullptr;
409     while (get_a_line(fp, buf, STRLEN))
410     {
411         srenew(model, nwm + 1);
412         snew(model[nwm], STRLEN);
413         sscanf(buf, "%s%n", model[nwm], &i);
414         if (i > 0)
415         {
416             ltrim(buf + i);
417             GMX_LOG(logger.info).asParagraph().appendTextFormatted("%2d: %s", nwm + 1, buf + i);
418             nwm++;
419         }
420         else
421         {
422             sfree(model[nwm]);
423         }
424     }
425     gmx_ffclose(fp);
426     GMX_LOG(logger.info).asParagraph().appendTextFormatted("%2d: %s", nwm + 1, "None");
427
428     sel = -1;
429     do
430     {
431         pret = fgets(buf, STRLEN, stdin);
432
433         if (pret != nullptr)
434         {
435             sel = strtol(buf, nullptr, 10);
436             sel--;
437         }
438     } while (pret == nullptr || sel < 0 || sel > nwm);
439
440     if (sel == nwm)
441     {
442         *watermodel = nullptr;
443     }
444     else
445     {
446         *watermodel = gmx_strdup(model[sel]);
447     }
448
449     for (i = 0; i < nwm; i++)
450     {
451         sfree(model[i]);
452     }
453     sfree(model);
454 }
455
456 static int name2type(t_atoms*                               at,
457                      int**                                  cgnr,
458                      gmx::ArrayRef<const PreprocessResidue> usedPpResidues,
459                      ResidueType*                           rt,
460                      const gmx::MDLogger&                   logger)
461 {
462     int    i, j, prevresind, i0, prevcg, cg, curcg;
463     char*  name;
464     bool   bNterm;
465     double qt;
466     int    nmissat;
467
468     nmissat = 0;
469
470     int resind = -1;
471     bNterm     = false;
472     i0         = 0;
473     snew(*cgnr, at->nr);
474     qt    = 0;
475     curcg = 0;
476     cg    = -1;
477
478     for (i = 0; (i < at->nr); i++)
479     {
480         prevresind = resind;
481         if (at->atom[i].resind != resind)
482         {
483             resind     = at->atom[i].resind;
484             bool bProt = rt->namedResidueHasType(*(at->resinfo[resind].name), "Protein");
485             bNterm     = bProt && (resind == 0);
486             if (resind > 0)
487             {
488                 nmissat += missing_atoms(&usedPpResidues[prevresind], prevresind, at, i0, i, logger);
489             }
490             i0 = i;
491         }
492         if (at->atom[i].m == 0)
493         {
494             qt               = 0;
495             prevcg           = cg;
496             name             = *(at->atomname[i]);
497             j                = search_jtype(usedPpResidues[resind], name, bNterm);
498             at->atom[i].type = usedPpResidues[resind].atom[j].type;
499             at->atom[i].q    = usedPpResidues[resind].atom[j].q;
500             at->atom[i].m    = usedPpResidues[resind].atom[j].m;
501             cg               = usedPpResidues[resind].cgnr[j];
502             /* A charge group number -1 signals a separate charge group
503              * for this atom.
504              */
505             if ((cg == -1) || (cg != prevcg) || (resind != prevresind))
506             {
507                 curcg++;
508             }
509         }
510         else
511         {
512             cg = -1;
513             if (is_int(qt))
514             {
515                 qt = 0;
516                 curcg++;
517             }
518             qt += at->atom[i].q;
519         }
520         (*cgnr)[i]        = curcg;
521         at->atom[i].typeB = at->atom[i].type;
522         at->atom[i].qB    = at->atom[i].q;
523         at->atom[i].mB    = at->atom[i].m;
524     }
525     nmissat += missing_atoms(&usedPpResidues[resind], resind, at, i0, i, logger);
526
527     return nmissat;
528 }
529
530 static void print_top_heavy_H(FILE* out, real mHmult)
531 {
532     if (mHmult == 2.0)
533     {
534         fprintf(out, "; Using deuterium instead of hydrogen\n\n");
535     }
536     else if (mHmult == 4.0)
537     {
538         fprintf(out, "#define HEAVY_H\n\n");
539     }
540     else if (mHmult != 1.0)
541     {
542         fprintf(stderr,
543                 "WARNING: unsupported proton mass multiplier (%g) "
544                 "in pdb2top\n",
545                 mHmult);
546     }
547 }
548
549 void print_top_comment(FILE* out, const char* filename, const char* ffdir, bool bITP)
550 {
551     char  ffdir_parent[STRLEN];
552     char* p;
553
554     try
555     {
556         gmx::TextWriter writer(out);
557         gmx::niceHeader(&writer, filename, ';');
558         writer.writeLine(gmx::formatString(";\tThis is a %s topology file", bITP ? "include" : "standalone"));
559         writer.writeLine(";");
560
561         gmx::BinaryInformationSettings settings;
562         settings.generatedByHeader(true);
563         settings.linePrefix(";\t");
564         gmx::printBinaryInformation(&writer, gmx::getProgramContext(), settings);
565     }
566     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
567
568     if (strchr(ffdir, '/') == nullptr)
569     {
570         fprintf(out, ";\tForce field was read from the standard GROMACS share directory.\n;\n\n");
571     }
572     else if (ffdir[0] == '.')
573     {
574         fprintf(out,
575                 ";\tForce field was read from current directory or a relative path - path "
576                 "added.\n;\n\n");
577     }
578     else
579     {
580         strncpy(ffdir_parent, ffdir, STRLEN - 1);
581         ffdir_parent[STRLEN - 1] = '\0'; /*make sure it is 0-terminated even for long string*/
582         p                        = strrchr(ffdir_parent, '/');
583
584         *p = '\0';
585
586         fprintf(out,
587                 ";\tForce field data was read from:\n"
588                 ";\t%s\n"
589                 ";\n"
590                 ";\tNote:\n"
591                 ";\tThis might be a non-standard force field location. When you use this topology, "
592                 "the\n"
593                 ";\tforce field must either be present in the current directory, or the location\n"
594                 ";\tspecified in the GMXLIB path variable or with the 'include' mdp file "
595                 "option.\n;\n\n",
596                 ffdir_parent);
597     }
598 }
599
600 void print_top_header(FILE* out, const char* filename, bool bITP, const char* ffdir, real mHmult)
601 {
602     const char* p;
603
604     print_top_comment(out, filename, ffdir, bITP);
605
606     print_top_heavy_H(out, mHmult);
607     fprintf(out, "; Include forcefield parameters\n");
608
609     p = strrchr(ffdir, '/');
610     p = (ffdir[0] == '.' || p == nullptr) ? ffdir : p + 1;
611
612     fprintf(out, "#include \"%s/%s\"\n\n", p, fflib_forcefield_itp());
613 }
614
615 static void print_top_posre(FILE* out, const char* pr)
616 {
617     fprintf(out, "; Include Position restraint file\n");
618     fprintf(out, "#ifdef POSRES\n");
619     fprintf(out, "#include \"%s\"\n", pr);
620     fprintf(out, "#endif\n\n");
621 }
622
623 static void print_top_water(FILE* out, const char* ffdir, const char* water)
624 {
625     const char* p;
626     char        buf[STRLEN];
627
628     fprintf(out, "; Include water topology\n");
629
630     p = strrchr(ffdir, '/');
631     p = (ffdir[0] == '.' || p == nullptr) ? ffdir : p + 1;
632     fprintf(out, "#include \"%s/%s.itp\"\n", p, water);
633
634     fprintf(out, "\n");
635     fprintf(out, "#ifdef POSRES_WATER\n");
636     fprintf(out, "; Position restraint for each water oxygen\n");
637     fprintf(out, "[ position_restraints ]\n");
638     fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
639     fprintf(out, "%4d %4d %10g %10g %10g\n", 1, 1, 1000.0, 1000.0, 1000.0);
640     fprintf(out, "#endif\n");
641     fprintf(out, "\n");
642
643     sprintf(buf, "%s/ions.itp", p);
644
645     if (fflib_fexist(buf))
646     {
647         fprintf(out, "; Include topology for ions\n");
648         fprintf(out, "#include \"%s\"\n", buf);
649         fprintf(out, "\n");
650     }
651 }
652
653 static void print_top_system(FILE* out, const char* title)
654 {
655     fprintf(out, "[ %s ]\n", dir2str(Directive::d_system));
656     fprintf(out, "; Name\n");
657     fprintf(out, "%s\n\n", title[0] ? title : "Protein");
658 }
659
660 void print_top_mols(FILE*                            out,
661                     const char*                      title,
662                     const char*                      ffdir,
663                     const char*                      water,
664                     gmx::ArrayRef<const std::string> incls,
665                     gmx::ArrayRef<const t_mols>      mols)
666 {
667
668     if (!incls.empty())
669     {
670         fprintf(out, "; Include chain topologies\n");
671         for (const auto& incl : incls)
672         {
673             fprintf(out, "#include \"%s\"\n", gmx::Path::getFilename(incl).c_str());
674         }
675         fprintf(out, "\n");
676     }
677
678     if (water)
679     {
680         print_top_water(out, ffdir, water);
681     }
682     print_top_system(out, title);
683
684     if (!mols.empty())
685     {
686         fprintf(out, "[ %s ]\n", dir2str(Directive::d_molecules));
687         fprintf(out, "; %-15s %5s\n", "Compound", "#mols");
688         for (const auto& mol : mols)
689         {
690             fprintf(out, "%-15s %5d\n", mol.name.c_str(), mol.nr);
691         }
692     }
693 }
694
695 void write_top(FILE*                                   out,
696                const char*                             pr,
697                const char*                             molname,
698                t_atoms*                                at,
699                bool                                    bRTPresname,
700                int                                     bts[],
701                gmx::ArrayRef<const InteractionsOfType> plist,
702                t_excls                                 excls[],
703                PreprocessingAtomTypes*                 atype,
704                int*                                    cgnr,
705                int                                     nrexcl)
706 /* NOTE: nrexcl is not the size of *excl! */
707 {
708     if (at && atype && cgnr)
709     {
710         fprintf(out, "[ %s ]\n", dir2str(Directive::d_moleculetype));
711         fprintf(out, "; %-15s %5s\n", "Name", "nrexcl");
712         fprintf(out, "%-15s %5d\n\n", molname ? molname : "Protein", nrexcl);
713
714         print_atoms(out, atype, at, cgnr, bRTPresname);
715         print_bondeds(out, at->nr, Directive::d_bonds, F_BONDS, bts[ebtsBONDS], plist);
716         print_bondeds(out, at->nr, Directive::d_constraints, F_CONSTR, 0, plist);
717         print_bondeds(out, at->nr, Directive::d_constraints, F_CONSTRNC, 0, plist);
718         print_bondeds(out, at->nr, Directive::d_pairs, F_LJ14, 0, plist);
719         print_excl(out, at->nr, excls);
720         print_bondeds(out, at->nr, Directive::d_angles, F_ANGLES, bts[ebtsANGLES], plist);
721         print_bondeds(out, at->nr, Directive::d_dihedrals, F_PDIHS, bts[ebtsPDIHS], plist);
722         print_bondeds(out, at->nr, Directive::d_dihedrals, F_IDIHS, bts[ebtsIDIHS], plist);
723         print_bondeds(out, at->nr, Directive::d_cmap, F_CMAP, bts[ebtsCMAP], plist);
724         print_bondeds(out, at->nr, Directive::d_polarization, F_POLARIZATION, 0, plist);
725         print_bondeds(out, at->nr, Directive::d_thole_polarization, F_THOLE_POL, 0, plist);
726         print_bondeds(out, at->nr, Directive::d_vsites2, F_VSITE2, 0, plist);
727         print_bondeds(out, at->nr, Directive::d_vsites3, F_VSITE3, 0, plist);
728         print_bondeds(out, at->nr, Directive::d_vsites3, F_VSITE3FD, 0, plist);
729         print_bondeds(out, at->nr, Directive::d_vsites3, F_VSITE3FAD, 0, plist);
730         print_bondeds(out, at->nr, Directive::d_vsites3, F_VSITE3OUT, 0, plist);
731         print_bondeds(out, at->nr, Directive::d_vsites4, F_VSITE4FD, 0, plist);
732         print_bondeds(out, at->nr, Directive::d_vsites4, F_VSITE4FDN, 0, plist);
733
734         if (pr)
735         {
736             print_top_posre(out, pr);
737         }
738     }
739 }
740
741
742 static void do_ssbonds(InteractionsOfType*                ps,
743                        t_atoms*                           atoms,
744                        gmx::ArrayRef<const DisulfideBond> ssbonds,
745                        bool                               bAllowMissing)
746 {
747     for (const auto& bond : ssbonds)
748     {
749         int ri = bond.firstResidue;
750         int rj = bond.secondResidue;
751         int ai = search_res_atom(bond.firstAtom.c_str(), ri, atoms, "special bond", bAllowMissing);
752         int aj = search_res_atom(bond.secondAtom.c_str(), rj, atoms, "special bond", bAllowMissing);
753         if ((ai == -1) || (aj == -1))
754         {
755             gmx_fatal(FARGS,
756                       "Trying to make impossible special bond (%s-%s)!",
757                       bond.firstAtom.c_str(),
758                       bond.secondAtom.c_str());
759         }
760         add_param(ps, ai, aj, {}, nullptr);
761     }
762 }
763
764 static void at2bonds(InteractionsOfType*                  psb,
765                      gmx::ArrayRef<MoleculePatchDatabase> globalPatches,
766                      t_atoms*                             atoms,
767                      gmx::ArrayRef<const gmx::RVec>       x,
768                      real                                 long_bond_dist,
769                      real                                 short_bond_dist,
770                      gmx::ArrayRef<const int>             cyclicBondsIndex,
771                      const gmx::MDLogger&                 logger)
772 {
773     real        long_bond_dist2, short_bond_dist2;
774     const char* ptr;
775
776     long_bond_dist2  = gmx::square(long_bond_dist);
777     short_bond_dist2 = gmx::square(short_bond_dist);
778
779     if (debug)
780     {
781         ptr = "bond";
782     }
783     else
784     {
785         ptr = "check";
786     }
787
788     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Making bonds...");
789     int i = 0;
790     for (int resind = 0; (resind < atoms->nres) && (i < atoms->nr); resind++)
791     {
792         /* add bonds from list of bonded interactions */
793         for (const auto& patch : globalPatches[resind].rb[ebtsBONDS].b)
794         {
795             /* Unfortunately we can not issue errors or warnings
796              * for missing atoms in bonds, as the hydrogens and terminal atoms
797              * have not been added yet.
798              */
799             int ai = search_atom(patch.ai().c_str(), i, atoms, ptr, TRUE, cyclicBondsIndex);
800             int aj = search_atom(patch.aj().c_str(), i, atoms, ptr, TRUE, cyclicBondsIndex);
801             if (ai != -1 && aj != -1)
802             {
803                 real dist2 = distance2(x[ai], x[aj]);
804                 if (dist2 > long_bond_dist2)
805
806                 {
807                     GMX_LOG(logger.warning)
808                             .asParagraph()
809                             .appendTextFormatted(
810                                     "Long Bond (%d-%d = %g nm)", ai + 1, aj + 1, std::sqrt(dist2));
811                 }
812                 else if (dist2 < short_bond_dist2)
813                 {
814                     GMX_LOG(logger.warning)
815                             .asParagraph()
816                             .appendTextFormatted(
817                                     "Short Bond (%d-%d = %g nm)", ai + 1, aj + 1, std::sqrt(dist2));
818                 }
819                 add_param(psb, ai, aj, {}, patch.s.c_str());
820             }
821         }
822         /* add bonds from list of hacks (each added atom gets a bond) */
823         while ((i < atoms->nr) && (atoms->atom[i].resind == resind))
824         {
825             for (const auto& patch : globalPatches[resind].hack)
826             {
827                 if ((patch.tp > 0 || patch.type() == MoleculePatchType::Add)
828                     && patch.a[0] == *(atoms->atomname[i]))
829                 {
830                     switch (patch.tp)
831                     {
832                         case 9:                                        /* COOH terminus */
833                             add_param(psb, i, i + 1, {}, nullptr);     /* C-O  */
834                             add_param(psb, i, i + 2, {}, nullptr);     /* C-OA */
835                             add_param(psb, i + 2, i + 3, {}, nullptr); /* OA-H */
836                             break;
837                         default:
838                             for (int k = 0; (k < patch.nr); k++)
839                             {
840                                 add_param(psb, i, i + k + 1, {}, nullptr);
841                             }
842                     }
843                 }
844             }
845             i++;
846         }
847         /* we're now at the start of the next residue */
848     }
849 }
850
851 static bool pcompar(const InteractionOfType& a, const InteractionOfType& b)
852 {
853     int d;
854
855     if (((d = a.ai() - b.ai()) != 0) || ((d = a.aj() - b.aj()) != 0))
856     {
857         return d < 0;
858     }
859     else
860     {
861         return a.interactionTypeName().length() > b.interactionTypeName().length();
862     }
863 }
864
865 static void clean_bonds(InteractionsOfType* ps, const gmx::MDLogger& logger)
866 {
867     if (ps->size() > 0)
868     {
869         /* Sort bonds */
870         for (auto& bond : ps->interactionTypes)
871         {
872             bond.sortAtomIds();
873         }
874         std::sort(ps->interactionTypes.begin(), ps->interactionTypes.end(), pcompar);
875
876         /* remove doubles, keep the first one always. */
877         int oldNumber = ps->size();
878         for (auto parm = ps->interactionTypes.begin() + 1; parm != ps->interactionTypes.end();)
879         {
880             auto prev = parm - 1;
881             if (parm->ai() == prev->ai() && parm->aj() == prev->aj())
882             {
883                 parm = ps->interactionTypes.erase(parm);
884             }
885             else
886             {
887                 ++parm;
888             }
889         }
890         GMX_LOG(logger.info)
891                 .asParagraph()
892                 .appendTextFormatted("Number of bonds was %d, now %zu", oldNumber, ps->size());
893     }
894     else
895     {
896         GMX_LOG(logger.info).asParagraph().appendTextFormatted("No bonds");
897     }
898 }
899
900 void print_sums(const t_atoms* atoms, bool bSystem, const gmx::MDLogger& logger)
901 {
902     double      m, qtot;
903     int         i;
904     const char* where;
905
906     if (bSystem)
907     {
908         where = " in system";
909     }
910     else
911     {
912         where = "";
913     }
914
915     m    = 0;
916     qtot = 0;
917     for (i = 0; (i < atoms->nr); i++)
918     {
919         m += atoms->atom[i].m;
920         qtot += atoms->atom[i].q;
921     }
922
923     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Total mass%s %.3f a.m.u.", where, m);
924     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Total charge%s %.3f e", where, qtot);
925 }
926
927 static void check_restp_type(const char* name, int t1, int t2)
928 {
929     if (t1 != t2)
930     {
931         gmx_fatal(FARGS, "Residues in one molecule have a different '%s' type: %d and %d", name, t1, t2);
932     }
933 }
934
935 static void check_restp_types(const PreprocessResidue& r0, const PreprocessResidue& r1)
936 {
937     check_restp_type("all dihedrals",
938                      static_cast<int>(r0.bKeepAllGeneratedDihedrals),
939                      static_cast<int>(r1.bKeepAllGeneratedDihedrals));
940     check_restp_type("nrexcl", r0.nrexcl, r1.nrexcl);
941     check_restp_type("HH14",
942                      static_cast<int>(r0.bGenerateHH14Interactions),
943                      static_cast<int>(r1.bGenerateHH14Interactions));
944     check_restp_type("remove dihedrals",
945                      static_cast<int>(r0.bRemoveDihedralIfWithImproper),
946                      static_cast<int>(r1.bRemoveDihedralIfWithImproper));
947
948     for (int i = 0; i < ebtsNR; i++)
949     {
950         check_restp_type(btsNames[i], r0.rb[i].type, r1.rb[i].type);
951     }
952 }
953
954 static void add_atom_to_restp(PreprocessResidue*   usedPpResidues,
955                               t_symtab*            symtab,
956                               int                  at_start,
957                               const MoleculePatch* patch)
958 {
959     /* now add them */
960     for (int k = 0; k < patch->nr; k++)
961     {
962         /* set counter in atomname */
963         std::string buf = patch->nname;
964         if (patch->nr > 1)
965         {
966             buf.append(gmx::formatString("%d", k + 1));
967         }
968         usedPpResidues->atomname.insert(usedPpResidues->atomname.begin() + at_start + 1 + k,
969                                         put_symtab(symtab, buf.c_str()));
970         usedPpResidues->atom.insert(usedPpResidues->atom.begin() + at_start + 1 + k, patch->atom.back());
971         if (patch->cgnr != NOTSET)
972         {
973             usedPpResidues->cgnr.insert(usedPpResidues->cgnr.begin() + at_start + 1 + k, patch->cgnr);
974         }
975         else
976         {
977             usedPpResidues->cgnr.insert(usedPpResidues->cgnr.begin() + at_start + 1 + k,
978                                         usedPpResidues->cgnr[at_start]);
979         }
980     }
981 }
982
983 void get_hackblocks_rtp(std::vector<MoleculePatchDatabase>*    globalPatches,
984                         std::vector<PreprocessResidue>*        usedPpResidues,
985                         gmx::ArrayRef<const PreprocessResidue> rtpFFDB,
986                         int                                    nres,
987                         t_resinfo*                             resinfo,
988                         int                                    nterpairs,
989                         t_symtab*                              symtab,
990                         gmx::ArrayRef<MoleculePatchDatabase*>  ntdb,
991                         gmx::ArrayRef<MoleculePatchDatabase*>  ctdb,
992                         gmx::ArrayRef<const int>               rn,
993                         gmx::ArrayRef<const int>               rc,
994                         bool                                   bAllowMissing,
995                         const gmx::MDLogger&                   logger)
996 {
997     char* key;
998     bool  bRM;
999
1000     globalPatches->resize(nres);
1001     usedPpResidues->clear();
1002     /* first the termini */
1003     for (int i = 0; i < nterpairs; i++)
1004     {
1005         if (rn[i] >= 0 && ntdb[i] != nullptr)
1006         {
1007             copyModificationBlocks(*ntdb[i], &globalPatches->at(rn[i]));
1008         }
1009         if (rc[i] >= 0 && ctdb[i] != nullptr)
1010         {
1011             mergeAtomAndBondModifications(*ctdb[i], &globalPatches->at(rc[i]));
1012         }
1013     }
1014
1015     /* then the whole rtp */
1016     for (int i = 0; i < nres; i++)
1017     {
1018         /* Here we allow a mismatch of one character when looking for the rtp entry.
1019          * For such a mismatch there should be only one mismatching name.
1020          * This is mainly useful for small molecules such as ions.
1021          * Note that this will usually not work for protein, DNA and RNA,
1022          * since there the residue names should be listed in residuetypes.dat
1023          * and an error will have been generated earlier in the process.
1024          */
1025         key = *resinfo[i].rtp;
1026
1027         resinfo[i].rtp = put_symtab(symtab, searchResidueDatabase(key, rtpFFDB, logger).c_str());
1028         auto res       = getDatabaseEntry(*resinfo[i].rtp, rtpFFDB);
1029         usedPpResidues->push_back(PreprocessResidue());
1030         PreprocessResidue* newentry = &usedPpResidues->back();
1031         copyPreprocessResidues(*res, newentry, symtab);
1032
1033         /* Check that we do not have different bonded types in one molecule */
1034         check_restp_types(usedPpResidues->front(), *newentry);
1035
1036         int tern = -1;
1037         for (int j = 0; j < nterpairs && tern == -1; j++)
1038         {
1039             if (i == rn[j])
1040             {
1041                 tern = j;
1042             }
1043         }
1044         int terc = -1;
1045         for (int j = 0; j < nterpairs && terc == -1; j++)
1046         {
1047             if (i == rc[j])
1048             {
1049                 terc = j;
1050             }
1051         }
1052         bRM = mergeBondedInteractionList(res->rb, globalPatches->at(i).rb, tern >= 0, terc >= 0);
1053
1054         if (bRM && ((tern >= 0 && ntdb[tern] == nullptr) || (terc >= 0 && ctdb[terc] == nullptr)))
1055         {
1056             const char* errString =
1057                     "There is a dangling bond at at least one of the terminal ends and the force "
1058                     "field does not provide terminal entries or files. Fix your terminal "
1059                     "residues so that they match the residue database (.rtp) entries, or provide "
1060                     "terminal database entries (.tdb).";
1061             if (bAllowMissing)
1062             {
1063                 GMX_LOG(logger.warning).asParagraph().appendTextFormatted("%s", errString);
1064             }
1065             else
1066             {
1067                 gmx_fatal(FARGS, "%s", errString);
1068             }
1069         }
1070         else if (bRM && ((tern >= 0 && ntdb[tern]->nhack() == 0) || (terc >= 0 && ctdb[terc]->nhack() == 0)))
1071         {
1072             const char* errString =
1073                     "There is a dangling bond at at least one of the terminal ends. Fix your "
1074                     "coordinate file, add a new terminal database entry (.tdb), or select the "
1075                     "proper existing terminal entry.";
1076             if (bAllowMissing)
1077             {
1078                 GMX_LOG(logger.warning).asParagraph().appendTextFormatted("%s", errString);
1079             }
1080             else
1081             {
1082                 gmx_fatal(FARGS, "%s", errString);
1083             }
1084         }
1085     }
1086
1087     /* Apply patchs to t_restp entries
1088        i.e. add's and deletes from termini database will be
1089        added to/removed from residue topology
1090        we'll do this on one big dirty loop, so it won't make easy reading! */
1091     for (auto modifiedResidue = globalPatches->begin(); modifiedResidue != globalPatches->end();
1092          modifiedResidue++)
1093     {
1094         const int          pos    = std::distance(globalPatches->begin(), modifiedResidue);
1095         PreprocessResidue* posres = &usedPpResidues->at(pos);
1096         for (auto patch = modifiedResidue->hack.begin(); patch != modifiedResidue->hack.end(); patch++)
1097         {
1098             if (patch->nr != 0)
1099             {
1100                 /* find atom in restp */
1101                 auto found = std::find_if(
1102                         posres->atomname.begin(), posres->atomname.end(), [&patch](char** name) {
1103                             return (patch->oname.empty() && patch->a[0] == *name)
1104                                    || (patch->oname == *name);
1105                         });
1106
1107                 if (found == posres->atomname.end())
1108                 {
1109                     /* If we are doing an atom rename only, we don't need
1110                      * to generate a fatal error if the old name is not found
1111                      * in the rtp.
1112                      */
1113                     /* Deleting can happen also only on the input atoms,
1114                      * not necessarily always on the rtp entry.
1115                      */
1116                     if (patch->type() == MoleculePatchType::Add)
1117                     {
1118                         gmx_fatal(FARGS,
1119                                   "atom %s not found in buiding block %d%s "
1120                                   "while combining tdb and rtp",
1121                                   patch->oname.empty() ? patch->a[0].c_str() : patch->oname.c_str(),
1122                                   pos + 1,
1123                                   *resinfo[pos].rtp);
1124                     }
1125                 }
1126                 else
1127                 {
1128                     int l = std::distance(posres->atomname.begin(), found);
1129                     switch (patch->type())
1130                     {
1131                         case MoleculePatchType::Add:
1132                         {
1133                             /* we're adding: */
1134                             add_atom_to_restp(posres, symtab, l, &(*patch));
1135                             break;
1136                         }
1137                         case MoleculePatchType::Delete:
1138                         { /* we're deleting */
1139                             posres->atom.erase(posres->atom.begin() + l);
1140                             posres->atomname.erase(posres->atomname.begin() + l);
1141                             posres->cgnr.erase(posres->cgnr.begin() + l);
1142                             break;
1143                         }
1144                         case MoleculePatchType::Replace:
1145                         {
1146                             /* we're replacing */
1147                             posres->atom[l]     = patch->atom.back();
1148                             posres->atomname[l] = put_symtab(symtab, patch->nname.c_str());
1149                             if (patch->cgnr != NOTSET)
1150                             {
1151                                 posres->cgnr[l] = patch->cgnr;
1152                             }
1153                             break;
1154                         }
1155                     }
1156                 }
1157             }
1158         }
1159     }
1160 }
1161
1162 static bool atomname_cmp_nr(const char* anm, const MoleculePatch* patch, int* nr)
1163 {
1164
1165     if (patch->nr == 1)
1166     {
1167         *nr = 0;
1168
1169         return (gmx::equalCaseInsensitive(anm, patch->nname));
1170     }
1171     else
1172     {
1173         if (isdigit(anm[strlen(anm) - 1]))
1174         {
1175             *nr = anm[strlen(anm) - 1] - '0';
1176         }
1177         else
1178         {
1179             *nr = 0;
1180         }
1181         if (*nr <= 0 || *nr > patch->nr)
1182         {
1183             return FALSE;
1184         }
1185         else
1186         {
1187             return (strlen(anm) == patch->nname.length() + 1
1188                     && gmx_strncasecmp(anm, patch->nname.c_str(), patch->nname.length()) == 0);
1189         }
1190     }
1191 }
1192
1193 static bool match_atomnames_with_rtp_atom(t_atoms*                     pdba,
1194                                           gmx::ArrayRef<gmx::RVec>     x,
1195                                           t_symtab*                    symtab,
1196                                           int                          atind,
1197                                           PreprocessResidue*           localPpResidue,
1198                                           const MoleculePatchDatabase& singlePatch,
1199                                           bool                         bVerbose,
1200                                           const gmx::MDLogger&         logger)
1201 {
1202     int   resnr;
1203     char* oldnm;
1204     int   anmnr;
1205     bool  bDeleted;
1206
1207     oldnm = *pdba->atomname[atind];
1208     resnr = pdba->resinfo[pdba->atom[atind].resind].nr;
1209
1210     bDeleted = FALSE;
1211     for (auto patch = singlePatch.hack.begin(); patch != singlePatch.hack.end(); patch++)
1212     {
1213         if (patch->type() == MoleculePatchType::Replace && gmx::equalCaseInsensitive(oldnm, patch->oname))
1214         {
1215             /* This is a replace entry. */
1216             /* Check if we are not replacing a replaced atom. */
1217             bool bReplaceReplace = false;
1218             for (auto selfPatch = singlePatch.hack.begin(); selfPatch != singlePatch.hack.end(); selfPatch++)
1219             {
1220                 if (patch != selfPatch && selfPatch->type() == MoleculePatchType::Replace
1221                     && gmx::equalCaseInsensitive(selfPatch->nname, patch->oname))
1222                 {
1223                     /* The replace in patch replaces an atom that
1224                      * was already replaced in selfPatch, we do not want
1225                      * second or higher level replaces at this stage.
1226                      */
1227                     bReplaceReplace = true;
1228                 }
1229             }
1230             if (bReplaceReplace)
1231             {
1232                 /* Skip this replace. */
1233                 continue;
1234             }
1235
1236             /* This atom still has the old name, rename it */
1237             std::string newnm = patch->nname;
1238             auto        found = std::find_if(
1239                     localPpResidue->atomname.begin(),
1240                     localPpResidue->atomname.end(),
1241                     [&newnm](char** name) { return gmx::equalCaseInsensitive(newnm, *name); });
1242             if (found == localPpResidue->atomname.end())
1243             {
1244                 /* The new name is not present in the rtp.
1245                  * We need to apply the replace also to the rtp entry.
1246                  */
1247
1248                 /* We need to find the add hack that can add this atom
1249                  * to find out after which atom it should be added.
1250                  */
1251                 bool bFoundInAdd = false;
1252                 for (auto rtpModification = singlePatch.hack.begin();
1253                      rtpModification != singlePatch.hack.end();
1254                      rtpModification++)
1255                 {
1256                     int         k = std::distance(localPpResidue->atomname.begin(), found);
1257                     std::string start_at;
1258                     if (rtpModification->type() == MoleculePatchType::Add
1259                         && atomname_cmp_nr(newnm.c_str(), &(*rtpModification), &anmnr))
1260                     {
1261                         if (anmnr <= 1)
1262                         {
1263                             start_at = singlePatch.hack[k].a[0];
1264                         }
1265                         else
1266                         {
1267                             start_at = gmx::formatString(
1268                                     "%s%d", singlePatch.hack[k].nname.c_str(), anmnr - 1);
1269                         }
1270                         auto found2 = std::find_if(localPpResidue->atomname.begin(),
1271                                                    localPpResidue->atomname.end(),
1272                                                    [&start_at](char** name) {
1273                                                        return gmx::equalCaseInsensitive(start_at, *name);
1274                                                    });
1275                         if (found2 == localPpResidue->atomname.end())
1276                         {
1277                             gmx_fatal(FARGS,
1278                                       "Could not find atom '%s' in residue building block '%s' to "
1279                                       "add atom '%s' to",
1280                                       start_at.c_str(),
1281                                       localPpResidue->resname.c_str(),
1282                                       newnm.c_str());
1283                         }
1284                         /* We can add the atom after atom start_nr */
1285                         add_atom_to_restp(localPpResidue,
1286                                           symtab,
1287                                           std::distance(localPpResidue->atomname.begin(), found2),
1288                                           &(*patch));
1289
1290                         bFoundInAdd = true;
1291                     }
1292                 }
1293
1294                 if (!bFoundInAdd)
1295                 {
1296                     gmx_fatal(FARGS,
1297                               "Could not find an 'add' entry for atom named '%s' corresponding to "
1298                               "the 'replace' entry from atom name '%s' to '%s' for tdb or hdb "
1299                               "database of residue type '%s'",
1300                               newnm.c_str(),
1301                               patch->oname.c_str(),
1302                               patch->nname.c_str(),
1303                               localPpResidue->resname.c_str());
1304                 }
1305             }
1306
1307             if (bVerbose)
1308             {
1309                 GMX_LOG(logger.info)
1310                         .asParagraph()
1311                         .appendTextFormatted("Renaming atom '%s' in residue '%s' %d to '%s'",
1312                                              oldnm,
1313                                              localPpResidue->resname.c_str(),
1314                                              resnr,
1315                                              newnm.c_str());
1316             }
1317             /* Rename the atom in pdba */
1318             pdba->atomname[atind] = put_symtab(symtab, newnm.c_str());
1319         }
1320         else if (patch->type() == MoleculePatchType::Delete
1321                  && gmx::equalCaseInsensitive(oldnm, patch->oname))
1322         {
1323             /* This is a delete entry, check if this atom is present
1324              * in the rtp entry of this residue.
1325              */
1326             auto found3 = std::find_if(
1327                     localPpResidue->atomname.begin(),
1328                     localPpResidue->atomname.end(),
1329                     [&oldnm](char** name) { return gmx::equalCaseInsensitive(oldnm, *name); });
1330             if (found3 == localPpResidue->atomname.end())
1331             {
1332                 /* This atom is not present in the rtp entry,
1333                  * delete is now from the input pdba.
1334                  */
1335                 if (bVerbose)
1336                 {
1337                     GMX_LOG(logger.info)
1338                             .asParagraph()
1339                             .appendTextFormatted("Deleting atom '%s' in residue '%s' %d",
1340                                                  oldnm,
1341                                                  localPpResidue->resname.c_str(),
1342                                                  resnr);
1343                 }
1344                 /* We should free the atom name,
1345                  * but it might be used multiple times in the symtab.
1346                  * sfree(pdba->atomname[atind]);
1347                  */
1348                 for (int k = atind + 1; k < pdba->nr; k++)
1349                 {
1350                     pdba->atom[k - 1]     = pdba->atom[k];
1351                     pdba->atomname[k - 1] = pdba->atomname[k];
1352                     copy_rvec(x[k], x[k - 1]);
1353                 }
1354                 pdba->nr--;
1355                 bDeleted = true;
1356             }
1357         }
1358     }
1359
1360     return bDeleted;
1361 }
1362
1363 void match_atomnames_with_rtp(gmx::ArrayRef<PreprocessResidue>     usedPpResidues,
1364                               gmx::ArrayRef<MoleculePatchDatabase> globalPatches,
1365                               t_atoms*                             pdba,
1366                               t_symtab*                            symtab,
1367                               gmx::ArrayRef<gmx::RVec>             x,
1368                               bool                                 bVerbose,
1369                               const gmx::MDLogger&                 logger)
1370 {
1371     for (int i = 0; i < pdba->nr; i++)
1372     {
1373         const char*        oldnm          = *pdba->atomname[i];
1374         PreprocessResidue* localPpResidue = &usedPpResidues[pdba->atom[i].resind];
1375         auto               found          = std::find_if(
1376                 localPpResidue->atomname.begin(), localPpResidue->atomname.end(), [&oldnm](char** name) {
1377                     return gmx::equalCaseInsensitive(oldnm, *name);
1378                 });
1379         if (found == localPpResidue->atomname.end())
1380         {
1381             /* Not found yet, check if we have to rename this atom */
1382             if (match_atomnames_with_rtp_atom(
1383                         pdba, x, symtab, i, localPpResidue, globalPatches[pdba->atom[i].resind], bVerbose, logger))
1384             {
1385                 /* We deleted this atom, decrease the atom counter by 1. */
1386                 i--;
1387             }
1388         }
1389     }
1390 }
1391
1392 #define NUM_CMAP_ATOMS 5
1393 static void gen_cmap(InteractionsOfType*                    psb,
1394                      gmx::ArrayRef<const PreprocessResidue> usedPpResidues,
1395                      t_atoms*                               atoms,
1396                      gmx::ArrayRef<const int>               cyclicBondsIndex,
1397                      const gmx::MDLogger&                   logger)
1398 {
1399     int         residx;
1400     const char* ptr;
1401     t_resinfo*  resinfo = atoms->resinfo;
1402     int         nres    = atoms->nres;
1403     int         cmap_atomid[NUM_CMAP_ATOMS];
1404     int         cmap_chainnum = -1;
1405
1406     if (debug)
1407     {
1408         ptr = "cmap";
1409     }
1410     else
1411     {
1412         ptr = "check";
1413     }
1414
1415     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Making cmap torsions...");
1416     int i = 0;
1417     /* Most cmap entries use the N atom from the next residue, so the last
1418      * residue should not have its CMAP entry in that case, but for things like
1419      * dipeptides we sometimes define a complete CMAP entry inside a residue,
1420      * and in this case we need to process everything through the last residue.
1421      */
1422     for (residx = 0; residx < nres; residx++)
1423     {
1424         /* Add CMAP terms from the list of CMAP interactions */
1425         for (const auto& b : usedPpResidues[residx].rb[ebtsCMAP].b)
1426         {
1427             bool bAddCMAP = true;
1428             /* Loop over atoms in a candidate CMAP interaction and
1429              * check that they exist, are from the same chain and are
1430              * from residues labelled as protein. */
1431             for (int k = 0; k < NUM_CMAP_ATOMS && bAddCMAP; k++)
1432             {
1433                 /* Assign the pointer to the name of the next reference atom.
1434                  * This can use -/+ labels to refer to previous/next residue.
1435                  */
1436                 const char* pname = b.a[k].c_str();
1437                 /* Skip this CMAP entry if it refers to residues before the
1438                  * first or after the last residue.
1439                  */
1440                 if ((cyclicBondsIndex.empty() && ((strchr(pname, '-') != nullptr) && (residx == 0)))
1441                     || ((strchr(pname, '+') != nullptr) && (residx == nres - 1)))
1442                 {
1443                     bAddCMAP = false;
1444                     break;
1445                 }
1446
1447                 cmap_atomid[k] = search_atom(pname, i, atoms, ptr, TRUE, cyclicBondsIndex);
1448                 bAddCMAP       = bAddCMAP && (cmap_atomid[k] != -1);
1449                 if (!bAddCMAP)
1450                 {
1451                     /* This break is necessary, because cmap_atomid[k]
1452                      * == -1 cannot be safely used as an index
1453                      * into the atom array. */
1454                     break;
1455                 }
1456                 int this_residue_index = atoms->atom[cmap_atomid[k]].resind;
1457                 if (0 == k)
1458                 {
1459                     cmap_chainnum = resinfo[this_residue_index].chainnum;
1460                 }
1461                 else
1462                 {
1463                     /* Does the residue for this atom have the same
1464                      * chain number as the residues for previous
1465                      * atoms? */
1466                     bAddCMAP = bAddCMAP && cmap_chainnum == resinfo[this_residue_index].chainnum;
1467                 }
1468                 /* Here we used to check that the residuetype was protein and
1469                  * disable bAddCMAP if that was not the case. However, some
1470                  * special residues (say, alanine dipeptides) might not adhere
1471                  * to standard naming, and if we start calling them normal
1472                  * protein residues the user will be bugged to select termini.
1473                  *
1474                  * Instead, I believe that the right course of action is to
1475                  * keep the CMAP interaction if it is present in the RTP file
1476                  * and we correctly identified all atoms (which is the case
1477                  * if we got here).
1478                  */
1479             }
1480
1481             if (bAddCMAP)
1482             {
1483                 add_cmap_param(psb,
1484                                cmap_atomid[0],
1485                                cmap_atomid[1],
1486                                cmap_atomid[2],
1487                                cmap_atomid[3],
1488                                cmap_atomid[4],
1489                                b.s.c_str());
1490             }
1491         }
1492
1493         if (residx < nres - 1)
1494         {
1495             while (atoms->atom[i].resind < residx + 1)
1496             {
1497                 i++;
1498             }
1499         }
1500     }
1501     /* Start the next residue */
1502 }
1503
1504 static void scrub_charge_groups(int* cgnr, int natoms)
1505 {
1506     int i;
1507
1508     for (i = 0; i < natoms; i++)
1509     {
1510         cgnr[i] = i + 1;
1511     }
1512 }
1513
1514
1515 void pdb2top(FILE*                                  top_file,
1516              const char*                            posre_fn,
1517              const char*                            molname,
1518              t_atoms*                               atoms,
1519              std::vector<gmx::RVec>*                x,
1520              PreprocessingAtomTypes*                atype,
1521              t_symtab*                              tab,
1522              gmx::ArrayRef<const PreprocessResidue> rtpFFDB,
1523              gmx::ArrayRef<PreprocessResidue>       usedPpResidues,
1524              gmx::ArrayRef<MoleculePatchDatabase>   globalPatches,
1525              bool                                   bAllowMissing,
1526              bool                                   bVsites,
1527              bool                                   bVsiteAromatics,
1528              const char*                            ffdir,
1529              real                                   mHmult,
1530              gmx::ArrayRef<const DisulfideBond>     ssbonds,
1531              real                                   long_bond_dist,
1532              real                                   short_bond_dist,
1533              bool                                   bDeuterate,
1534              bool                                   bChargeGroups,
1535              bool                                   bCmap,
1536              bool                                   bRenumRes,
1537              bool                                   bRTPresname,
1538              gmx::ArrayRef<const int>               cyclicBondsIndex,
1539              const gmx::MDLogger&                   logger)
1540 {
1541     std::array<InteractionsOfType, F_NRE> plist;
1542     t_excls*                              excls;
1543     int*                                  cgnr;
1544     int*                                  vsite_type;
1545     int                                   i, nmissat;
1546     int                                   bts[ebtsNR];
1547
1548     ResidueType rt;
1549
1550     /* Make bonds */
1551     at2bonds(&(plist[F_BONDS]), globalPatches, atoms, *x, long_bond_dist, short_bond_dist, cyclicBondsIndex, logger);
1552
1553     /* specbonds: disulphide bonds & heme-his */
1554     do_ssbonds(&(plist[F_BONDS]), atoms, ssbonds, bAllowMissing);
1555
1556     nmissat = name2type(atoms, &cgnr, usedPpResidues, &rt, logger);
1557     if (nmissat)
1558     {
1559         if (bAllowMissing)
1560         {
1561             GMX_LOG(logger.warning)
1562                     .asParagraph()
1563                     .appendTextFormatted("There were %d missing atoms in molecule %s", nmissat, molname);
1564         }
1565         else
1566         {
1567             gmx_fatal(FARGS,
1568                       "There were %d missing atoms in molecule %s, if you want to use this "
1569                       "incomplete topology anyhow, use the option -missing",
1570                       nmissat,
1571                       molname);
1572         }
1573     }
1574
1575     /* Cleanup bonds (sort and rm doubles) */
1576     clean_bonds(&(plist[F_BONDS]), logger);
1577
1578     snew(vsite_type, atoms->nr);
1579     for (i = 0; i < atoms->nr; i++)
1580     {
1581         vsite_type[i] = NOTSET;
1582     }
1583     if (bVsites)
1584     {
1585         if (bVsiteAromatics)
1586         {
1587             GMX_LOG(logger.info)
1588                     .asParagraph()
1589                     .appendTextFormatted(
1590                             "The conversion of aromatic rings into virtual sites is deprecated "
1591                             "and may be removed in a future version of GROMACS");
1592         }
1593         /* determine which atoms will be vsites and add dummy masses
1594            also renumber atom numbers in plist[0..F_NRE]! */
1595         do_vsites(rtpFFDB, atype, atoms, tab, x, plist, &vsite_type, &cgnr, mHmult, bVsiteAromatics, ffdir);
1596     }
1597
1598     /* Make Angles and Dihedrals */
1599     GMX_LOG(logger.info)
1600             .asParagraph()
1601             .appendTextFormatted("Generating angles, dihedrals and pairs...");
1602     snew(excls, atoms->nr);
1603     gen_pad(atoms, usedPpResidues, plist, excls, globalPatches, bAllowMissing, cyclicBondsIndex);
1604
1605     /* Make CMAP */
1606     if (bCmap)
1607     {
1608         gen_cmap(&(plist[F_CMAP]), usedPpResidues, atoms, cyclicBondsIndex, logger);
1609         if (plist[F_CMAP].size() > 0)
1610         {
1611             GMX_LOG(logger.info)
1612                     .asParagraph()
1613                     .appendTextFormatted("There are %4zu cmap torsion pairs", plist[F_CMAP].size());
1614         }
1615     }
1616
1617     /* set mass of all remaining hydrogen atoms */
1618     if (mHmult != 1.0)
1619     {
1620         do_h_mass(&(plist[F_BONDS]), vsite_type, atoms, mHmult, bDeuterate);
1621     }
1622     sfree(vsite_type);
1623
1624     /* Cleanup bonds (sort and rm doubles) */
1625     /* clean_bonds(&(plist[F_BONDS]));*/
1626
1627     GMX_LOG(logger.info)
1628             .asParagraph()
1629             .appendTextFormatted(
1630                     "There are %4zu dihedrals, %4zu impropers, %4zu angles\n"
1631                     "          %4zu pairs,     %4zu bonds and  %4zu virtual sites",
1632                     plist[F_PDIHS].size(),
1633                     plist[F_IDIHS].size(),
1634                     plist[F_ANGLES].size(),
1635                     plist[F_LJ14].size(),
1636                     plist[F_BONDS].size(),
1637                     plist[F_VSITE2].size() + plist[F_VSITE3].size() + plist[F_VSITE3FD].size()
1638                             + plist[F_VSITE3FAD].size() + plist[F_VSITE3OUT].size()
1639                             + plist[F_VSITE4FD].size() + plist[F_VSITE4FDN].size());
1640
1641     print_sums(atoms, FALSE, logger);
1642
1643     if (!bChargeGroups)
1644     {
1645         scrub_charge_groups(cgnr, atoms->nr);
1646     }
1647
1648     if (bRenumRes)
1649     {
1650         for (i = 0; i < atoms->nres; i++)
1651         {
1652             atoms->resinfo[i].nr = i + 1;
1653             atoms->resinfo[i].ic = ' ';
1654         }
1655     }
1656
1657     if (top_file)
1658     {
1659         GMX_LOG(logger.info).asParagraph().appendTextFormatted("Writing topology");
1660         /* We can copy the bonded types from the first restp,
1661          * since the types have to be identical for all residues in one molecule.
1662          */
1663         for (i = 0; i < ebtsNR; i++)
1664         {
1665             bts[i] = usedPpResidues[0].rb[i].type;
1666         }
1667         write_top(top_file,
1668                   posre_fn,
1669                   molname,
1670                   atoms,
1671                   bRTPresname,
1672                   bts,
1673                   plist,
1674                   excls,
1675                   atype,
1676                   cgnr,
1677                   usedPpResidues[0].nrexcl);
1678     }
1679
1680
1681     /* we should clean up hb and restp here, but that is a *L*O*T* of work! */
1682     sfree(cgnr);
1683     for (i = 0; i < atoms->nr; i++)
1684     {
1685         sfree(excls[i].e);
1686     }
1687     sfree(excls);
1688 }