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