b8567d2a5328d33cc1a0fcf94b424ee5f18fdf67
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / grompp.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 "grompp.h"
41
42 #include <cerrno>
43 #include <climits>
44 #include <cmath>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <memory>
49 #include <vector>
50
51 #include <sys/types.h>
52
53 #include "gromacs/applied_forces/awh/read_params.h"
54 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/ewald/ewald_utils.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/fft/calcgrid.h"
58 #include "gromacs/fileio/confio.h"
59 #include "gromacs/fileio/enxio.h"
60 #include "gromacs/fileio/tpxio.h"
61 #include "gromacs/fileio/trxio.h"
62 #include "gromacs/fileio/warninp.h"
63 #include "gromacs/gmxpreprocess/add_par.h"
64 #include "gromacs/gmxpreprocess/convparm.h"
65 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
66 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
67 #include "gromacs/gmxpreprocess/grompp_impl.h"
68 #include "gromacs/gmxpreprocess/notset.h"
69 #include "gromacs/gmxpreprocess/readir.h"
70 #include "gromacs/gmxpreprocess/tomorse.h"
71 #include "gromacs/gmxpreprocess/topio.h"
72 #include "gromacs/gmxpreprocess/toputil.h"
73 #include "gromacs/gmxpreprocess/vsite_parm.h"
74 #include "gromacs/imd/imd.h"
75 #include "gromacs/math/functions.h"
76 #include "gromacs/math/invertmatrix.h"
77 #include "gromacs/math/units.h"
78 #include "gromacs/math/vec.h"
79 #include "gromacs/mdlib/calc_verletbuf.h"
80 #include "gromacs/mdlib/compute_io.h"
81 #include "gromacs/mdlib/constr.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/perf_est.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdrun/mdmodules.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/mdtypes/md_enums.h"
88 #include "gromacs/mdtypes/nblist.h"
89 #include "gromacs/mdtypes/state.h"
90 #include "gromacs/pbcutil/boxutilities.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/pulling/pull.h"
93 #include "gromacs/random/seed.h"
94 #include "gromacs/topology/ifunc.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/symtab.h"
97 #include "gromacs/topology/topology.h"
98 #include "gromacs/trajectory/trajectoryframe.h"
99 #include "gromacs/utility/arraysize.h"
100 #include "gromacs/utility/cstringutil.h"
101 #include "gromacs/utility/exceptions.h"
102 #include "gromacs/utility/fatalerror.h"
103 #include "gromacs/utility/filestream.h"
104 #include "gromacs/utility/futil.h"
105 #include "gromacs/utility/gmxassert.h"
106 #include "gromacs/utility/keyvaluetreebuilder.h"
107 #include "gromacs/utility/listoflists.h"
108 #include "gromacs/utility/logger.h"
109 #include "gromacs/utility/loggerbuilder.h"
110 #include "gromacs/utility/mdmodulenotification.h"
111 #include "gromacs/utility/smalloc.h"
112 #include "gromacs/utility/snprintf.h"
113
114 /* TODO The implementation details should move to their own source file. */
115 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int>  atoms,
116                                      gmx::ArrayRef<const real> params,
117                                      const std::string&        name) :
118     atoms_(atoms.begin(), atoms.end()),
119     interactionTypeName_(name)
120 {
121     GMX_RELEASE_ASSERT(
122             params.size() <= forceParam_.size(),
123             gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM)
124                     .c_str());
125     auto forceParamIt = forceParam_.begin();
126     for (const auto param : params)
127     {
128         *forceParamIt++ = param;
129     }
130     while (forceParamIt != forceParam_.end())
131     {
132         *forceParamIt++ = NOTSET;
133     }
134 }
135
136 const int& InteractionOfType::ai() const
137 {
138     GMX_RELEASE_ASSERT(!atoms_.empty(), "Need to have at least one atom set");
139     return atoms_[0];
140 }
141
142 const int& InteractionOfType::aj() const
143 {
144     GMX_RELEASE_ASSERT(atoms_.size() > 1, "Need to have at least two atoms set");
145     return atoms_[1];
146 }
147
148 const int& InteractionOfType::ak() const
149 {
150     GMX_RELEASE_ASSERT(atoms_.size() > 2, "Need to have at least three atoms set");
151     return atoms_[2];
152 }
153
154 const int& InteractionOfType::al() const
155 {
156     GMX_RELEASE_ASSERT(atoms_.size() > 3, "Need to have at least four atoms set");
157     return atoms_[3];
158 }
159
160 const int& InteractionOfType::am() const
161 {
162     GMX_RELEASE_ASSERT(atoms_.size() > 4, "Need to have at least five atoms set");
163     return atoms_[4];
164 }
165
166 const real& InteractionOfType::c0() const
167 {
168     return forceParam_[0];
169 }
170
171 const real& InteractionOfType::c1() const
172 {
173     return forceParam_[1];
174 }
175
176 const real& InteractionOfType::c2() const
177 {
178     return forceParam_[2];
179 }
180
181 const std::string& InteractionOfType::interactionTypeName() const
182 {
183     return interactionTypeName_;
184 }
185
186 void InteractionOfType::sortBondAtomIds()
187 {
188     if (aj() < ai())
189     {
190         std::swap(atoms_[0], atoms_[1]);
191     }
192 }
193
194 void InteractionOfType::sortAngleAtomIds()
195 {
196     if (ak() < ai())
197     {
198         std::swap(atoms_[0], atoms_[2]);
199     }
200 }
201
202 void InteractionOfType::sortDihedralAtomIds()
203 {
204     if (al() < ai())
205     {
206         std::swap(atoms_[0], atoms_[3]);
207         std::swap(atoms_[1], atoms_[2]);
208     }
209 }
210
211 void InteractionOfType::sortAtomIds()
212 {
213     if (isBond())
214     {
215         sortBondAtomIds();
216     }
217     else if (isAngle())
218     {
219         sortAngleAtomIds();
220     }
221     else if (isDihedral())
222     {
223         sortDihedralAtomIds();
224     }
225     else
226     {
227         GMX_THROW(gmx::InternalError(
228                 "Can not sort parameters other than bonds, angles or dihedrals"));
229     }
230 };
231
232 void InteractionOfType::setForceParameter(int pos, real value)
233 {
234     GMX_RELEASE_ASSERT(pos < MAXFORCEPARAM,
235                        "Can't set parameter beyond the maximum number of parameters");
236     forceParam_[pos] = value;
237 }
238
239 void MoleculeInformation::initMolInfo()
240 {
241     init_block(&mols);
242     excls.clear();
243     init_t_atoms(&atoms, 0, FALSE);
244 }
245
246 void MoleculeInformation::partialCleanUp()
247 {
248     done_block(&mols);
249 }
250
251 void MoleculeInformation::fullCleanUp()
252 {
253     done_atom(&atoms);
254     done_block(&mols);
255 }
256
257 static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
258 {
259     int n = 0;
260     /* For all the molecule types */
261     for (auto& mol : mols)
262     {
263         n += mol.interactions[ifunc].size();
264         mol.interactions[ifunc].interactionTypes.clear();
265     }
266     return n;
267 }
268
269 static int check_atom_names(const char*          fn1,
270                             const char*          fn2,
271                             gmx_mtop_t*          mtop,
272                             const t_atoms*       at,
273                             const gmx::MDLogger& logger)
274 {
275     int      m, i, j, nmismatch;
276     t_atoms* tat;
277
278     constexpr int c_maxNumberOfMismatches = 20;
279
280     if (mtop->natoms != at->nr)
281     {
282         gmx_incons("comparing atom names");
283     }
284
285     nmismatch = 0;
286     i         = 0;
287     for (const gmx_molblock_t& molb : mtop->molblock)
288     {
289         tat = &mtop->moltype[molb.type].atoms;
290         for (m = 0; m < molb.nmol; m++)
291         {
292             for (j = 0; j < tat->nr; j++)
293             {
294                 if (strcmp(*(tat->atomname[j]), *(at->atomname[i])) != 0)
295                 {
296                     if (nmismatch < c_maxNumberOfMismatches)
297                     {
298                         GMX_LOG(logger.warning)
299                                 .asParagraph()
300                                 .appendTextFormatted(
301                                         "atom name %d in %s and %s does not match (%s - %s)",
302                                         i + 1,
303                                         fn1,
304                                         fn2,
305                                         *(tat->atomname[j]),
306                                         *(at->atomname[i]));
307                     }
308                     else if (nmismatch == c_maxNumberOfMismatches)
309                     {
310                         GMX_LOG(logger.warning)
311                                 .asParagraph()
312                                 .appendTextFormatted("(more than %d non-matching atom names)",
313                                                      c_maxNumberOfMismatches);
314                     }
315                     nmismatch++;
316                 }
317                 i++;
318             }
319         }
320     }
321
322     return nmismatch;
323 }
324
325 static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
326 {
327     /* This check is not intended to ensure accurate integration,
328      * rather it is to signal mistakes in the mdp settings.
329      * A common mistake is to forget to turn on constraints
330      * for MD after energy minimization with flexible bonds.
331      * This check can also detect too large time steps for flexible water
332      * models, but such errors will often be masked by the constraints
333      * mdp options, which turns flexible water into water with bond constraints,
334      * but without an angle constraint. Unfortunately such incorrect use
335      * of water models can not easily be detected without checking
336      * for specific model names.
337      *
338      * The stability limit of leap-frog or velocity verlet is 4.44 steps
339      * per oscillational period.
340      * But accurate bonds distributions are lost far before that limit.
341      * To allow relatively common schemes (although not common with Gromacs)
342      * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
343      * we set the note limit to 10.
344      */
345     int  min_steps_warn = 5;
346     int  min_steps_note = 10;
347     int  ftype;
348     int  i, a1, a2, w_a1, w_a2, j;
349     real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
350     bool bFound, bWater, bWarn;
351
352     /* Get the interaction parameters */
353     gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
354
355     twopi2 = gmx::square(2 * M_PI);
356
357     limit2 = gmx::square(min_steps_note * dt);
358
359     w_a1 = w_a2 = -1;
360     w_period2   = -1.0;
361
362     const gmx_moltype_t* w_moltype = nullptr;
363     for (const gmx_moltype_t& moltype : mtop->moltype)
364     {
365         const t_atom*           atom  = moltype.atoms.atom;
366         const InteractionLists& ilist = moltype.ilist;
367         const InteractionList&  ilc   = ilist[F_CONSTR];
368         const InteractionList&  ils   = ilist[F_SETTLE];
369         for (ftype = 0; ftype < F_NRE; ftype++)
370         {
371             if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
372             {
373                 continue;
374             }
375
376             const InteractionList& ilb = ilist[ftype];
377             for (i = 0; i < ilb.size(); i += 3)
378             {
379                 fc = ip[ilb.iatoms[i]].harmonic.krA;
380                 re = ip[ilb.iatoms[i]].harmonic.rA;
381                 if (ftype == F_G96BONDS)
382                 {
383                     /* Convert squared sqaure fc to harmonic fc */
384                     fc = 2 * fc * re;
385                 }
386                 a1 = ilb.iatoms[i + 1];
387                 a2 = ilb.iatoms[i + 2];
388                 m1 = atom[a1].m;
389                 m2 = atom[a2].m;
390                 if (fc > 0 && m1 > 0 && m2 > 0)
391                 {
392                     period2 = twopi2 * m1 * m2 / ((m1 + m2) * fc);
393                 }
394                 else
395                 {
396                     period2 = GMX_FLOAT_MAX;
397                 }
398                 if (debug)
399                 {
400                     fprintf(debug, "fc %g m1 %g m2 %g period %g\n", fc, m1, m2, std::sqrt(period2));
401                 }
402                 if (period2 < limit2)
403                 {
404                     bFound = false;
405                     for (j = 0; j < ilc.size(); j += 3)
406                     {
407                         if ((ilc.iatoms[j + 1] == a1 && ilc.iatoms[j + 2] == a2)
408                             || (ilc.iatoms[j + 1] == a2 && ilc.iatoms[j + 2] == a1))
409                         {
410                             bFound = true;
411                         }
412                     }
413                     for (j = 0; j < ils.size(); j += 4)
414                     {
415                         if ((a1 == ils.iatoms[j + 1] || a1 == ils.iatoms[j + 2] || a1 == ils.iatoms[j + 3])
416                             && (a2 == ils.iatoms[j + 1] || a2 == ils.iatoms[j + 2]
417                                 || a2 == ils.iatoms[j + 3]))
418                         {
419                             bFound = true;
420                         }
421                     }
422                     if (!bFound && (w_moltype == nullptr || period2 < w_period2))
423                     {
424                         w_moltype = &moltype;
425                         w_a1      = a1;
426                         w_a2      = a2;
427                         w_period2 = period2;
428                     }
429                 }
430             }
431         }
432     }
433
434     if (w_moltype != nullptr)
435     {
436         bWarn = (w_period2 < gmx::square(min_steps_warn * dt));
437         /* A check that would recognize most water models */
438         bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' && w_moltype->atoms.nr <= 5);
439         std::string warningMessage = gmx::formatString(
440                 "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated "
441                 "oscillational period of %.1e ps, which is less than %d times the time step of "
442                 "%.1e ps.\n"
443                 "%s",
444                 *w_moltype->name,
445                 w_a1 + 1,
446                 *w_moltype->atoms.atomname[w_a1],
447                 w_a2 + 1,
448                 *w_moltype->atoms.atomname[w_a2],
449                 std::sqrt(w_period2),
450                 bWarn ? min_steps_warn : min_steps_note,
451                 dt,
452                 bWater ? "Maybe you asked for fexible water."
453                        : "Maybe you forgot to change the constraints mdp option.");
454         if (bWarn)
455         {
456             warning(wi, warningMessage.c_str());
457         }
458         else
459         {
460             warning_note(wi, warningMessage.c_str());
461         }
462     }
463 }
464
465 static void check_vel(gmx_mtop_t* mtop, rvec v[])
466 {
467     for (const AtomProxy atomP : AtomRange(*mtop))
468     {
469         const t_atom& local = atomP.atom();
470         int           i     = atomP.globalAtomNumber();
471         if (local.ptype == eptShell || local.ptype == eptBond || local.ptype == eptVSite)
472         {
473             clear_rvec(v[i]);
474         }
475     }
476 }
477
478 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
479 {
480     int nshells = 0;
481
482     for (const AtomProxy atomP : AtomRange(*mtop))
483     {
484         const t_atom& local = atomP.atom();
485         if (local.ptype == eptShell || local.ptype == eptBond)
486         {
487             nshells++;
488         }
489     }
490     if ((nshells > 0) && (ir->nstcalcenergy != 1))
491     {
492         set_warning_line(wi, "unknown", -1);
493         std::string warningMessage = gmx::formatString(
494                 "There are %d shells, changing nstcalcenergy from %d to 1", nshells, ir->nstcalcenergy);
495         ir->nstcalcenergy = 1;
496         warning(wi, warningMessage.c_str());
497     }
498 }
499
500 /* TODO Decide whether this function can be consolidated with
501  * gmx_mtop_ftype_count */
502 static int nint_ftype(gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
503 {
504     int nint = 0;
505     for (const gmx_molblock_t& molb : mtop->molblock)
506     {
507         nint += molb.nmol * mi[molb.type].interactions[ftype].size();
508     }
509
510     return nint;
511 }
512
513 /* This routine reorders the molecule type array
514  * in the order of use in the molblocks,
515  * unused molecule types are deleted.
516  */
517 static void renumber_moltypes(gmx_mtop_t* sys, std::vector<MoleculeInformation>* molinfo)
518 {
519
520     std::vector<int> order;
521     for (gmx_molblock_t& molblock : sys->molblock)
522     {
523         const auto found = std::find_if(order.begin(), order.end(), [&molblock](const auto& entry) {
524             return molblock.type == entry;
525         });
526         if (found == order.end())
527         {
528             /* This type did not occur yet, add it */
529             order.push_back(molblock.type);
530             molblock.type = order.size() - 1;
531         }
532         else
533         {
534             molblock.type = std::distance(order.begin(), found);
535         }
536     }
537
538     /* We still need to reorder the molinfo structs */
539     std::vector<MoleculeInformation> minew(order.size());
540     int                              index = 0;
541     for (auto& mi : *molinfo)
542     {
543         const auto found = std::find(order.begin(), order.end(), index);
544         if (found != order.end())
545         {
546             int pos    = std::distance(order.begin(), found);
547             minew[pos] = mi;
548         }
549         else
550         {
551             // Need to manually clean up memory ....
552             mi.fullCleanUp();
553         }
554         index++;
555     }
556
557     *molinfo = minew;
558 }
559
560 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t* mtop)
561 {
562     mtop->moltype.resize(mi.size());
563     int pos = 0;
564     for (const auto& mol : mi)
565     {
566         gmx_moltype_t& molt = mtop->moltype[pos];
567         molt.name           = mol.name;
568         molt.atoms          = mol.atoms;
569         /* ilists are copied later */
570         molt.excls = mol.excls;
571         pos++;
572     }
573 }
574
575 static void new_status(const char*                           topfile,
576                        const char*                           topppfile,
577                        const char*                           confin,
578                        t_gromppopts*                         opts,
579                        t_inputrec*                           ir,
580                        gmx_bool                              bZero,
581                        bool                                  bGenVel,
582                        bool                                  bVerbose,
583                        t_state*                              state,
584                        PreprocessingAtomTypes*               atypes,
585                        gmx_mtop_t*                           sys,
586                        std::vector<MoleculeInformation>*     mi,
587                        std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
588                        gmx::ArrayRef<InteractionsOfType>     interactions,
589                        CombinationRule*                      comb,
590                        double*                               reppow,
591                        real*                                 fudgeQQ,
592                        gmx_bool                              bMorse,
593                        warninp*                              wi,
594                        const gmx::MDLogger&                  logger)
595 {
596     std::vector<gmx_molblock_t> molblock;
597     int                         i, nmismatch;
598     bool                        ffParametrizedWithHBondConstraints;
599
600     /* TOPOLOGY processing */
601     sys->name = do_top(bVerbose,
602                        topfile,
603                        topppfile,
604                        opts,
605                        bZero,
606                        &(sys->symtab),
607                        interactions,
608                        comb,
609                        reppow,
610                        fudgeQQ,
611                        atypes,
612                        mi,
613                        intermolecular_interactions,
614                        ir,
615                        &molblock,
616                        &ffParametrizedWithHBondConstraints,
617                        wi,
618                        logger);
619
620     sys->molblock.clear();
621
622     sys->natoms = 0;
623     for (const gmx_molblock_t& molb : molblock)
624     {
625         if (!sys->molblock.empty() && molb.type == sys->molblock.back().type)
626         {
627             /* Merge consecutive blocks with the same molecule type */
628             sys->molblock.back().nmol += molb.nmol;
629         }
630         else if (molb.nmol > 0)
631         {
632             /* Add a new molblock to the topology */
633             sys->molblock.push_back(molb);
634         }
635         sys->natoms += molb.nmol * (*mi)[sys->molblock.back().type].atoms.nr;
636     }
637     if (sys->molblock.empty())
638     {
639         gmx_fatal(FARGS, "No molecules were defined in the system");
640     }
641
642     renumber_moltypes(sys, mi);
643
644     if (bMorse)
645     {
646         convert_harmonics(*mi, atypes);
647     }
648
649     if (ir->eDisre == DistanceRestraintRefinement::None)
650     {
651         i = rm_interactions(F_DISRES, *mi);
652         if (i > 0)
653         {
654             set_warning_line(wi, "unknown", -1);
655             std::string warningMessage =
656                     gmx::formatString("disre = no, removed %d distance restraints", i);
657             warning_note(wi, warningMessage.c_str());
658         }
659     }
660     if (!opts->bOrire)
661     {
662         i = rm_interactions(F_ORIRES, *mi);
663         if (i > 0)
664         {
665             set_warning_line(wi, "unknown", -1);
666             std::string warningMessage =
667                     gmx::formatString("orire = no, removed %d orientation restraints", i);
668             warning_note(wi, warningMessage.c_str());
669         }
670     }
671
672     /* Copy structures from msys to sys */
673     molinfo2mtop(*mi, sys);
674
675     sys->finalize();
676
677     /* Discourage using the, previously common, setup of applying constraints
678      * to all bonds with force fields that have been parametrized with H-bond
679      * constraints only. Do not print note with large timesteps or vsites.
680      */
681     if (opts->nshake == eshALLBONDS && ffParametrizedWithHBondConstraints && ir->delta_t < 0.0026
682         && gmx_mtop_ftype_count(sys, F_VSITE3FD) == 0)
683     {
684         set_warning_line(wi, "unknown", -1);
685         warning_note(wi,
686                      "You are using constraints on all bonds, whereas the forcefield "
687                      "has been parametrized only with constraints involving hydrogen atoms. "
688                      "We suggest using constraints = h-bonds instead, this will also improve "
689                      "performance.");
690     }
691
692     /* COORDINATE file processing */
693     if (bVerbose)
694     {
695         GMX_LOG(logger.info).asParagraph().appendTextFormatted("processing coordinates...");
696     }
697
698     t_topology* conftop;
699     rvec*       x = nullptr;
700     rvec*       v = nullptr;
701     snew(conftop, 1);
702     read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
703     state->natoms = conftop->atoms.nr;
704     if (state->natoms != sys->natoms)
705     {
706         gmx_fatal(FARGS,
707                   "number of coordinates in coordinate file (%s, %d)\n"
708                   "             does not match topology (%s, %d)",
709                   confin,
710                   state->natoms,
711                   topfile,
712                   sys->natoms);
713     }
714     /* It would be nice to get rid of the copies below, but we don't know
715      * a priori if the number of atoms in confin matches what we expect.
716      */
717     state->flags |= enumValueToBitMask(StateEntry::X);
718     if (v != nullptr)
719     {
720         state->flags |= enumValueToBitMask(StateEntry::V);
721     }
722     state_change_natoms(state, state->natoms);
723     std::copy(x, x + state->natoms, state->x.data());
724     sfree(x);
725     if (v != nullptr)
726     {
727         std::copy(v, v + state->natoms, state->v.data());
728         sfree(v);
729     }
730     /* This call fixes the box shape for runs with pressure scaling */
731     set_box_rel(ir, state);
732
733     nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms, logger);
734     done_top(conftop);
735     sfree(conftop);
736
737     if (nmismatch)
738     {
739         std::string warningMessage = gmx::formatString(
740                 "%d non-matching atom name%s\n"
741                 "atom names from %s will be used\n"
742                 "atom names from %s will be ignored\n",
743                 nmismatch,
744                 (nmismatch == 1) ? "" : "s",
745                 topfile,
746                 confin);
747         warning(wi, warningMessage.c_str());
748     }
749
750     /* Do more checks, mostly related to constraints */
751     if (bVerbose)
752     {
753         GMX_LOG(logger.info)
754                 .asParagraph()
755                 .appendTextFormatted("double-checking input for internal consistency...");
756     }
757     {
758         bool bHasNormalConstraints =
759                 0 < (nint_ftype(sys, *mi, F_CONSTR) + nint_ftype(sys, *mi, F_CONSTRNC));
760         bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
761         double_check(ir, state->box, bHasNormalConstraints, bHasAnyConstraints, wi);
762     }
763
764     if (bGenVel)
765     {
766         real* mass;
767
768         snew(mass, state->natoms);
769         for (const AtomProxy atomP : AtomRange(*sys))
770         {
771             const t_atom& local = atomP.atom();
772             int           i     = atomP.globalAtomNumber();
773             mass[i]             = local.m;
774         }
775
776         if (opts->seed == -1)
777         {
778             opts->seed = static_cast<int>(gmx::makeRandomSeed());
779             GMX_LOG(logger.info).asParagraph().appendTextFormatted("Setting gen_seed to %d", opts->seed);
780         }
781         state->flags |= enumValueToBitMask(StateEntry::V);
782         maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array(), logger);
783
784         stop_cm(logger, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
785         sfree(mass);
786     }
787 }
788
789 static void copy_state(const char* slog, t_trxframe* fr, bool bReadVel, t_state* state, double* use_time)
790 {
791     if (fr->not_ok & FRAME_NOT_OK)
792     {
793         gmx_fatal(FARGS, "Can not start from an incomplete frame");
794     }
795     if (!fr->bX)
796     {
797         gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s", slog);
798     }
799
800     std::copy(fr->x, fr->x + state->natoms, state->x.data());
801     if (bReadVel)
802     {
803         if (!fr->bV)
804         {
805             gmx_incons("Trajecory frame unexpectedly does not contain velocities");
806         }
807         std::copy(fr->v, fr->v + state->natoms, state->v.data());
808     }
809     if (fr->bBox)
810     {
811         copy_mat(fr->box, state->box);
812     }
813
814     *use_time = fr->time;
815 }
816
817 static void cont_status(const char*             slog,
818                         const char*             ener,
819                         bool                    bNeedVel,
820                         bool                    bGenVel,
821                         real                    fr_time,
822                         t_inputrec*             ir,
823                         t_state*                state,
824                         gmx_mtop_t*             sys,
825                         const gmx_output_env_t* oenv,
826                         const gmx::MDLogger&    logger)
827 /* If fr_time == -1 read the last frame available which is complete */
828 {
829     bool         bReadVel;
830     t_trxframe   fr;
831     t_trxstatus* fp;
832     double       use_time;
833
834     bReadVel = (bNeedVel && !bGenVel);
835
836     GMX_LOG(logger.info)
837             .asParagraph()
838             .appendTextFormatted("Reading Coordinates%s and Box size from old trajectory",
839                                  bReadVel ? ", Velocities" : "");
840     if (fr_time == -1)
841     {
842         GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read whole trajectory");
843     }
844     else
845     {
846         GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read till time %g", fr_time);
847     }
848     if (!bReadVel)
849     {
850         if (bGenVel)
851         {
852             GMX_LOG(logger.info)
853                     .asParagraph()
854                     .appendTextFormatted(
855                             "Velocities generated: "
856                             "ignoring velocities in input trajectory");
857         }
858         read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
859     }
860     else
861     {
862         read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
863
864         if (!fr.bV)
865         {
866             GMX_LOG(logger.warning)
867                     .asParagraph()
868                     .appendTextFormatted(
869                             "WARNING: Did not find a frame with velocities in file %s,\n"
870                             "         all velocities will be set to zero!",
871                             slog);
872             for (auto& vi : makeArrayRef(state->v))
873             {
874                 vi = { 0, 0, 0 };
875             }
876             close_trx(fp);
877             /* Search for a frame without velocities */
878             bReadVel = false;
879             read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
880         }
881     }
882
883     state->natoms = fr.natoms;
884
885     if (sys->natoms != state->natoms)
886     {
887         gmx_fatal(FARGS,
888                   "Number of atoms in Topology "
889                   "is not the same as in Trajectory");
890     }
891     copy_state(slog, &fr, bReadVel, state, &use_time);
892
893     /* Find the appropriate frame */
894     while ((fr_time == -1 || fr.time < fr_time) && read_next_frame(oenv, fp, &fr))
895     {
896         copy_state(slog, &fr, bReadVel, state, &use_time);
897     }
898
899     close_trx(fp);
900
901     /* Set the relative box lengths for preserving the box shape.
902      * Note that this call can lead to differences in the last bit
903      * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
904      */
905     set_box_rel(ir, state);
906
907     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Using frame at t = %g ps", use_time);
908     GMX_LOG(logger.info).asParagraph().appendTextFormatted("Starting time for run is %g ps", ir->init_t);
909
910     if ((ir->epc != PressureCoupling::No || ir->etc == TemperatureCoupling::NoseHoover) && ener)
911     {
912         get_enx_state(ener, use_time, sys->groups, ir, state);
913         preserve_box_shape(ir, state->box_rel, state->boxv);
914     }
915 }
916
917 static void read_posres(gmx_mtop_t*                              mtop,
918                         gmx::ArrayRef<const MoleculeInformation> molinfo,
919                         gmx_bool                                 bTopB,
920                         const char*                              fn,
921                         RefCoordScaling                          rc_scaling,
922                         PbcType                                  pbcType,
923                         rvec                                     com,
924                         warninp*                                 wi,
925                         const gmx::MDLogger&                     logger)
926 {
927     gmx_bool*   hadAtom;
928     rvec *      x, *v;
929     dvec        sum;
930     double      totmass;
931     t_topology* top;
932     matrix      box, invbox;
933     int         natoms, npbcdim = 0;
934     int         a, nat_molb;
935     t_atom*     atom;
936
937     snew(top, 1);
938     read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
939     natoms = top->atoms.nr;
940     done_top(top);
941     sfree(top);
942     if (natoms != mtop->natoms)
943     {
944         std::string warningMessage = gmx::formatString(
945                 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
946                 "(%d). Will assume that the first %d atoms in the topology and %s match.",
947                 fn,
948                 natoms,
949                 mtop->natoms,
950                 std::min(mtop->natoms, natoms),
951                 fn);
952         warning(wi, warningMessage.c_str());
953     }
954
955     npbcdim = numPbcDimensions(pbcType);
956     GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
957     clear_rvec(com);
958     if (rc_scaling != RefCoordScaling::No)
959     {
960         copy_mat(box, invbox);
961         for (int j = npbcdim; j < DIM; j++)
962         {
963             clear_rvec(invbox[j]);
964             invbox[j][j] = 1;
965         }
966         gmx::invertBoxMatrix(invbox, invbox);
967     }
968
969     /* Copy the reference coordinates to mtop */
970     clear_dvec(sum);
971     totmass = 0;
972     a       = 0;
973     snew(hadAtom, natoms);
974     for (gmx_molblock_t& molb : mtop->molblock)
975     {
976         nat_molb                       = molb.nmol * mtop->moltype[molb.type].atoms.nr;
977         const InteractionsOfType* pr   = &(molinfo[molb.type].interactions[F_POSRES]);
978         const InteractionsOfType* prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
979         if (pr->size() > 0 || prfb->size() > 0)
980         {
981             atom = mtop->moltype[molb.type].atoms.atom;
982             for (const auto& restraint : pr->interactionTypes)
983             {
984                 int ai = restraint.ai();
985                 if (ai >= natoms)
986                 {
987                     gmx_fatal(FARGS,
988                               "Position restraint atom index (%d) in moltype '%s' is larger than "
989                               "number of atoms in %s (%d).\n",
990                               ai + 1,
991                               *molinfo[molb.type].name,
992                               fn,
993                               natoms);
994                 }
995                 hadAtom[ai] = TRUE;
996                 if (rc_scaling == RefCoordScaling::Com)
997                 {
998                     /* Determine the center of mass of the posres reference coordinates */
999                     for (int j = 0; j < npbcdim; j++)
1000                     {
1001                         sum[j] += atom[ai].m * x[a + ai][j];
1002                     }
1003                     totmass += atom[ai].m;
1004                 }
1005             }
1006             /* Same for flat-bottomed posres, but do not count an atom twice for COM */
1007             for (const auto& restraint : prfb->interactionTypes)
1008             {
1009                 int ai = restraint.ai();
1010                 if (ai >= natoms)
1011                 {
1012                     gmx_fatal(FARGS,
1013                               "Position restraint atom index (%d) in moltype '%s' is larger than "
1014                               "number of atoms in %s (%d).\n",
1015                               ai + 1,
1016                               *molinfo[molb.type].name,
1017                               fn,
1018                               natoms);
1019                 }
1020                 if (rc_scaling == RefCoordScaling::Com && !hadAtom[ai])
1021                 {
1022                     /* Determine the center of mass of the posres reference coordinates */
1023                     for (int j = 0; j < npbcdim; j++)
1024                     {
1025                         sum[j] += atom[ai].m * x[a + ai][j];
1026                     }
1027                     totmass += atom[ai].m;
1028                 }
1029             }
1030             if (!bTopB)
1031             {
1032                 molb.posres_xA.resize(nat_molb);
1033                 for (int i = 0; i < nat_molb; i++)
1034                 {
1035                     copy_rvec(x[a + i], molb.posres_xA[i]);
1036                 }
1037             }
1038             else
1039             {
1040                 molb.posres_xB.resize(nat_molb);
1041                 for (int i = 0; i < nat_molb; i++)
1042                 {
1043                     copy_rvec(x[a + i], molb.posres_xB[i]);
1044                 }
1045             }
1046         }
1047         a += nat_molb;
1048     }
1049     if (rc_scaling == RefCoordScaling::Com)
1050     {
1051         if (totmass == 0)
1052         {
1053             gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
1054         }
1055         for (int j = 0; j < npbcdim; j++)
1056         {
1057             com[j] = sum[j] / totmass;
1058         }
1059         GMX_LOG(logger.info)
1060                 .asParagraph()
1061                 .appendTextFormatted(
1062                         "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f",
1063                         com[XX],
1064                         com[YY],
1065                         com[ZZ]);
1066     }
1067
1068     if (rc_scaling != RefCoordScaling::No)
1069     {
1070         GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
1071
1072         for (gmx_molblock_t& molb : mtop->molblock)
1073         {
1074             nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
1075             if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1076             {
1077                 std::vector<gmx::RVec>& xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1078                 for (int i = 0; i < nat_molb; i++)
1079                 {
1080                     for (int j = 0; j < npbcdim; j++)
1081                     {
1082                         if (rc_scaling == RefCoordScaling::All)
1083                         {
1084                             /* Convert from Cartesian to crystal coordinates */
1085                             xp[i][j] *= invbox[j][j];
1086                             for (int k = j + 1; k < npbcdim; k++)
1087                             {
1088                                 xp[i][j] += invbox[k][j] * xp[i][k];
1089                             }
1090                         }
1091                         else if (rc_scaling == RefCoordScaling::Com)
1092                         {
1093                             /* Subtract the center of mass */
1094                             xp[i][j] -= com[j];
1095                         }
1096                     }
1097                 }
1098             }
1099         }
1100
1101         if (rc_scaling == RefCoordScaling::Com)
1102         {
1103             /* Convert the COM from Cartesian to crystal coordinates */
1104             for (int j = 0; j < npbcdim; j++)
1105             {
1106                 com[j] *= invbox[j][j];
1107                 for (int k = j + 1; k < npbcdim; k++)
1108                 {
1109                     com[j] += invbox[k][j] * com[k];
1110                 }
1111             }
1112         }
1113     }
1114
1115     sfree(x);
1116     sfree(v);
1117     sfree(hadAtom);
1118 }
1119
1120 static void gen_posres(gmx_mtop_t*                              mtop,
1121                        gmx::ArrayRef<const MoleculeInformation> mi,
1122                        const char*                              fnA,
1123                        const char*                              fnB,
1124                        RefCoordScaling                          rc_scaling,
1125                        PbcType                                  pbcType,
1126                        rvec                                     com,
1127                        rvec                                     comB,
1128                        warninp*                                 wi,
1129                        const gmx::MDLogger&                     logger)
1130 {
1131     read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi, logger);
1132     /* It is safer to simply read the b-state posres rather than trying
1133      * to be smart and copy the positions.
1134      */
1135     read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi, logger);
1136 }
1137
1138 static void set_wall_atomtype(PreprocessingAtomTypes* at,
1139                               t_gromppopts*           opts,
1140                               t_inputrec*             ir,
1141                               warninp*                wi,
1142                               const gmx::MDLogger&    logger)
1143 {
1144     int i;
1145
1146     if (ir->nwall > 0)
1147     {
1148         GMX_LOG(logger.info).asParagraph().appendTextFormatted("Searching the wall atom type(s)");
1149     }
1150     for (i = 0; i < ir->nwall; i++)
1151     {
1152         auto atomType = at->atomTypeFromName(opts->wall_atomtype[i]);
1153         if (!atomType.has_value())
1154         {
1155             std::string warningMessage = gmx::formatString(
1156                     "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1157             warning_error(wi, warningMessage.c_str());
1158         }
1159         else
1160         {
1161             ir->wall_atomtype[i] = *atomType;
1162         }
1163     }
1164 }
1165
1166 static int nrdf_internal(const t_atoms* atoms)
1167 {
1168     int i, nmass, nrdf;
1169
1170     nmass = 0;
1171     for (i = 0; i < atoms->nr; i++)
1172     {
1173         /* Vsite ptype might not be set here yet, so also check the mass */
1174         if ((atoms->atom[i].ptype == eptAtom || atoms->atom[i].ptype == eptNucleus) && atoms->atom[i].m > 0)
1175         {
1176             nmass++;
1177         }
1178     }
1179     switch (nmass)
1180     {
1181         case 0: // Fall through intended
1182         case 1: nrdf = 0; break;
1183         case 2: nrdf = 1; break;
1184         default: nrdf = nmass * 3 - 6; break;
1185     }
1186
1187     return nrdf;
1188 }
1189
1190 static void spline1d(double dx, const double* y, int n, double* u, double* y2)
1191 {
1192     int    i;
1193     double p, q;
1194
1195     y2[0] = 0.0;
1196     u[0]  = 0.0;
1197
1198     for (i = 1; i < n - 1; i++)
1199     {
1200         p     = 0.5 * y2[i - 1] + 2.0;
1201         y2[i] = -0.5 / p;
1202         q     = (y[i + 1] - 2.0 * y[i] + y[i - 1]) / dx;
1203         u[i]  = (3.0 * q / dx - 0.5 * u[i - 1]) / p;
1204     }
1205
1206     y2[n - 1] = 0.0;
1207
1208     for (i = n - 2; i >= 0; i--)
1209     {
1210         y2[i] = y2[i] * y2[i + 1] + u[i];
1211     }
1212 }
1213
1214
1215 static void
1216 interpolate1d(double xmin, double dx, const double* ya, const double* y2a, double x, double* y, double* y1)
1217 {
1218     int    ix;
1219     double a, b;
1220
1221     ix = static_cast<int>((x - xmin) / dx);
1222
1223     a = (xmin + (ix + 1) * dx - x) / dx;
1224     b = (x - xmin - ix * dx) / dx;
1225
1226     *y = a * ya[ix] + b * ya[ix + 1]
1227          + ((a * a * a - a) * y2a[ix] + (b * b * b - b) * y2a[ix + 1]) * (dx * dx) / 6.0;
1228     *y1 = (ya[ix + 1] - ya[ix]) / dx - (3.0 * a * a - 1.0) / 6.0 * dx * y2a[ix]
1229           + (3.0 * b * b - 1.0) / 6.0 * dx * y2a[ix + 1];
1230 }
1231
1232
1233 static void setup_cmap(int grid_spacing, int nc, gmx::ArrayRef<const real> grid, gmx_cmap_t* cmap_grid)
1234 {
1235     int    i, j, k, ii, jj, kk, idx;
1236     int    offset;
1237     double dx, xmin, v, v1, v2, v12;
1238     double phi, psi;
1239
1240     std::vector<double> tmp_u(2 * grid_spacing, 0.0);
1241     std::vector<double> tmp_u2(2 * grid_spacing, 0.0);
1242     std::vector<double> tmp_yy(2 * grid_spacing, 0.0);
1243     std::vector<double> tmp_y1(2 * grid_spacing, 0.0);
1244     std::vector<double> tmp_t2(2 * grid_spacing * 2 * grid_spacing, 0.0);
1245     std::vector<double> tmp_grid(2 * grid_spacing * 2 * grid_spacing, 0.0);
1246
1247     dx   = 360.0 / grid_spacing;
1248     xmin = -180.0 - dx * grid_spacing / 2;
1249
1250     for (kk = 0; kk < nc; kk++)
1251     {
1252         /* Compute an offset depending on which cmap we are using
1253          * Offset will be the map number multiplied with the
1254          * grid_spacing * grid_spacing * 2
1255          */
1256         offset = kk * grid_spacing * grid_spacing * 2;
1257
1258         for (i = 0; i < 2 * grid_spacing; i++)
1259         {
1260             ii = (i + grid_spacing - grid_spacing / 2) % grid_spacing;
1261
1262             for (j = 0; j < 2 * grid_spacing; j++)
1263             {
1264                 jj = (j + grid_spacing - grid_spacing / 2) % grid_spacing;
1265                 tmp_grid[i * grid_spacing * 2 + j] = grid[offset + ii * grid_spacing + jj];
1266             }
1267         }
1268
1269         for (i = 0; i < 2 * grid_spacing; i++)
1270         {
1271             spline1d(dx,
1272                      &(tmp_grid[2 * grid_spacing * i]),
1273                      2 * grid_spacing,
1274                      tmp_u.data(),
1275                      &(tmp_t2[2 * grid_spacing * i]));
1276         }
1277
1278         for (i = grid_spacing / 2; i < grid_spacing + grid_spacing / 2; i++)
1279         {
1280             ii  = i - grid_spacing / 2;
1281             phi = ii * dx - 180.0;
1282
1283             for (j = grid_spacing / 2; j < grid_spacing + grid_spacing / 2; j++)
1284             {
1285                 jj  = j - grid_spacing / 2;
1286                 psi = jj * dx - 180.0;
1287
1288                 for (k = 0; k < 2 * grid_spacing; k++)
1289                 {
1290                     interpolate1d(xmin,
1291                                   dx,
1292                                   &(tmp_grid[2 * grid_spacing * k]),
1293                                   &(tmp_t2[2 * grid_spacing * k]),
1294                                   psi,
1295                                   &tmp_yy[k],
1296                                   &tmp_y1[k]);
1297                 }
1298
1299                 spline1d(dx, tmp_yy.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1300                 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1301                 spline1d(dx, tmp_y1.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1302                 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1303
1304                 idx                                       = ii * grid_spacing + jj;
1305                 cmap_grid->cmapdata[kk].cmap[idx * 4]     = grid[offset + ii * grid_spacing + jj];
1306                 cmap_grid->cmapdata[kk].cmap[idx * 4 + 1] = v1;
1307                 cmap_grid->cmapdata[kk].cmap[idx * 4 + 2] = v2;
1308                 cmap_grid->cmapdata[kk].cmap[idx * 4 + 3] = v12;
1309             }
1310         }
1311     }
1312 }
1313
1314 static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
1315 {
1316     int i, nelem;
1317
1318     cmap_grid->grid_spacing = grid_spacing;
1319     nelem                   = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
1320
1321     cmap_grid->cmapdata.resize(ngrid);
1322
1323     for (i = 0; i < ngrid; i++)
1324     {
1325         cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
1326     }
1327 }
1328
1329
1330 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
1331 {
1332     int count, count_mol;
1333
1334     count = 0;
1335     for (const gmx_molblock_t& molb : mtop->molblock)
1336     {
1337         count_mol                                            = 0;
1338         gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1339
1340         for (int i = 0; i < F_NRE; i++)
1341         {
1342             if (i == F_SETTLE)
1343             {
1344                 count_mol += 3 * interactions[i].size();
1345             }
1346             else if (interaction_function[i].flags & IF_CONSTRAINT)
1347             {
1348                 count_mol += interactions[i].size();
1349             }
1350         }
1351
1352         if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1353         {
1354             std::string warningMessage = gmx::formatString(
1355                     "Molecule type '%s' has %d constraints.\n"
1356                     "For stability and efficiency there should not be more constraints than "
1357                     "internal number of degrees of freedom: %d.\n",
1358                     *mi[molb.type].name,
1359                     count_mol,
1360                     nrdf_internal(&mi[molb.type].atoms));
1361             warning(wi, warningMessage.c_str());
1362         }
1363         count += molb.nmol * count_mol;
1364     }
1365
1366     return count;
1367 }
1368
1369 static real calc_temp(const gmx_mtop_t* mtop, const t_inputrec* ir, rvec* v)
1370 {
1371     double sum_mv2 = 0;
1372     for (const AtomProxy atomP : AtomRange(*mtop))
1373     {
1374         const t_atom& local = atomP.atom();
1375         int           i     = atomP.globalAtomNumber();
1376         sum_mv2 += local.m * norm2(v[i]);
1377     }
1378
1379     double nrdf = 0;
1380     for (int g = 0; g < ir->opts.ngtc; g++)
1381     {
1382         nrdf += ir->opts.nrdf[g];
1383     }
1384
1385     return sum_mv2 / (nrdf * BOLTZ);
1386 }
1387
1388 static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
1389 {
1390     real ref_t;
1391     int  i;
1392     bool bNoCoupl;
1393
1394     ref_t    = 0;
1395     bNoCoupl = false;
1396     for (i = 0; i < ir->opts.ngtc; i++)
1397     {
1398         if (ir->opts.tau_t[i] < 0)
1399         {
1400             bNoCoupl = true;
1401         }
1402         else
1403         {
1404             ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1405         }
1406     }
1407
1408     if (bNoCoupl)
1409     {
1410         std::string warningMessage = gmx::formatString(
1411                 "Some temperature coupling groups do not use temperature coupling. We will assume "
1412                 "their temperature is not more than %.3f K. If their temperature is higher, the "
1413                 "energy error and the Verlet buffer might be underestimated.",
1414                 ref_t);
1415         warning(wi, warningMessage.c_str());
1416     }
1417
1418     return ref_t;
1419 }
1420
1421 /* Checks if there are unbound atoms in moleculetype molt.
1422  * Prints a note for each unbound atoms and a warning if any is present.
1423  */
1424 static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1425 {
1426     const t_atoms* atoms = &molt->atoms;
1427
1428     if (atoms->nr == 1)
1429     {
1430         /* Only one atom, there can't be unbound atoms */
1431         return;
1432     }
1433
1434     std::vector<int> count(atoms->nr, 0);
1435
1436     for (int ftype = 0; ftype < F_NRE; ftype++)
1437     {
1438         if (((interaction_function[ftype].flags & IF_BOND) && NRAL(ftype) == 2 && ftype != F_CONNBONDS)
1439             || (interaction_function[ftype].flags & IF_CONSTRAINT) || ftype == F_SETTLE)
1440         {
1441             const InteractionList& il   = molt->ilist[ftype];
1442             const int              nral = NRAL(ftype);
1443
1444             for (int i = 0; i < il.size(); i += 1 + nral)
1445             {
1446                 for (int j = 0; j < nral; j++)
1447                 {
1448                     count[il.iatoms[i + 1 + j]]++;
1449                 }
1450             }
1451         }
1452     }
1453
1454     int numDanglingAtoms = 0;
1455     for (int a = 0; a < atoms->nr; a++)
1456     {
1457         if (atoms->atom[a].ptype != eptVSite && count[a] == 0)
1458         {
1459             if (bVerbose)
1460             {
1461                 GMX_LOG(logger.warning)
1462                         .asParagraph()
1463                         .appendTextFormatted(
1464                                 "Atom %d '%s' in moleculetype '%s' is not bound by a potential or "
1465                                 "constraint to any other atom in the same moleculetype.",
1466                                 a + 1,
1467                                 *atoms->atomname[a],
1468                                 *molt->name);
1469             }
1470             numDanglingAtoms++;
1471         }
1472     }
1473
1474     if (numDanglingAtoms > 0)
1475     {
1476         std::string warningMessage = gmx::formatString(
1477                 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
1478                 "other atom in the same moleculetype. Although technically this might not cause "
1479                 "issues in a simulation, this often means that the user forgot to add a "
1480                 "bond/potential/constraint or put multiple molecules in the same moleculetype "
1481                 "definition by mistake. Run with -v to get information for each atom.",
1482                 *molt->name,
1483                 numDanglingAtoms);
1484         warning_note(wi, warningMessage.c_str());
1485     }
1486 }
1487
1488 /* Checks all moleculetypes for unbound atoms */
1489 static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1490 {
1491     for (const gmx_moltype_t& molt : mtop->moltype)
1492     {
1493         checkForUnboundAtoms(&molt, bVerbose, wi, logger);
1494     }
1495 }
1496
1497 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1498  *
1499  * The specific decoupled modes this routine check for are angle modes
1500  * where the two bonds are constrained and the atoms a both ends are only
1501  * involved in a single constraint; the mass of the two atoms needs to
1502  * differ by more than \p massFactorThreshold.
1503  */
1504 static bool haveDecoupledModeInMol(const gmx_moltype_t&           molt,
1505                                    gmx::ArrayRef<const t_iparams> iparams,
1506                                    real                           massFactorThreshold)
1507 {
1508     if (molt.ilist[F_CONSTR].empty() && molt.ilist[F_CONSTRNC].empty())
1509     {
1510         return false;
1511     }
1512
1513     const t_atom* atom = molt.atoms.atom;
1514
1515     const auto atomToConstraints =
1516             gmx::make_at2con(molt, iparams, gmx::FlexibleConstraintTreatment::Exclude);
1517
1518     bool haveDecoupledMode = false;
1519     for (int ftype = 0; ftype < F_NRE; ftype++)
1520     {
1521         if (interaction_function[ftype].flags & IF_ATYPE)
1522         {
1523             const int              nral = NRAL(ftype);
1524             const InteractionList& il   = molt.ilist[ftype];
1525             for (int i = 0; i < il.size(); i += 1 + nral)
1526             {
1527                 /* Here we check for the mass difference between the atoms
1528                  * at both ends of the angle, that the atoms at the ends have
1529                  * 1 contraint and the atom in the middle at least 3; we check
1530                  * that the 3 atoms are linked by constraints below.
1531                  * We check for at least three constraints for the middle atom,
1532                  * since with only the two bonds in the angle, we have 3-atom
1533                  * molecule, which has much more thermal exhange in this single
1534                  * angle mode than molecules with more atoms.
1535                  * Note that this check also catches molecules where more atoms
1536                  * are connected to one or more atoms in the angle, but by
1537                  * bond potentials instead of angles. But such cases will not
1538                  * occur in "normal" molecules and it doesn't hurt running
1539                  * those with higher accuracy settings as well.
1540                  */
1541                 int a0 = il.iatoms[1 + i];
1542                 int a1 = il.iatoms[1 + i + 1];
1543                 int a2 = il.iatoms[1 + i + 2];
1544                 if ((atom[a0].m > atom[a2].m * massFactorThreshold || atom[a2].m > atom[a0].m * massFactorThreshold)
1545                     && atomToConstraints[a0].ssize() == 1 && atomToConstraints[a2].ssize() == 1
1546                     && atomToConstraints[a1].ssize() >= 3)
1547                 {
1548                     int constraint0 = atomToConstraints[a0][0];
1549                     int constraint2 = atomToConstraints[a2][0];
1550
1551                     bool foundAtom0 = false;
1552                     bool foundAtom2 = false;
1553                     for (const int constraint : atomToConstraints[a1])
1554                     {
1555                         if (constraint == constraint0)
1556                         {
1557                             foundAtom0 = true;
1558                         }
1559                         if (constraint == constraint2)
1560                         {
1561                             foundAtom2 = true;
1562                         }
1563                     }
1564                     if (foundAtom0 && foundAtom2)
1565                     {
1566                         haveDecoupledMode = true;
1567                     }
1568                 }
1569             }
1570         }
1571     }
1572
1573     return haveDecoupledMode;
1574 }
1575
1576 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1577  *
1578  * When decoupled modes are present and the accuray in insufficient,
1579  * this routine issues a warning if the accuracy is insufficient.
1580  */
1581 static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec* ir, warninp* wi)
1582 {
1583     /* We only have issues with decoupled modes with normal MD.
1584      * With stochastic dynamics equipartitioning is enforced strongly.
1585      */
1586     if (!EI_MD(ir->eI))
1587     {
1588         return;
1589     }
1590
1591     /* When atoms of very different mass are involved in an angle potential
1592      * and both bonds in the angle are constrained, the dynamic modes in such
1593      * angles have very different periods and significant energy exchange
1594      * takes several nanoseconds. Thus even a small amount of error in
1595      * different algorithms can lead to violation of equipartitioning.
1596      * The parameters below are mainly based on an all-atom chloroform model
1597      * with all bonds constrained. Equipartitioning is off by more than 1%
1598      * (need to run 5-10 ns) when the difference in mass between H and Cl
1599      * is more than a factor 13 and the accuracy is less than the thresholds
1600      * given below. This has been verified on some other molecules.
1601      *
1602      * Note that the buffer and shake parameters have unit length and
1603      * energy/time, respectively, so they will "only" work correctly
1604      * for atomistic force fields using MD units.
1605      */
1606     const real massFactorThreshold      = 13.0;
1607     const real bufferToleranceThreshold = 1e-4;
1608     const int  lincsIterationThreshold  = 2;
1609     const int  lincsOrderThreshold      = 4;
1610     const real shakeToleranceThreshold  = 0.005 * ir->delta_t;
1611
1612     bool lincsWithSufficientTolerance =
1613             (ir->eConstrAlg == ConstraintAlgorithm::Lincs
1614              && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1615     bool shakeWithSufficientTolerance = (ir->eConstrAlg == ConstraintAlgorithm::Shake
1616                                          && ir->shake_tol <= 1.1 * shakeToleranceThreshold);
1617     if (ir->cutoff_scheme == CutoffScheme::Verlet && ir->verletbuf_tol <= 1.1 * bufferToleranceThreshold
1618         && (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1619     {
1620         return;
1621     }
1622
1623     bool haveDecoupledMode = false;
1624     for (const gmx_moltype_t& molt : mtop->moltype)
1625     {
1626         if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams, massFactorThreshold))
1627         {
1628             haveDecoupledMode = true;
1629         }
1630     }
1631
1632     if (haveDecoupledMode)
1633     {
1634         std::string message = gmx::formatString(
1635                 "There are atoms at both ends of an angle, connected by constraints "
1636                 "and with masses that differ by more than a factor of %g. This means "
1637                 "that there are likely dynamic modes that are only very weakly coupled.",
1638                 massFactorThreshold);
1639         if (ir->cutoff_scheme == CutoffScheme::Verlet)
1640         {
1641             message += gmx::formatString(
1642                     " To ensure good equipartitioning, you need to either not use "
1643                     "constraints on all bonds (but, if possible, only on bonds involving "
1644                     "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1645                     "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1646                     ">= %d or SHAKE tolerance <= %g",
1647                     enumValueToString(IntegrationAlgorithm::SD1),
1648                     bufferToleranceThreshold,
1649                     lincsIterationThreshold,
1650                     lincsOrderThreshold,
1651                     shakeToleranceThreshold);
1652         }
1653         else
1654         {
1655             message += gmx::formatString(
1656                     " To ensure good equipartitioning, we suggest to switch to the %s "
1657                     "cutoff-scheme, since that allows for better control over the Verlet "
1658                     "buffer size and thus over the energy drift.",
1659                     enumValueToString(CutoffScheme::Verlet));
1660         }
1661         warning(wi, message);
1662     }
1663 }
1664
1665 static void set_verlet_buffer(const gmx_mtop_t*    mtop,
1666                               t_inputrec*          ir,
1667                               real                 buffer_temp,
1668                               matrix               box,
1669                               warninp*             wi,
1670                               const gmx::MDLogger& logger)
1671 {
1672     GMX_LOG(logger.info)
1673             .asParagraph()
1674             .appendTextFormatted(
1675                     "Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K",
1676                     ir->verletbuf_tol,
1677                     buffer_temp);
1678
1679     /* Calculate the buffer size for simple atom vs atoms list */
1680     VerletbufListSetup listSetup1x1;
1681     listSetup1x1.cluster_size_i = 1;
1682     listSetup1x1.cluster_size_j = 1;
1683     const real rlist_1x1        = calcVerletBufferSize(
1684             *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup1x1);
1685
1686     /* Set the pair-list buffer size in ir */
1687     VerletbufListSetup listSetup4x4 = verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1688     ir->rlist                       = calcVerletBufferSize(
1689             *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup4x4);
1690
1691     const int n_nonlin_vsite = gmx::countNonlinearVsites(*mtop);
1692     if (n_nonlin_vsite > 0)
1693     {
1694         std::string warningMessage = gmx::formatString(
1695                 "There are %d non-linear virtual site constructions. Their contribution to the "
1696                 "energy error is approximated. In most cases this does not affect the error "
1697                 "significantly.",
1698                 n_nonlin_vsite);
1699         warning_note(wi, warningMessage);
1700     }
1701
1702     GMX_LOG(logger.info)
1703             .asParagraph()
1704             .appendTextFormatted(
1705                     "Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm",
1706                     1,
1707                     1,
1708                     rlist_1x1,
1709                     rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
1710
1711     GMX_LOG(logger.info)
1712             .asParagraph()
1713             .appendTextFormatted(
1714                     "Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm",
1715                     listSetup4x4.cluster_size_i,
1716                     listSetup4x4.cluster_size_j,
1717                     ir->rlist,
1718                     ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
1719
1720     GMX_LOG(logger.info)
1721             .asParagraph()
1722             .appendTextFormatted(
1723                     "Note that mdrun will redetermine rlist based on the actual pair-list setup");
1724
1725     if (gmx::square(ir->rlist) >= max_cutoff2(ir->pbcType, box))
1726     {
1727         gmx_fatal(FARGS,
1728                   "The pair-list cut-off (%g nm) is longer than half the shortest box vector or "
1729                   "longer than the smallest box diagonal element (%g nm). Increase the box size or "
1730                   "decrease nstlist or increase verlet-buffer-tolerance.",
1731                   ir->rlist,
1732                   std::sqrt(max_cutoff2(ir->pbcType, box)));
1733     }
1734 }
1735
1736 int gmx_grompp(int argc, char* argv[])
1737 {
1738     const char* desc[] = {
1739         "[THISMODULE] (the gromacs preprocessor)",
1740         "reads a molecular topology file, checks the validity of the",
1741         "file, expands the topology from a molecular description to an atomic",
1742         "description. The topology file contains information about",
1743         "molecule types and the number of molecules, the preprocessor",
1744         "copies each molecule as needed. ",
1745         "There is no limitation on the number of molecule types. ",
1746         "Bonds and bond-angles can be converted into constraints, separately",
1747         "for hydrogens and heavy atoms.",
1748         "Then a coordinate file is read and velocities can be generated",
1749         "from a Maxwellian distribution if requested.",
1750         "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1751         "(eg. number of MD steps, time step, cut-off).",
1752         "Eventually a binary file is produced that can serve as the sole input",
1753         "file for the MD program.[PAR]",
1754
1755         "[THISMODULE] uses the atom names from the topology file. The atom names",
1756         "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1757         "warnings when they do not match the atom names in the topology.",
1758         "Note that the atom names are irrelevant for the simulation as",
1759         "only the atom types are used for generating interaction parameters.[PAR]",
1760
1761         "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1762         "etc. The preprocessor supports the following keywords::",
1763         "",
1764         "    #ifdef VARIABLE",
1765         "    #ifndef VARIABLE",
1766         "    #else",
1767         "    #endif",
1768         "    #define VARIABLE",
1769         "    #undef VARIABLE",
1770         "    #include \"filename\"",
1771         "    #include <filename>",
1772         "",
1773         "The functioning of these statements in your topology may be modulated by",
1774         "using the following two flags in your [REF].mdp[ref] file::",
1775         "",
1776         "    define = -DVARIABLE1 -DVARIABLE2",
1777         "    include = -I/home/john/doe",
1778         "",
1779         "For further information a C-programming textbook may help you out.",
1780         "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1781         "topology file written out so that you can verify its contents.[PAR]",
1782
1783         "When using position restraints, a file with restraint coordinates",
1784         "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1785         "for [TT]-c[tt]). For free energy calculations, separate reference",
1786         "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1787         "otherwise they will be equal to those of the A topology.[PAR]",
1788
1789         "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1790         "The last frame with coordinates and velocities will be read,",
1791         "unless the [TT]-time[tt] option is used. Only if this information",
1792         "is absent will the coordinates in the [TT]-c[tt] file be used.",
1793         "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1794         "in your [REF].mdp[ref] file. An energy file can be supplied with",
1795         "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1796         "variables.[PAR]",
1797
1798         "[THISMODULE] can be used to restart simulations (preserving",
1799         "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1800         "However, for simply changing the number of run steps to extend",
1801         "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1802         "You then supply the old checkpoint file directly to [gmx-mdrun]",
1803         "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1804         "like output frequency, then supplying the checkpoint file to",
1805         "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1806         "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1807         "the ensemble (if possible) still requires passing the checkpoint",
1808         "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1809
1810         "By default, all bonded interactions which have constant energy due to",
1811         "virtual site constructions will be removed. If this constant energy is",
1812         "not zero, this will result in a shift in the total energy. All bonded",
1813         "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1814         "all constraints for distances which will be constant anyway because",
1815         "of virtual site constructions will be removed. If any constraints remain",
1816         "which involve virtual sites, a fatal error will result.[PAR]",
1817
1818         "To verify your run input file, please take note of all warnings",
1819         "on the screen, and correct where necessary. Do also look at the contents",
1820         "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1821         "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1822         "with the [TT]-debug[tt] option which will give you more information",
1823         "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1824         "can see the contents of the run input file with the [gmx-dump]",
1825         "program. [gmx-check] can be used to compare the contents of two",
1826         "run input files.[PAR]",
1827
1828         "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1829         "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1830         "harmless, but usually they are not. The user is advised to carefully",
1831         "interpret the output messages before attempting to bypass them with",
1832         "this option."
1833     };
1834     std::vector<MoleculeInformation>     mi;
1835     std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1836     int                                  nvsite;
1837     CombinationRule                      comb;
1838     real                                 fudgeQQ;
1839     double                               reppow;
1840     const char*                          mdparin;
1841     bool                                 bNeedVel, bGenVel;
1842     gmx_output_env_t*                    oenv;
1843     gmx_bool                             bVerbose = FALSE;
1844     warninp*                             wi;
1845
1846     t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
1847                        { efMDP, "-po", "mdout", ffWRITE },
1848                        { efSTX, "-c", nullptr, ffREAD },
1849                        { efSTX, "-r", "restraint", ffOPTRD },
1850                        { efSTX, "-rb", "restraint", ffOPTRD },
1851                        { efNDX, nullptr, nullptr, ffOPTRD },
1852                        { efTOP, nullptr, nullptr, ffREAD },
1853                        { efTOP, "-pp", "processed", ffOPTWR },
1854                        { efTPR, "-o", nullptr, ffWRITE },
1855                        { efTRN, "-t", nullptr, ffOPTRD },
1856                        { efEDR, "-e", nullptr, ffOPTRD },
1857                        /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1858                        { efGRO, "-imd", "imdgroup", ffOPTWR },
1859                        { efTRN, "-ref", "rotref", ffOPTRW | ffALLOW_MISSING } };
1860 #define NFILE asize(fnm)
1861
1862     /* Command line options */
1863     gmx_bool bRenum   = TRUE;
1864     gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1865     int      i, maxwarn             = 0;
1866     real     fr_time = -1;
1867     t_pargs  pa[]    = {
1868         { "-v", FALSE, etBOOL, { &bVerbose }, "Be loud and noisy" },
1869         { "-time", FALSE, etREAL, { &fr_time }, "Take frame at or first after this time." },
1870         { "-rmvsbds",
1871           FALSE,
1872           etBOOL,
1873           { &bRmVSBds },
1874           "Remove constant bonded interactions with virtual sites" },
1875         { "-maxwarn",
1876           FALSE,
1877           etINT,
1878           { &maxwarn },
1879           "Number of allowed warnings during input processing. Not for normal use and may "
1880           "generate unstable systems" },
1881         { "-zero",
1882           FALSE,
1883           etBOOL,
1884           { &bZero },
1885           "Set parameters for bonded interactions without defaults to zero instead of "
1886           "generating an error" },
1887         { "-renum",
1888           FALSE,
1889           etBOOL,
1890           { &bRenum },
1891           "Renumber atomtypes and minimize number of atomtypes" }
1892     };
1893
1894     /* Parse the command line */
1895     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
1896     {
1897         return 0;
1898     }
1899
1900     /* Initiate some variables */
1901     gmx::MDModules mdModules;
1902     t_inputrec     irInstance;
1903     t_inputrec*    ir = &irInstance;
1904     t_gromppopts   optsInstance;
1905     t_gromppopts*  opts = &optsInstance;
1906     snew(opts->include, STRLEN);
1907     snew(opts->define, STRLEN);
1908
1909     gmx::LoggerBuilder builder;
1910     builder.addTargetStream(gmx::MDLogger::LogLevel::Info, &gmx::TextOutputFile::standardOutput());
1911     builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
1912     gmx::LoggerOwner    logOwner(builder.build());
1913     const gmx::MDLogger logger(logOwner.logger());
1914
1915
1916     wi = init_warning(TRUE, maxwarn);
1917
1918     /* PARAMETER file processing */
1919     mdparin = opt2fn("-f", NFILE, fnm);
1920     set_warning_line(wi, mdparin, -1);
1921     try
1922     {
1923         get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1924     }
1925     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1926
1927     // Now that the MdModules have their options assigned from get_ir, subscribe
1928     // to eventual notifications during pre-processing their data
1929     mdModules.subscribeToPreProcessingNotifications();
1930
1931
1932     if (bVerbose)
1933     {
1934         GMX_LOG(logger.info)
1935                 .asParagraph()
1936                 .appendTextFormatted("checking input for internal consistency...");
1937     }
1938     check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
1939
1940     if (ir->ld_seed == -1)
1941     {
1942         ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1943         GMX_LOG(logger.info)
1944                 .asParagraph()
1945                 .appendTextFormatted("Setting the LD random seed to %" PRId64 "", ir->ld_seed);
1946     }
1947
1948     if (ir->expandedvals->lmc_seed == -1)
1949     {
1950         ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1951         GMX_LOG(logger.info)
1952                 .asParagraph()
1953                 .appendTextFormatted("Setting the lambda MC random seed to %d", ir->expandedvals->lmc_seed);
1954     }
1955
1956     bNeedVel = EI_STATE_VELOCITY(ir->eI);
1957     bGenVel  = (bNeedVel && opts->bGenVel);
1958     if (bGenVel && ir->bContinuation)
1959     {
1960         std::string warningMessage = gmx::formatString(
1961                 "Generating velocities is inconsistent with attempting "
1962                 "to continue a previous run. Choose only one of "
1963                 "gen-vel = yes and continuation = yes.");
1964         warning_error(wi, warningMessage);
1965     }
1966
1967     std::array<InteractionsOfType, F_NRE> interactions;
1968     gmx_mtop_t                            sys;
1969     PreprocessingAtomTypes                atypes;
1970     if (debug)
1971     {
1972         pr_symtab(debug, 0, "Just opened", &sys.symtab);
1973     }
1974
1975     const char* fn = ftp2fn(efTOP, NFILE, fnm);
1976     if (!gmx_fexist(fn))
1977     {
1978         gmx_fatal(FARGS, "%s does not exist", fn);
1979     }
1980
1981     t_state state;
1982     new_status(fn,
1983                opt2fn_null("-pp", NFILE, fnm),
1984                opt2fn("-c", NFILE, fnm),
1985                opts,
1986                ir,
1987                bZero,
1988                bGenVel,
1989                bVerbose,
1990                &state,
1991                &atypes,
1992                &sys,
1993                &mi,
1994                &intermolecular_interactions,
1995                interactions,
1996                &comb,
1997                &reppow,
1998                &fudgeQQ,
1999                opts->bMorse,
2000                wi,
2001                logger);
2002
2003     if (debug)
2004     {
2005         pr_symtab(debug, 0, "After new_status", &sys.symtab);
2006     }
2007
2008     nvsite = 0;
2009     /* set parameters for virtual site construction (not for vsiten) */
2010     for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2011     {
2012         nvsite += set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions, logger);
2013     }
2014     /* now throw away all obsolete bonds, angles and dihedrals: */
2015     /* note: constraints are ALWAYS removed */
2016     if (nvsite)
2017     {
2018         for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2019         {
2020             clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds, logger);
2021         }
2022     }
2023
2024     if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == ConstraintAlgorithm::Shake))
2025     {
2026         if (ir->eI == IntegrationAlgorithm::CG || ir->eI == IntegrationAlgorithm::LBFGS)
2027         {
2028             std::string warningMessage =
2029                     gmx::formatString("Can not do %s with %s, use %s",
2030                                       enumValueToString(ir->eI),
2031                                       enumValueToString(ConstraintAlgorithm::Shake),
2032                                       enumValueToString(ConstraintAlgorithm::Lincs));
2033             warning_error(wi, warningMessage);
2034         }
2035         if (ir->bPeriodicMols)
2036         {
2037             std::string warningMessage =
2038                     gmx::formatString("Can not do periodic molecules with %s, use %s",
2039                                       enumValueToString(ConstraintAlgorithm::Shake),
2040                                       enumValueToString(ConstraintAlgorithm::Lincs));
2041             warning_error(wi, warningMessage);
2042         }
2043     }
2044
2045     if (EI_SD(ir->eI) && ir->etc != TemperatureCoupling::No)
2046     {
2047         warning_note(wi, "Temperature coupling is ignored with SD integrators.");
2048     }
2049
2050     /* Check for errors in the input now, since they might cause problems
2051      * during processing further down.
2052      */
2053     check_warning_error(wi, FARGS);
2054
2055     if (nint_ftype(&sys, mi, F_POSRES) > 0 || nint_ftype(&sys, mi, F_FBPOSRES) > 0)
2056     {
2057         if (ir->epc == PressureCoupling::ParrinelloRahman || ir->epc == PressureCoupling::Mttk)
2058         {
2059             std::string warningMessage = gmx::formatString(
2060                     "You are combining position restraints with %s pressure coupling, which can "
2061                     "lead to instabilities. If you really want to combine position restraints with "
2062                     "pressure coupling, we suggest to use %s pressure coupling instead.",
2063                     enumValueToString(ir->epc),
2064                     enumValueToString(PressureCoupling::Berendsen));
2065             warning_note(wi, warningMessage);
2066         }
2067
2068         const char* fn = opt2fn("-r", NFILE, fnm);
2069         const char* fnB;
2070
2071         if (!gmx_fexist(fn))
2072         {
2073             gmx_fatal(FARGS,
2074                       "Cannot find position restraint file %s (option -r).\n"
2075                       "From GROMACS-2018, you need to specify the position restraint "
2076                       "coordinate files explicitly to avoid mistakes, although you can "
2077                       "still use the same file as you specify for the -c option.",
2078                       fn);
2079         }
2080
2081         if (opt2bSet("-rb", NFILE, fnm))
2082         {
2083             fnB = opt2fn("-rb", NFILE, fnm);
2084             if (!gmx_fexist(fnB))
2085             {
2086                 gmx_fatal(FARGS,
2087                           "Cannot find B-state position restraint file %s (option -rb).\n"
2088                           "From GROMACS-2018, you need to specify the position restraint "
2089                           "coordinate files explicitly to avoid mistakes, although you can "
2090                           "still use the same file as you specify for the -c option.",
2091                           fn);
2092             }
2093         }
2094         else
2095         {
2096             fnB = fn;
2097         }
2098
2099         if (bVerbose)
2100         {
2101             std::string message = gmx::formatString("Reading position restraint coords from %s", fn);
2102             if (strcmp(fn, fnB) != 0)
2103             {
2104                 message += gmx::formatString(" and %s", fnB);
2105             }
2106             GMX_LOG(logger.info).asParagraph().appendText(message);
2107         }
2108         gen_posres(&sys, mi, fn, fnB, ir->refcoord_scaling, ir->pbcType, ir->posres_com, ir->posres_comB, wi, logger);
2109     }
2110
2111     /* If we are using CMAP, setup the pre-interpolation grid */
2112     if (interactions[F_CMAP].ncmap() > 0)
2113     {
2114         init_cmap_grid(&sys.ffparams.cmap_grid,
2115                        interactions[F_CMAP].cmapAngles,
2116                        interactions[F_CMAP].cmakeGridSpacing);
2117         setup_cmap(interactions[F_CMAP].cmakeGridSpacing,
2118                    interactions[F_CMAP].cmapAngles,
2119                    interactions[F_CMAP].cmap,
2120                    &sys.ffparams.cmap_grid);
2121     }
2122
2123     set_wall_atomtype(&atypes, opts, ir, wi, logger);
2124     if (bRenum)
2125     {
2126         atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
2127     }
2128
2129     if (ir->implicit_solvent)
2130     {
2131         gmx_fatal(FARGS, "Implicit solvation is no longer supported");
2132     }
2133
2134     /* PELA: Copy the atomtype data to the topology atomtype list */
2135     atypes.copyTot_atomtypes(&(sys.atomtypes));
2136
2137     if (debug)
2138     {
2139         pr_symtab(debug, 0, "After atype.renumberTypes", &sys.symtab);
2140     }
2141
2142     if (bVerbose)
2143     {
2144         GMX_LOG(logger.info).asParagraph().appendTextFormatted("converting bonded parameters...");
2145     }
2146
2147     const int ntype = atypes.size();
2148     convertInteractionsOfType(
2149             ntype, interactions, mi, intermolecular_interactions.get(), comb, reppow, fudgeQQ, &sys);
2150
2151     if (debug)
2152     {
2153         pr_symtab(debug, 0, "After converInteractionsOfType", &sys.symtab);
2154     }
2155
2156     /* set ptype to VSite for virtual sites */
2157     for (gmx_moltype_t& moltype : sys.moltype)
2158     {
2159         set_vsites_ptype(FALSE, &moltype, logger);
2160     }
2161     if (debug)
2162     {
2163         pr_symtab(debug, 0, "After virtual sites", &sys.symtab);
2164     }
2165     /* Check velocity for virtual sites and shells */
2166     if (bGenVel)
2167     {
2168         check_vel(&sys, state.v.rvec_array());
2169     }
2170
2171     /* check for shells and inpurecs */
2172     check_shells_inputrec(&sys, ir, wi);
2173
2174     /* check masses */
2175     check_mol(&sys, wi);
2176
2177     checkForUnboundAtoms(&sys, bVerbose, wi, logger);
2178
2179     if (EI_DYNAMICS(ir->eI) && ir->eI != IntegrationAlgorithm::BD)
2180     {
2181         check_bonds_timestep(&sys, ir->delta_t, wi);
2182     }
2183
2184     checkDecoupledModeAccuracy(&sys, ir, wi);
2185
2186     if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2187     {
2188         warning_note(
2189                 wi,
2190                 "Zero-step energy minimization will alter the coordinates before calculating the "
2191                 "energy. If you just want the energy of a single point, try zero-step MD (with "
2192                 "unconstrained_start = yes). To do multiple single-point energy evaluations of "
2193                 "different configurations of the same topology, use mdrun -rerun.");
2194     }
2195
2196     check_warning_error(wi, FARGS);
2197
2198     if (bVerbose)
2199     {
2200         GMX_LOG(logger.info).asParagraph().appendTextFormatted("initialising group options...");
2201     }
2202     do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifier(), ir, wi);
2203
2204     if (ir->cutoff_scheme == CutoffScheme::Verlet && ir->verletbuf_tol > 0)
2205     {
2206         if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2207         {
2208             real buffer_temp;
2209
2210             if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No)
2211             {
2212                 if (bGenVel)
2213                 {
2214                     buffer_temp = opts->tempi;
2215                 }
2216                 else
2217                 {
2218                     buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2219                 }
2220                 if (buffer_temp > 0)
2221                 {
2222                     std::string warningMessage = gmx::formatString(
2223                             "NVE simulation: will use the initial temperature of %.3f K for "
2224                             "determining the Verlet buffer size",
2225                             buffer_temp);
2226                     warning_note(wi, warningMessage);
2227                 }
2228                 else
2229                 {
2230                     std::string warningMessage = gmx::formatString(
2231                             "NVE simulation with an initial temperature of zero: will use a Verlet "
2232                             "buffer of %d%%. Check your energy drift!",
2233                             gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
2234                     warning_note(wi, warningMessage);
2235                 }
2236             }
2237             else
2238             {
2239                 buffer_temp = get_max_reference_temp(ir, wi);
2240             }
2241
2242             if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No && buffer_temp == 0)
2243             {
2244                 /* NVE with initial T=0: we add a fixed ratio to rlist.
2245                  * Since we don't actually use verletbuf_tol, we set it to -1
2246                  * so it can't be misused later.
2247                  */
2248                 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2249                 ir->verletbuf_tol = -1;
2250             }
2251             else
2252             {
2253                 /* We warn for NVE simulations with a drift tolerance that
2254                  * might result in a 1(.1)% drift over the total run-time.
2255                  * Note that we can't warn when nsteps=0, since we don't
2256                  * know how many steps the user intends to run.
2257                  */
2258                 if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No && ir->nstlist > 1 && ir->nsteps > 0)
2259                 {
2260                     const real driftTolerance = 0.01;
2261                     /* We use 2 DOF per atom = 2kT pot+kin energy,
2262                      * to be on the safe side with constraints.
2263                      */
2264                     const real totalEnergyDriftPerAtomPerPicosecond =
2265                             2 * BOLTZ * buffer_temp / (ir->nsteps * ir->delta_t);
2266
2267                     if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
2268                     {
2269                         std::string warningMessage = gmx::formatString(
2270                                 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
2271                                 "NVE simulation of length %g ps, which can give a final drift of "
2272                                 "%d%%. For conserving energy to %d%% when using constraints, you "
2273                                 "might need to set verlet-buffer-tolerance to %.1e.",
2274                                 ir->verletbuf_tol,
2275                                 ir->nsteps * ir->delta_t,
2276                                 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
2277                                 gmx::roundToInt(100 * driftTolerance),
2278                                 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
2279                         warning_note(wi, warningMessage);
2280                     }
2281                 }
2282
2283                 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi, logger);
2284             }
2285         }
2286     }
2287
2288     /* Init the temperature coupling state */
2289     init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2290
2291     if (debug)
2292     {
2293         pr_symtab(debug, 0, "After index", &sys.symtab);
2294     }
2295
2296     triple_check(mdparin, ir, &sys, wi);
2297     close_symtab(&sys.symtab);
2298     if (debug)
2299     {
2300         pr_symtab(debug, 0, "After close", &sys.symtab);
2301     }
2302
2303     if (ir->eI == IntegrationAlgorithm::Mimic)
2304     {
2305         generate_qmexcl(&sys, ir, logger);
2306     }
2307
2308     if (ftp2bSet(efTRN, NFILE, fnm))
2309     {
2310         if (bVerbose)
2311         {
2312             GMX_LOG(logger.info)
2313                     .asParagraph()
2314                     .appendTextFormatted("getting data from old trajectory ...");
2315         }
2316         cont_status(ftp2fn(efTRN, NFILE, fnm),
2317                     ftp2fn_null(efEDR, NFILE, fnm),
2318                     bNeedVel,
2319                     bGenVel,
2320                     fr_time,
2321                     ir,
2322                     &state,
2323                     &sys,
2324                     oenv,
2325                     logger);
2326     }
2327
2328     if (ir->pbcType == PbcType::XY && ir->nwall != 2)
2329     {
2330         clear_rvec(state.box[ZZ]);
2331     }
2332
2333     if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2334     {
2335         /* Calculate the optimal grid dimensions */
2336         matrix          scaledBox;
2337         EwaldBoxZScaler boxScaler(*ir);
2338         boxScaler.scaleBox(state.box, scaledBox);
2339
2340         if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2341         {
2342             /* Mark fourier_spacing as not used */
2343             ir->fourier_spacing = 0;
2344         }
2345         else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2346         {
2347             set_warning_line(wi, mdparin, -1);
2348             warning_error(
2349                     wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2350         }
2351         const int minGridSize = minimalPmeGridSize(ir->pme_order);
2352         calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize, &(ir->nkx), &(ir->nky), &(ir->nkz));
2353         if (ir->nkx < minGridSize || ir->nky < minGridSize || ir->nkz < minGridSize)
2354         {
2355             warning_error(wi,
2356                           "The PME grid size should be >= 2*(pme-order - 1); either manually "
2357                           "increase the grid size or decrease pme-order");
2358         }
2359     }
2360
2361     /* MRS: eventually figure out better logic for initializing the fep
2362        values that makes declaring the lambda and declaring the state not
2363        potentially conflict if not handled correctly. */
2364     if (ir->efep != FreeEnergyPerturbationType::No)
2365     {
2366         state.fep_state = ir->fepvals->init_fep_state;
2367         for (i = 0; i < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); i++)
2368         {
2369             /* init_lambda trumps state definitions*/
2370             if (ir->fepvals->init_lambda >= 0)
2371             {
2372                 state.lambda[i] = ir->fepvals->init_lambda;
2373             }
2374             else
2375             {
2376                 if (ir->fepvals->all_lambda[i].empty())
2377                 {
2378                     gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2379                 }
2380                 else
2381                 {
2382                     state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2383                 }
2384             }
2385         }
2386     }
2387
2388     pull_t* pull = nullptr;
2389
2390     if (ir->bPull)
2391     {
2392         pull = set_pull_init(ir,
2393                              &sys,
2394                              state.x.rvec_array(),
2395                              state.box,
2396                              state.lambda[FreeEnergyPerturbationCouplingType::Mass],
2397                              wi);
2398     }
2399
2400     /* Modules that supply external potential for pull coordinates
2401      * should register those potentials here. finish_pull() will check
2402      * that providers have been registerd for all external potentials.
2403      */
2404
2405     if (ir->bDoAwh)
2406     {
2407         tensor compressibility = { { 0 } };
2408         if (ir->epc != PressureCoupling::No)
2409         {
2410             copy_mat(ir->compress, compressibility);
2411         }
2412         setStateDependentAwhParams(
2413                 ir->awhParams,
2414                 *ir->pull,
2415                 pull,
2416                 state.box,
2417                 ir->pbcType,
2418                 compressibility,
2419                 &ir->opts,
2420                 ir->efep != FreeEnergyPerturbationType::No
2421                         ? ir->fepvals->all_lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Fep)]
2422                                                  [ir->fepvals->init_fep_state]
2423                         : 0,
2424                 sys,
2425                 wi);
2426     }
2427
2428     if (ir->bPull)
2429     {
2430         finish_pull(pull);
2431     }
2432
2433     if (ir->bRot)
2434     {
2435         set_reference_positions(ir->rot,
2436                                 state.x.rvec_array(),
2437                                 state.box,
2438                                 opt2fn("-ref", NFILE, fnm),
2439                                 opt2bSet("-ref", NFILE, fnm),
2440                                 wi);
2441     }
2442
2443     /*  reset_multinr(sys); */
2444
2445     if (EEL_PME(ir->coulombtype))
2446     {
2447         float ratio = pme_load_estimate(sys, *ir, state.box);
2448         GMX_LOG(logger.info)
2449                 .asParagraph()
2450                 .appendTextFormatted(
2451                         "Estimate for the relative computational load of the PME mesh part: %.2f", ratio);
2452         /* With free energy we might need to do PME both for the A and B state
2453          * charges. This will double the cost, but the optimal performance will
2454          * then probably be at a slightly larger cut-off and grid spacing.
2455          */
2456         if ((ir->efep == FreeEnergyPerturbationType::No && ratio > 1.0 / 2.0)
2457             || (ir->efep != FreeEnergyPerturbationType::No && ratio > 2.0 / 3.0))
2458         {
2459             warning_note(
2460                     wi,
2461                     "The optimal PME mesh load for parallel simulations is below 0.5\n"
2462                     "and for highly parallel simulations between 0.25 and 0.33,\n"
2463                     "for higher performance, increase the cut-off and the PME grid spacing.\n");
2464             if (ir->efep != FreeEnergyPerturbationType::No)
2465             {
2466                 warning_note(wi,
2467                              "For free energy simulations, the optimal load limit increases from "
2468                              "0.5 to 0.667\n");
2469             }
2470         }
2471     }
2472
2473     {
2474         double      cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2475         std::string warningMessage =
2476                 gmx::formatString("This run will generate roughly %.0f Mb of data", cio);
2477         if (cio > 2000)
2478         {
2479             set_warning_line(wi, mdparin, -1);
2480             warning_note(wi, warningMessage);
2481         }
2482         else
2483         {
2484             GMX_LOG(logger.info).asParagraph().appendText(warningMessage);
2485         }
2486     }
2487
2488     // Add the md modules internal parameters that are not mdp options
2489     // e.g., atom indices
2490
2491     {
2492         gmx::KeyValueTreeBuilder internalParameterBuilder;
2493         mdModules.notifier().preProcessingNotifications_.notify(internalParameterBuilder.rootObject());
2494         ir->internalParameters =
2495                 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2496     }
2497
2498     if (ir->comm_mode != ComRemovalAlgorithm::No)
2499     {
2500         const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
2501         if (ir->nstcomm % nstglobalcomm != 0)
2502         {
2503             warning_note(
2504                     wi,
2505                     gmx::formatString(
2506                             "COM removal frequency is set to (%d).\n"
2507                             "Other settings require a global communication frequency of %d.\n"
2508                             "Note that this will require additional global communication steps,\n"
2509                             "which will reduce performance when using multiple ranks.\n"
2510                             "Consider setting nstcomm to a multiple of %d.",
2511                             ir->nstcomm,
2512                             nstglobalcomm,
2513                             nstglobalcomm));
2514         }
2515     }
2516
2517     if (bVerbose)
2518     {
2519         GMX_LOG(logger.info).asParagraph().appendTextFormatted("writing run input file...");
2520     }
2521
2522     done_warning(wi, FARGS);
2523     write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, &sys);
2524
2525     /* Output IMD group, if bIMD is TRUE */
2526     gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, &sys, NFILE, fnm);
2527
2528     sfree(opts->define);
2529     sfree(opts->wall_atomtype[0]);
2530     sfree(opts->wall_atomtype[1]);
2531     sfree(opts->include);
2532     for (auto& mol : mi)
2533     {
2534         // Some of the contents of molinfo have been stolen, so
2535         // fullCleanUp can't be called.
2536         mol.partialCleanUp();
2537     }
2538     done_inputrec_strings();
2539     output_env_done(oenv);
2540
2541     return 0;
2542 }