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