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