Remove unused function generate_excls and make clean_excls static
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.cpp
blob11755e631f6caab7c614e3de4cec3d228b7653d0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "grompp.h"
41 #include <cerrno>
42 #include <climits>
43 #include <cmath>
44 #include <cstring>
46 #include <algorithm>
47 #include <memory>
48 #include <vector>
50 #include <sys/types.h>
52 #include "gromacs/awh/read_params.h"
53 #include "gromacs/commandline/pargs.h"
54 #include "gromacs/ewald/ewald_utils.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/fft/calcgrid.h"
57 #include "gromacs/fileio/confio.h"
58 #include "gromacs/fileio/enxio.h"
59 #include "gromacs/fileio/tpxio.h"
60 #include "gromacs/fileio/trxio.h"
61 #include "gromacs/fileio/warninp.h"
62 #include "gromacs/gmxpreprocess/add_par.h"
63 #include "gromacs/gmxpreprocess/convparm.h"
64 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
65 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
66 #include "gromacs/gmxpreprocess/grompp_impl.h"
67 #include "gromacs/gmxpreprocess/notset.h"
68 #include "gromacs/gmxpreprocess/readir.h"
69 #include "gromacs/gmxpreprocess/tomorse.h"
70 #include "gromacs/gmxpreprocess/topio.h"
71 #include "gromacs/gmxpreprocess/toputil.h"
72 #include "gromacs/gmxpreprocess/vsite_parm.h"
73 #include "gromacs/imd/imd.h"
74 #include "gromacs/math/functions.h"
75 #include "gromacs/math/invertmatrix.h"
76 #include "gromacs/math/units.h"
77 #include "gromacs/math/vec.h"
78 #include "gromacs/mdlib/calc_verletbuf.h"
79 #include "gromacs/mdlib/compute_io.h"
80 #include "gromacs/mdlib/constr.h"
81 #include "gromacs/mdlib/perf_est.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/vsite.h"
84 #include "gromacs/mdrun/mdmodules.h"
85 #include "gromacs/mdtypes/inputrec.h"
86 #include "gromacs/mdtypes/md_enums.h"
87 #include "gromacs/mdtypes/nblist.h"
88 #include "gromacs/mdtypes/state.h"
89 #include "gromacs/pbcutil/boxutilities.h"
90 #include "gromacs/pbcutil/pbc.h"
91 #include "gromacs/pulling/pull.h"
92 #include "gromacs/random/seed.h"
93 #include "gromacs/topology/ifunc.h"
94 #include "gromacs/topology/mtop_util.h"
95 #include "gromacs/topology/symtab.h"
96 #include "gromacs/topology/topology.h"
97 #include "gromacs/trajectory/trajectoryframe.h"
98 #include "gromacs/utility/arraysize.h"
99 #include "gromacs/utility/cstringutil.h"
100 #include "gromacs/utility/exceptions.h"
101 #include "gromacs/utility/fatalerror.h"
102 #include "gromacs/utility/futil.h"
103 #include "gromacs/utility/gmxassert.h"
104 #include "gromacs/utility/keyvaluetreebuilder.h"
105 #include "gromacs/utility/mdmodulenotification.h"
106 #include "gromacs/utility/smalloc.h"
107 #include "gromacs/utility/snprintf.h"
109 /* TODO The implementation details should move to their own source file. */
110 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int> atoms,
111 gmx::ArrayRef<const real> params,
112 const std::string &name)
113 : atoms_(atoms.begin(), atoms.end()),
114 interactionTypeName_(name)
116 GMX_RELEASE_ASSERT(params.size() <= forceParam_.size(),
117 gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM).c_str());
118 auto forceParamIt = forceParam_.begin();
119 for (const auto param : params)
121 *forceParamIt++ = param;
123 while (forceParamIt != forceParam_.end())
125 *forceParamIt++ = NOTSET;
129 const int &InteractionOfType::ai() const
131 GMX_RELEASE_ASSERT(!atoms_.empty(), "Need to have at least one atom set");
132 return atoms_[0];
135 const int &InteractionOfType::aj() const
137 GMX_RELEASE_ASSERT(atoms_.size() > 1, "Need to have at least two atoms set");
138 return atoms_[1];
141 const int &InteractionOfType::ak() const
143 GMX_RELEASE_ASSERT(atoms_.size() > 2, "Need to have at least three atoms set");
144 return atoms_[2];
147 const int &InteractionOfType::al() const
149 GMX_RELEASE_ASSERT(atoms_.size() > 3, "Need to have at least four atoms set");
150 return atoms_[3];
153 const int &InteractionOfType::am() const
155 GMX_RELEASE_ASSERT(atoms_.size() > 4, "Need to have at least five atoms set");
156 return atoms_[4];
159 const real &InteractionOfType::c0() const
161 return forceParam_[0];
164 const real &InteractionOfType::c1() const
166 return forceParam_[1];
169 const real &InteractionOfType::c2() const
171 return forceParam_[2];
174 const std::string &InteractionOfType::interactionTypeName() const
176 return interactionTypeName_;
179 void InteractionOfType::sortBondAtomIds()
181 if (aj() < ai())
183 std::swap(atoms_[0], atoms_[1]);
187 void InteractionOfType::sortAngleAtomIds()
189 if (ak() < ai())
191 std::swap(atoms_[0], atoms_[2]);
195 void InteractionOfType::sortDihedralAtomIds()
197 if (al() < ai())
199 std::swap(atoms_[0], atoms_[3]);
200 std::swap(atoms_[1], atoms_[2]);
204 void InteractionOfType::sortAtomIds()
206 if (isBond())
208 sortBondAtomIds();
210 else if (isAngle())
212 sortAngleAtomIds();
214 else if (isDihedral())
216 sortDihedralAtomIds();
218 else
220 GMX_THROW(gmx::InternalError("Can not sort parameters other than bonds, angles or dihedrals"));
224 void InteractionOfType::setForceParameter(int pos, real value)
226 GMX_RELEASE_ASSERT(pos < MAXFORCEPARAM, "Can't set parameter beyond the maximum number of parameters");
227 forceParam_[pos] = value;
230 void MoleculeInformation::initMolInfo()
232 init_block(&mols);
233 init_blocka(&excls);
234 init_t_atoms(&atoms, 0, FALSE);
237 void MoleculeInformation::partialCleanUp()
239 done_block(&mols);
242 void MoleculeInformation::fullCleanUp()
244 done_atom (&atoms);
245 done_block(&mols);
248 static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
250 int n = 0;
251 /* For all the molecule types */
252 for (auto &mol : mols)
254 n += mol.interactions[ifunc].size();
255 mol.interactions[ifunc].interactionTypes.clear();
257 return n;
260 static int check_atom_names(const char *fn1, const char *fn2,
261 gmx_mtop_t *mtop, const t_atoms *at)
263 int m, i, j, nmismatch;
264 t_atoms *tat;
265 #define MAXMISMATCH 20
267 if (mtop->natoms != at->nr)
269 gmx_incons("comparing atom names");
272 nmismatch = 0;
273 i = 0;
274 for (const gmx_molblock_t &molb : mtop->molblock)
276 tat = &mtop->moltype[molb.type].atoms;
277 for (m = 0; m < molb.nmol; m++)
279 for (j = 0; j < tat->nr; j++)
281 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
283 if (nmismatch < MAXMISMATCH)
285 fprintf(stderr,
286 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
287 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
289 else if (nmismatch == MAXMISMATCH)
291 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
293 nmismatch++;
295 i++;
300 return nmismatch;
303 static void check_bonds_timestep(const gmx_mtop_t *mtop, double dt, warninp *wi)
305 /* This check is not intended to ensure accurate integration,
306 * rather it is to signal mistakes in the mdp settings.
307 * A common mistake is to forget to turn on constraints
308 * for MD after energy minimization with flexible bonds.
309 * This check can also detect too large time steps for flexible water
310 * models, but such errors will often be masked by the constraints
311 * mdp options, which turns flexible water into water with bond constraints,
312 * but without an angle constraint. Unfortunately such incorrect use
313 * of water models can not easily be detected without checking
314 * for specific model names.
316 * The stability limit of leap-frog or velocity verlet is 4.44 steps
317 * per oscillational period.
318 * But accurate bonds distributions are lost far before that limit.
319 * To allow relatively common schemes (although not common with Gromacs)
320 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
321 * we set the note limit to 10.
323 int min_steps_warn = 5;
324 int min_steps_note = 10;
325 int ftype;
326 int i, a1, a2, w_a1, w_a2, j;
327 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
328 bool bFound, bWater, bWarn;
329 char warn_buf[STRLEN];
331 /* Get the interaction parameters */
332 gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
334 twopi2 = gmx::square(2*M_PI);
336 limit2 = gmx::square(min_steps_note*dt);
338 w_a1 = w_a2 = -1;
339 w_period2 = -1.0;
341 const gmx_moltype_t *w_moltype = nullptr;
342 for (const gmx_moltype_t &moltype : mtop->moltype)
344 const t_atom *atom = moltype.atoms.atom;
345 const InteractionLists &ilist = moltype.ilist;
346 const InteractionList &ilc = ilist[F_CONSTR];
347 const InteractionList &ils = ilist[F_SETTLE];
348 for (ftype = 0; ftype < F_NRE; ftype++)
350 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
352 continue;
355 const InteractionList &ilb = ilist[ftype];
356 for (i = 0; i < ilb.size(); i += 3)
358 fc = ip[ilb.iatoms[i]].harmonic.krA;
359 re = ip[ilb.iatoms[i]].harmonic.rA;
360 if (ftype == F_G96BONDS)
362 /* Convert squared sqaure fc to harmonic fc */
363 fc = 2*fc*re;
365 a1 = ilb.iatoms[i+1];
366 a2 = ilb.iatoms[i+2];
367 m1 = atom[a1].m;
368 m2 = atom[a2].m;
369 if (fc > 0 && m1 > 0 && m2 > 0)
371 period2 = twopi2*m1*m2/((m1 + m2)*fc);
373 else
375 period2 = GMX_FLOAT_MAX;
377 if (debug)
379 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
380 fc, m1, m2, std::sqrt(period2));
382 if (period2 < limit2)
384 bFound = false;
385 for (j = 0; j < ilc.size(); j += 3)
387 if ((ilc.iatoms[j+1] == a1 && ilc.iatoms[j+2] == a2) ||
388 (ilc.iatoms[j+1] == a2 && ilc.iatoms[j+2] == a1))
390 bFound = true;
393 for (j = 0; j < ils.size(); j += 4)
395 if ((a1 == ils.iatoms[j+1] || a1 == ils.iatoms[j+2] || a1 == ils.iatoms[j+3]) &&
396 (a2 == ils.iatoms[j+1] || a2 == ils.iatoms[j+2] || a2 == ils.iatoms[j+3]))
398 bFound = true;
401 if (!bFound &&
402 (w_moltype == nullptr || period2 < w_period2))
404 w_moltype = &moltype;
405 w_a1 = a1;
406 w_a2 = a2;
407 w_period2 = period2;
414 if (w_moltype != nullptr)
416 bWarn = (w_period2 < gmx::square(min_steps_warn*dt));
417 /* A check that would recognize most water models */
418 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
419 w_moltype->atoms.nr <= 5);
420 sprintf(warn_buf, "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated oscillational period of %.1e ps, which is less than %d times the time step of %.1e ps.\n"
421 "%s",
422 *w_moltype->name,
423 w_a1+1, *w_moltype->atoms.atomname[w_a1],
424 w_a2+1, *w_moltype->atoms.atomname[w_a2],
425 std::sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
426 bWater ?
427 "Maybe you asked for fexible water." :
428 "Maybe you forgot to change the constraints mdp option.");
429 if (bWarn)
431 warning(wi, warn_buf);
433 else
435 warning_note(wi, warn_buf);
440 static void check_vel(gmx_mtop_t *mtop, rvec v[])
442 for (const AtomProxy atomP : AtomRange(*mtop))
444 const t_atom &local = atomP.atom();
445 int i = atomP.globalAtomNumber();
446 if (local.ptype == eptShell ||
447 local.ptype == eptBond ||
448 local.ptype == eptVSite)
450 clear_rvec(v[i]);
455 static void check_shells_inputrec(gmx_mtop_t *mtop,
456 t_inputrec *ir,
457 warninp *wi)
459 int nshells = 0;
460 char warn_buf[STRLEN];
462 for (const AtomProxy atomP : AtomRange(*mtop))
464 const t_atom &local = atomP.atom();
465 if (local.ptype == eptShell ||
466 local.ptype == eptBond)
468 nshells++;
471 if ((nshells > 0) && (ir->nstcalcenergy != 1))
473 set_warning_line(wi, "unknown", -1);
474 snprintf(warn_buf, STRLEN,
475 "There are %d shells, changing nstcalcenergy from %d to 1",
476 nshells, ir->nstcalcenergy);
477 ir->nstcalcenergy = 1;
478 warning(wi, warn_buf);
482 /* TODO Decide whether this function can be consolidated with
483 * gmx_mtop_ftype_count */
484 static int nint_ftype(gmx_mtop_t *mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
486 int nint = 0;
487 for (const gmx_molblock_t &molb : mtop->molblock)
489 nint += molb.nmol*mi[molb.type].interactions[ftype].size();
492 return nint;
495 /* This routine reorders the molecule type array
496 * in the order of use in the molblocks,
497 * unused molecule types are deleted.
499 static void renumber_moltypes(gmx_mtop_t *sys,
500 std::vector<MoleculeInformation> *molinfo)
503 std::vector<int> order;
504 for (gmx_molblock_t &molblock : sys->molblock)
506 const auto found = std::find_if(order.begin(), order.end(),
507 [&molblock](const auto &entry)
508 { return molblock.type == entry; });
509 if (found == order.end())
511 /* This type did not occur yet, add it */
512 order.push_back(molblock.type);
513 molblock.type = order.size() - 1;
515 else
517 molblock.type = std::distance(order.begin(), found);
521 /* We still need to reorder the molinfo structs */
522 std::vector<MoleculeInformation> minew(order.size());
523 int index = 0;
524 for (auto &mi : *molinfo)
526 const auto found = std::find(order.begin(), order.end(), index);
527 if (found != order.end())
529 int pos = std::distance(order.begin(), found);
530 minew[pos] = mi;
532 else
534 // Need to manually clean up memory ....
535 mi.fullCleanUp();
537 index++;
540 *molinfo = minew;
543 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t *mtop)
545 mtop->moltype.resize(mi.size());
546 int pos = 0;
547 for (const auto &mol : mi)
549 gmx_moltype_t &molt = mtop->moltype[pos];
550 molt.name = mol.name;
551 molt.atoms = mol.atoms;
552 /* ilists are copied later */
553 molt.excls = mol.excls;
554 pos++;
558 static void
559 new_status(const char *topfile, const char *topppfile, const char *confin,
560 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
561 bool bGenVel, bool bVerbose, t_state *state,
562 PreprocessingAtomTypes *atypes, gmx_mtop_t *sys,
563 std::vector<MoleculeInformation> *mi,
564 std::unique_ptr<MoleculeInformation> *intermolecular_interactions,
565 gmx::ArrayRef<InteractionsOfType> interactions,
566 int *comb, double *reppow, real *fudgeQQ,
567 gmx_bool bMorse,
568 warninp *wi)
570 std::vector<gmx_molblock_t> molblock;
571 int i, nmismatch;
572 bool ffParametrizedWithHBondConstraints;
573 char buf[STRLEN];
574 char warn_buf[STRLEN];
576 /* TOPOLOGY processing */
577 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
578 interactions, comb, reppow, fudgeQQ,
579 atypes, mi, intermolecular_interactions,
581 &molblock,
582 &ffParametrizedWithHBondConstraints,
583 wi);
585 sys->molblock.clear();
587 sys->natoms = 0;
588 for (const gmx_molblock_t &molb : molblock)
590 if (!sys->molblock.empty() &&
591 molb.type == sys->molblock.back().type)
593 /* Merge consecutive blocks with the same molecule type */
594 sys->molblock.back().nmol += molb.nmol;
596 else if (molb.nmol > 0)
598 /* Add a new molblock to the topology */
599 sys->molblock.push_back(molb);
601 sys->natoms += molb.nmol*(*mi)[sys->molblock.back().type].atoms.nr;
603 if (sys->molblock.empty())
605 gmx_fatal(FARGS, "No molecules were defined in the system");
608 renumber_moltypes(sys, mi);
610 if (bMorse)
612 convert_harmonics(*mi, atypes);
615 if (ir->eDisre == edrNone)
617 i = rm_interactions(F_DISRES, *mi);
618 if (i > 0)
620 set_warning_line(wi, "unknown", -1);
621 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
622 warning_note(wi, warn_buf);
625 if (!opts->bOrire)
627 i = rm_interactions(F_ORIRES, *mi);
628 if (i > 0)
630 set_warning_line(wi, "unknown", -1);
631 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
632 warning_note(wi, warn_buf);
636 /* Copy structures from msys to sys */
637 molinfo2mtop(*mi, sys);
639 gmx_mtop_finalize(sys);
641 /* Discourage using the, previously common, setup of applying constraints
642 * to all bonds with force fields that have been parametrized with H-bond
643 * constraints only. Do not print note with large timesteps or vsites.
645 if (opts->nshake == eshALLBONDS &&
646 ffParametrizedWithHBondConstraints &&
647 ir->delta_t < 0.0026 &&
648 gmx_mtop_ftype_count(sys, F_VSITE3FD) == 0)
650 set_warning_line(wi, "unknown", -1);
651 warning_note(wi, "You are using constraints on all bonds, whereas the forcefield "
652 "has been parametrized only with constraints involving hydrogen atoms. "
653 "We suggest using constraints = h-bonds instead, this will also improve performance.");
656 /* COORDINATE file processing */
657 if (bVerbose)
659 fprintf(stderr, "processing coordinates...\n");
662 t_topology *conftop;
663 rvec *x = nullptr;
664 rvec *v = nullptr;
665 snew(conftop, 1);
666 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
667 state->natoms = conftop->atoms.nr;
668 if (state->natoms != sys->natoms)
670 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
671 " does not match topology (%s, %d)",
672 confin, state->natoms, topfile, sys->natoms);
674 /* It would be nice to get rid of the copies below, but we don't know
675 * a priori if the number of atoms in confin matches what we expect.
677 state->flags |= (1 << estX);
678 if (v != nullptr)
680 state->flags |= (1 << estV);
682 state_change_natoms(state, state->natoms);
683 std::copy(x, x+state->natoms, state->x.data());
684 sfree(x);
685 if (v != nullptr)
687 std::copy(v, v+state->natoms, state->v.data());
688 sfree(v);
690 /* This call fixes the box shape for runs with pressure scaling */
691 set_box_rel(ir, state);
693 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
694 done_top(conftop);
695 sfree(conftop);
697 if (nmismatch)
699 sprintf(buf, "%d non-matching atom name%s\n"
700 "atom names from %s will be used\n"
701 "atom names from %s will be ignored\n",
702 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
703 warning(wi, buf);
706 /* Do more checks, mostly related to constraints */
707 if (bVerbose)
709 fprintf(stderr, "double-checking input for internal consistency...\n");
712 bool bHasNormalConstraints = 0 < (nint_ftype(sys, *mi, F_CONSTR) +
713 nint_ftype(sys, *mi, F_CONSTRNC));
714 bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
715 double_check(ir, state->box,
716 bHasNormalConstraints,
717 bHasAnyConstraints,
718 wi);
721 if (bGenVel)
723 real *mass;
725 snew(mass, state->natoms);
726 for (const AtomProxy atomP : AtomRange(*sys))
728 const t_atom &local = atomP.atom();
729 int i = atomP.globalAtomNumber();
730 mass[i] = local.m;
733 if (opts->seed == -1)
735 opts->seed = static_cast<int>(gmx::makeRandomSeed());
736 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
738 state->flags |= (1 << estV);
739 maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array());
741 stop_cm(stdout, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
742 sfree(mass);
746 static void copy_state(const char *slog, t_trxframe *fr,
747 bool bReadVel, t_state *state,
748 double *use_time)
750 if (fr->not_ok & FRAME_NOT_OK)
752 gmx_fatal(FARGS, "Can not start from an incomplete frame");
754 if (!fr->bX)
756 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
757 slog);
760 std::copy(fr->x, fr->x+state->natoms, state->x.data());
761 if (bReadVel)
763 if (!fr->bV)
765 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
767 std::copy(fr->v, fr->v+state->natoms, state->v.data());
769 if (fr->bBox)
771 copy_mat(fr->box, state->box);
774 *use_time = fr->time;
777 static void cont_status(const char *slog, const char *ener,
778 bool bNeedVel, bool bGenVel, real fr_time,
779 t_inputrec *ir, t_state *state,
780 gmx_mtop_t *sys,
781 const gmx_output_env_t *oenv)
782 /* If fr_time == -1 read the last frame available which is complete */
784 bool bReadVel;
785 t_trxframe fr;
786 t_trxstatus *fp;
787 double use_time;
789 bReadVel = (bNeedVel && !bGenVel);
791 fprintf(stderr,
792 "Reading Coordinates%s and Box size from old trajectory\n",
793 bReadVel ? ", Velocities" : "");
794 if (fr_time == -1)
796 fprintf(stderr, "Will read whole trajectory\n");
798 else
800 fprintf(stderr, "Will read till time %g\n", fr_time);
802 if (!bReadVel)
804 if (bGenVel)
806 fprintf(stderr, "Velocities generated: "
807 "ignoring velocities in input trajectory\n");
809 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
811 else
813 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
815 if (!fr.bV)
817 fprintf(stderr,
818 "\n"
819 "WARNING: Did not find a frame with velocities in file %s,\n"
820 " all velocities will be set to zero!\n\n", slog);
821 for (auto &vi : makeArrayRef(state->v))
823 vi = {0, 0, 0};
825 close_trx(fp);
826 /* Search for a frame without velocities */
827 bReadVel = false;
828 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
832 state->natoms = fr.natoms;
834 if (sys->natoms != state->natoms)
836 gmx_fatal(FARGS, "Number of atoms in Topology "
837 "is not the same as in Trajectory");
839 copy_state(slog, &fr, bReadVel, state, &use_time);
841 /* Find the appropriate frame */
842 while ((fr_time == -1 || fr.time < fr_time) &&
843 read_next_frame(oenv, fp, &fr))
845 copy_state(slog, &fr, bReadVel, state, &use_time);
848 close_trx(fp);
850 /* Set the relative box lengths for preserving the box shape.
851 * Note that this call can lead to differences in the last bit
852 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
854 set_box_rel(ir, state);
856 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
857 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
859 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
861 get_enx_state(ener, use_time, sys->groups, ir, state);
862 preserve_box_shape(ir, state->box_rel, state->boxv);
866 static void read_posres(gmx_mtop_t *mtop,
867 gmx::ArrayRef<const MoleculeInformation> molinfo,
868 gmx_bool bTopB,
869 const char *fn,
870 int rc_scaling, int ePBC,
871 rvec com,
872 warninp *wi)
874 gmx_bool *hadAtom;
875 rvec *x, *v;
876 dvec sum;
877 double totmass;
878 t_topology *top;
879 matrix box, invbox;
880 int natoms, npbcdim = 0;
881 char warn_buf[STRLEN];
882 int a, nat_molb;
883 t_atom *atom;
885 snew(top, 1);
886 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
887 natoms = top->atoms.nr;
888 done_top(top);
889 sfree(top);
890 if (natoms != mtop->natoms)
892 sprintf(warn_buf, "The number of atoms in %s (%d) does not match the number of atoms in the topology (%d). Will assume that the first %d atoms in the topology and %s match.", fn, natoms, mtop->natoms, std::min(mtop->natoms, natoms), fn);
893 warning(wi, warn_buf);
896 npbcdim = ePBC2npbcdim(ePBC);
897 GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
898 clear_rvec(com);
899 if (rc_scaling != erscNO)
901 copy_mat(box, invbox);
902 for (int j = npbcdim; j < DIM; j++)
904 clear_rvec(invbox[j]);
905 invbox[j][j] = 1;
907 gmx::invertBoxMatrix(invbox, invbox);
910 /* Copy the reference coordinates to mtop */
911 clear_dvec(sum);
912 totmass = 0;
913 a = 0;
914 snew(hadAtom, natoms);
915 for (gmx_molblock_t &molb : mtop->molblock)
917 nat_molb = molb.nmol*mtop->moltype[molb.type].atoms.nr;
918 const InteractionsOfType *pr = &(molinfo[molb.type].interactions[F_POSRES]);
919 const InteractionsOfType *prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
920 if (pr->size() > 0 || prfb->size() > 0)
922 atom = mtop->moltype[molb.type].atoms.atom;
923 for (const auto &restraint : pr->interactionTypes)
925 int ai = restraint.ai();
926 if (ai >= natoms)
928 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
929 ai+1, *molinfo[molb.type].name, fn, natoms);
931 hadAtom[ai] = TRUE;
932 if (rc_scaling == erscCOM)
934 /* Determine the center of mass of the posres reference coordinates */
935 for (int j = 0; j < npbcdim; j++)
937 sum[j] += atom[ai].m*x[a+ai][j];
939 totmass += atom[ai].m;
942 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
943 for (const auto &restraint : prfb->interactionTypes)
945 int ai = restraint.ai();
946 if (ai >= natoms)
948 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
949 ai+1, *molinfo[molb.type].name, fn, natoms);
951 if (rc_scaling == erscCOM && !hadAtom[ai])
953 /* Determine the center of mass of the posres reference coordinates */
954 for (int j = 0; j < npbcdim; j++)
956 sum[j] += atom[ai].m*x[a+ai][j];
958 totmass += atom[ai].m;
961 if (!bTopB)
963 molb.posres_xA.resize(nat_molb);
964 for (int i = 0; i < nat_molb; i++)
966 copy_rvec(x[a+i], molb.posres_xA[i]);
969 else
971 molb.posres_xB.resize(nat_molb);
972 for (int i = 0; i < nat_molb; i++)
974 copy_rvec(x[a+i], molb.posres_xB[i]);
978 a += nat_molb;
980 if (rc_scaling == erscCOM)
982 if (totmass == 0)
984 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
986 for (int j = 0; j < npbcdim; j++)
988 com[j] = sum[j]/totmass;
990 fprintf(stderr, "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n", com[XX], com[YY], com[ZZ]);
993 if (rc_scaling != erscNO)
995 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
997 for (gmx_molblock_t &molb : mtop->molblock)
999 nat_molb = molb.nmol*mtop->moltype[molb.type].atoms.nr;
1000 if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1002 std::vector<gmx::RVec> &xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1003 for (int i = 0; i < nat_molb; i++)
1005 for (int j = 0; j < npbcdim; j++)
1007 if (rc_scaling == erscALL)
1009 /* Convert from Cartesian to crystal coordinates */
1010 xp[i][j] *= invbox[j][j];
1011 for (int k = j+1; k < npbcdim; k++)
1013 xp[i][j] += invbox[k][j]*xp[i][k];
1016 else if (rc_scaling == erscCOM)
1018 /* Subtract the center of mass */
1019 xp[i][j] -= com[j];
1026 if (rc_scaling == erscCOM)
1028 /* Convert the COM from Cartesian to crystal coordinates */
1029 for (int j = 0; j < npbcdim; j++)
1031 com[j] *= invbox[j][j];
1032 for (int k = j+1; k < npbcdim; k++)
1034 com[j] += invbox[k][j]*com[k];
1040 sfree(x);
1041 sfree(v);
1042 sfree(hadAtom);
1045 static void gen_posres(gmx_mtop_t *mtop,
1046 gmx::ArrayRef<const MoleculeInformation> mi,
1047 const char *fnA, const char *fnB,
1048 int rc_scaling, int ePBC,
1049 rvec com, rvec comB,
1050 warninp *wi)
1052 read_posres(mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1053 /* It is safer to simply read the b-state posres rather than trying
1054 * to be smart and copy the positions.
1056 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1059 static void set_wall_atomtype(PreprocessingAtomTypes *at, t_gromppopts *opts,
1060 t_inputrec *ir, warninp *wi)
1062 int i;
1063 char warn_buf[STRLEN];
1065 if (ir->nwall > 0)
1067 fprintf(stderr, "Searching the wall atom type(s)\n");
1069 for (i = 0; i < ir->nwall; i++)
1071 ir->wall_atomtype[i] = at->atomTypeFromName(opts->wall_atomtype[i]);
1072 if (ir->wall_atomtype[i] == NOTSET)
1074 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1075 warning_error(wi, warn_buf);
1080 static int nrdf_internal(const t_atoms *atoms)
1082 int i, nmass, nrdf;
1084 nmass = 0;
1085 for (i = 0; i < atoms->nr; i++)
1087 /* Vsite ptype might not be set here yet, so also check the mass */
1088 if ((atoms->atom[i].ptype == eptAtom ||
1089 atoms->atom[i].ptype == eptNucleus)
1090 && atoms->atom[i].m > 0)
1092 nmass++;
1095 switch (nmass)
1097 case 0: nrdf = 0; break;
1098 case 1: nrdf = 0; break;
1099 case 2: nrdf = 1; break;
1100 default: nrdf = nmass*3 - 6; break;
1103 return nrdf;
1106 static void
1107 spline1d( double dx,
1108 const double * y,
1109 int n,
1110 double * u,
1111 double * y2 )
1113 int i;
1114 double p, q;
1116 y2[0] = 0.0;
1117 u[0] = 0.0;
1119 for (i = 1; i < n-1; i++)
1121 p = 0.5*y2[i-1]+2.0;
1122 y2[i] = -0.5/p;
1123 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1124 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1127 y2[n-1] = 0.0;
1129 for (i = n-2; i >= 0; i--)
1131 y2[i] = y2[i]*y2[i+1]+u[i];
1136 static void
1137 interpolate1d( double xmin,
1138 double dx,
1139 const double * ya,
1140 const double * y2a,
1141 double x,
1142 double * y,
1143 double * y1)
1145 int ix;
1146 double a, b;
1148 ix = static_cast<int>((x-xmin)/dx);
1150 a = (xmin+(ix+1)*dx-x)/dx;
1151 b = (x-xmin-ix*dx)/dx;
1153 *y = a*ya[ix]+b*ya[ix+1]+((a*a*a-a)*y2a[ix]+(b*b*b-b)*y2a[ix+1])*(dx*dx)/6.0;
1154 *y1 = (ya[ix+1]-ya[ix])/dx-(3.0*a*a-1.0)/6.0*dx*y2a[ix]+(3.0*b*b-1.0)/6.0*dx*y2a[ix+1];
1158 static void
1159 setup_cmap (int grid_spacing,
1160 int nc,
1161 gmx::ArrayRef<const real> grid,
1162 gmx_cmap_t * cmap_grid)
1164 int i, j, k, ii, jj, kk, idx;
1165 int offset;
1166 double dx, xmin, v, v1, v2, v12;
1167 double phi, psi;
1169 std::vector<double> tmp_u(2*grid_spacing, 0.0);
1170 std::vector<double> tmp_u2(2*grid_spacing, 0.0);
1171 std::vector<double> tmp_yy(2*grid_spacing, 0.0);
1172 std::vector<double> tmp_y1(2*grid_spacing, 0.0);
1173 std::vector<double> tmp_t2(2*grid_spacing*2*grid_spacing, 0.0);
1174 std::vector<double> tmp_grid(2*grid_spacing*2*grid_spacing, 0.0);
1176 dx = 360.0/grid_spacing;
1177 xmin = -180.0-dx*grid_spacing/2;
1179 for (kk = 0; kk < nc; kk++)
1181 /* Compute an offset depending on which cmap we are using
1182 * Offset will be the map number multiplied with the
1183 * grid_spacing * grid_spacing * 2
1185 offset = kk * grid_spacing * grid_spacing * 2;
1187 for (i = 0; i < 2*grid_spacing; i++)
1189 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1191 for (j = 0; j < 2*grid_spacing; j++)
1193 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1194 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1198 for (i = 0; i < 2*grid_spacing; i++)
1200 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u.data(), &(tmp_t2[2*grid_spacing*i]));
1203 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1205 ii = i-grid_spacing/2;
1206 phi = ii*dx-180.0;
1208 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1210 jj = j-grid_spacing/2;
1211 psi = jj*dx-180.0;
1213 for (k = 0; k < 2*grid_spacing; k++)
1215 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1216 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1219 spline1d(dx, tmp_yy.data(), 2*grid_spacing, tmp_u.data(), tmp_u2.data());
1220 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1221 spline1d(dx, tmp_y1.data(), 2*grid_spacing, tmp_u.data(), tmp_u2.data());
1222 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1224 idx = ii*grid_spacing+jj;
1225 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1226 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1227 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1228 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1234 static void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1236 int i, nelem;
1238 cmap_grid->grid_spacing = grid_spacing;
1239 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1241 cmap_grid->cmapdata.resize(ngrid);
1243 for (i = 0; i < ngrid; i++)
1245 cmap_grid->cmapdata[i].cmap.resize(4*nelem);
1250 static int count_constraints(const gmx_mtop_t *mtop,
1251 gmx::ArrayRef<const MoleculeInformation> mi,
1252 warninp *wi)
1254 int count, count_mol;
1255 char buf[STRLEN];
1257 count = 0;
1258 for (const gmx_molblock_t &molb : mtop->molblock)
1260 count_mol = 0;
1261 gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1263 for (int i = 0; i < F_NRE; i++)
1265 if (i == F_SETTLE)
1267 count_mol += 3*interactions[i].size();
1269 else if (interaction_function[i].flags & IF_CONSTRAINT)
1271 count_mol += interactions[i].size();
1275 if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1277 sprintf(buf,
1278 "Molecule type '%s' has %d constraints.\n"
1279 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1280 *mi[molb.type].name, count_mol,
1281 nrdf_internal(&mi[molb.type].atoms));
1282 warning(wi, buf);
1284 count += molb.nmol*count_mol;
1287 return count;
1290 static real calc_temp(const gmx_mtop_t *mtop,
1291 const t_inputrec *ir,
1292 rvec *v)
1294 double sum_mv2 = 0;
1295 for (const AtomProxy atomP : AtomRange(*mtop))
1297 const t_atom &local = atomP.atom();
1298 int i = atomP.globalAtomNumber();
1299 sum_mv2 += local.m*norm2(v[i]);
1302 double nrdf = 0;
1303 for (int g = 0; g < ir->opts.ngtc; g++)
1305 nrdf += ir->opts.nrdf[g];
1308 return sum_mv2/(nrdf*BOLTZ);
1311 static real get_max_reference_temp(const t_inputrec *ir,
1312 warninp *wi)
1314 real ref_t;
1315 int i;
1316 bool bNoCoupl;
1318 ref_t = 0;
1319 bNoCoupl = false;
1320 for (i = 0; i < ir->opts.ngtc; i++)
1322 if (ir->opts.tau_t[i] < 0)
1324 bNoCoupl = true;
1326 else
1328 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1332 if (bNoCoupl)
1334 char buf[STRLEN];
1336 sprintf(buf, "Some temperature coupling groups do not use temperature coupling. We will assume their temperature is not more than %.3f K. If their temperature is higher, the energy error and the Verlet buffer might be underestimated.",
1337 ref_t);
1338 warning(wi, buf);
1341 return ref_t;
1344 /* Checks if there are unbound atoms in moleculetype molt.
1345 * Prints a note for each unbound atoms and a warning if any is present.
1347 static void checkForUnboundAtoms(const gmx_moltype_t *molt,
1348 gmx_bool bVerbose,
1349 warninp *wi)
1351 const t_atoms *atoms = &molt->atoms;
1353 if (atoms->nr == 1)
1355 /* Only one atom, there can't be unbound atoms */
1356 return;
1359 std::vector<int> count(atoms->nr, 0);
1361 for (int ftype = 0; ftype < F_NRE; ftype++)
1363 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS) ||
1364 (interaction_function[ftype].flags & IF_CONSTRAINT) ||
1365 ftype == F_SETTLE)
1367 const InteractionList &il = molt->ilist[ftype];
1368 const int nral = NRAL(ftype);
1370 for (int i = 0; i < il.size(); i += 1 + nral)
1372 for (int j = 0; j < nral; j++)
1374 count[il.iatoms[i + 1 + j]]++;
1380 int numDanglingAtoms = 0;
1381 for (int a = 0; a < atoms->nr; a++)
1383 if (atoms->atom[a].ptype != eptVSite &&
1384 count[a] == 0)
1386 if (bVerbose)
1388 fprintf(stderr, "\nAtom %d '%s' in moleculetype '%s' is not bound by a potential or constraint to any other atom in the same moleculetype.\n",
1389 a + 1, *atoms->atomname[a], *molt->name);
1391 numDanglingAtoms++;
1395 if (numDanglingAtoms > 0)
1397 char buf[STRLEN];
1398 sprintf(buf, "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any other atom in the same moleculetype. Although technically this might not cause issues in a simulation, this often means that the user forgot to add a bond/potential/constraint or put multiple molecules in the same moleculetype definition by mistake. Run with -v to get information for each atom.",
1399 *molt->name, numDanglingAtoms);
1400 warning_note(wi, buf);
1404 /* Checks all moleculetypes for unbound atoms */
1405 static void checkForUnboundAtoms(const gmx_mtop_t *mtop,
1406 gmx_bool bVerbose,
1407 warninp *wi)
1409 for (const gmx_moltype_t &molt : mtop->moltype)
1411 checkForUnboundAtoms(&molt, bVerbose, wi);
1415 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1417 * The specific decoupled modes this routine check for are angle modes
1418 * where the two bonds are constrained and the atoms a both ends are only
1419 * involved in a single constraint; the mass of the two atoms needs to
1420 * differ by more than \p massFactorThreshold.
1422 static bool haveDecoupledModeInMol(const gmx_moltype_t &molt,
1423 gmx::ArrayRef<const t_iparams> iparams,
1424 real massFactorThreshold)
1426 if (molt.ilist[F_CONSTR].size() == 0 &&
1427 molt.ilist[F_CONSTRNC].size() == 0)
1429 return false;
1432 const t_atom * atom = molt.atoms.atom;
1434 t_blocka atomToConstraints =
1435 gmx::make_at2con(molt, iparams,
1436 gmx::FlexibleConstraintTreatment::Exclude);
1438 bool haveDecoupledMode = false;
1439 for (int ftype = 0; ftype < F_NRE; ftype++)
1441 if (interaction_function[ftype].flags & IF_ATYPE)
1443 const int nral = NRAL(ftype);
1444 const InteractionList &il = molt.ilist[ftype];
1445 for (int i = 0; i < il.size(); i += 1 + nral)
1447 /* Here we check for the mass difference between the atoms
1448 * at both ends of the angle, that the atoms at the ends have
1449 * 1 contraint and the atom in the middle at least 3; we check
1450 * that the 3 atoms are linked by constraints below.
1451 * We check for at least three constraints for the middle atom,
1452 * since with only the two bonds in the angle, we have 3-atom
1453 * molecule, which has much more thermal exhange in this single
1454 * angle mode than molecules with more atoms.
1455 * Note that this check also catches molecules where more atoms
1456 * are connected to one or more atoms in the angle, but by
1457 * bond potentials instead of angles. But such cases will not
1458 * occur in "normal" molecules and it doesn't hurt running
1459 * those with higher accuracy settings as well.
1461 int a0 = il.iatoms[1 + i];
1462 int a1 = il.iatoms[1 + i + 1];
1463 int a2 = il.iatoms[1 + i + 2];
1464 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
1465 atom[a2].m > atom[a0].m*massFactorThreshold) &&
1466 atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
1467 atomToConstraints.index[a2 + 1] - atomToConstraints.index[a2] == 1 &&
1468 atomToConstraints.index[a1 + 1] - atomToConstraints.index[a1] >= 3)
1470 int constraint0 = atomToConstraints.a[atomToConstraints.index[a0]];
1471 int constraint2 = atomToConstraints.a[atomToConstraints.index[a2]];
1473 bool foundAtom0 = false;
1474 bool foundAtom2 = false;
1475 for (int conIndex = atomToConstraints.index[a1]; conIndex < atomToConstraints.index[a1 + 1]; conIndex++)
1477 if (atomToConstraints.a[conIndex] == constraint0)
1479 foundAtom0 = true;
1481 if (atomToConstraints.a[conIndex] == constraint2)
1483 foundAtom2 = true;
1486 if (foundAtom0 && foundAtom2)
1488 haveDecoupledMode = true;
1495 done_blocka(&atomToConstraints);
1497 return haveDecoupledMode;
1500 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1502 * When decoupled modes are present and the accuray in insufficient,
1503 * this routine issues a warning if the accuracy is insufficient.
1505 static void checkDecoupledModeAccuracy(const gmx_mtop_t *mtop,
1506 const t_inputrec *ir,
1507 warninp *wi)
1509 /* We only have issues with decoupled modes with normal MD.
1510 * With stochastic dynamics equipartitioning is enforced strongly.
1512 if (!EI_MD(ir->eI))
1514 return;
1517 /* When atoms of very different mass are involved in an angle potential
1518 * and both bonds in the angle are constrained, the dynamic modes in such
1519 * angles have very different periods and significant energy exchange
1520 * takes several nanoseconds. Thus even a small amount of error in
1521 * different algorithms can lead to violation of equipartitioning.
1522 * The parameters below are mainly based on an all-atom chloroform model
1523 * with all bonds constrained. Equipartitioning is off by more than 1%
1524 * (need to run 5-10 ns) when the difference in mass between H and Cl
1525 * is more than a factor 13 and the accuracy is less than the thresholds
1526 * given below. This has been verified on some other molecules.
1528 * Note that the buffer and shake parameters have unit length and
1529 * energy/time, respectively, so they will "only" work correctly
1530 * for atomistic force fields using MD units.
1532 const real massFactorThreshold = 13.0;
1533 const real bufferToleranceThreshold = 1e-4;
1534 const int lincsIterationThreshold = 2;
1535 const int lincsOrderThreshold = 4;
1536 const real shakeToleranceThreshold = 0.005*ir->delta_t;
1538 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1539 bool shakeWithSufficientTolerance = (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1*shakeToleranceThreshold);
1540 if (ir->cutoff_scheme == ecutsVERLET &&
1541 ir->verletbuf_tol <= 1.1*bufferToleranceThreshold &&
1542 (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1544 return;
1547 bool haveDecoupledMode = false;
1548 for (const gmx_moltype_t &molt : mtop->moltype)
1550 if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams,
1551 massFactorThreshold))
1553 haveDecoupledMode = true;
1557 if (haveDecoupledMode)
1559 std::string message = gmx::formatString(
1560 "There are atoms at both ends of an angle, connected by constraints "
1561 "and with masses that differ by more than a factor of %g. This means "
1562 "that there are likely dynamic modes that are only very weakly coupled.",
1563 massFactorThreshold);
1564 if (ir->cutoff_scheme == ecutsVERLET)
1566 message += gmx::formatString(
1567 " To ensure good equipartitioning, you need to either not use "
1568 "constraints on all bonds (but, if possible, only on bonds involving "
1569 "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1570 "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1571 ">= %d or SHAKE tolerance <= %g",
1572 ei_names[eiSD1],
1573 bufferToleranceThreshold,
1574 lincsIterationThreshold, lincsOrderThreshold,
1575 shakeToleranceThreshold);
1577 else
1579 message += gmx::formatString(
1580 " To ensure good equipartitioning, we suggest to switch to the %s "
1581 "cutoff-scheme, since that allows for better control over the Verlet "
1582 "buffer size and thus over the energy drift.",
1583 ecutscheme_names[ecutsVERLET]);
1585 warning(wi, message);
1589 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1590 t_inputrec *ir,
1591 real buffer_temp,
1592 matrix box,
1593 warninp *wi)
1595 char warn_buf[STRLEN];
1597 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1599 /* Calculate the buffer size for simple atom vs atoms list */
1600 VerletbufListSetup listSetup1x1;
1601 listSetup1x1.cluster_size_i = 1;
1602 listSetup1x1.cluster_size_j = 1;
1603 const real rlist_1x1 =
1604 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1605 buffer_temp, listSetup1x1);
1607 /* Set the pair-list buffer size in ir */
1608 VerletbufListSetup listSetup4x4 =
1609 verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1610 ir->rlist =
1611 calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1612 buffer_temp, listSetup4x4);
1614 const int n_nonlin_vsite = countNonlinearVsites(*mtop);
1615 if (n_nonlin_vsite > 0)
1617 sprintf(warn_buf, "There are %d non-linear virtual site constructions. Their contribution to the energy error is approximated. In most cases this does not affect the error significantly.", n_nonlin_vsite);
1618 warning_note(wi, warn_buf);
1621 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1622 1, 1, rlist_1x1, rlist_1x1-std::max(ir->rvdw, ir->rcoulomb));
1624 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1625 listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j,
1626 ir->rlist, ir->rlist-std::max(ir->rvdw, ir->rcoulomb));
1628 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1630 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1632 gmx_fatal(FARGS, "The pair-list cut-off (%g nm) is longer than half the shortest box vector or longer than the smallest box diagonal element (%g nm). Increase the box size or decrease nstlist or increase verlet-buffer-tolerance.", ir->rlist, std::sqrt(max_cutoff2(ir->ePBC, box)));
1636 int gmx_grompp(int argc, char *argv[])
1638 const char *desc[] = {
1639 "[THISMODULE] (the gromacs preprocessor)",
1640 "reads a molecular topology file, checks the validity of the",
1641 "file, expands the topology from a molecular description to an atomic",
1642 "description. The topology file contains information about",
1643 "molecule types and the number of molecules, the preprocessor",
1644 "copies each molecule as needed. ",
1645 "There is no limitation on the number of molecule types. ",
1646 "Bonds and bond-angles can be converted into constraints, separately",
1647 "for hydrogens and heavy atoms.",
1648 "Then a coordinate file is read and velocities can be generated",
1649 "from a Maxwellian distribution if requested.",
1650 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1651 "(eg. number of MD steps, time step, cut-off), and others such as",
1652 "NEMD parameters, which are corrected so that the net acceleration",
1653 "is zero.",
1654 "Eventually a binary file is produced that can serve as the sole input",
1655 "file for the MD program.[PAR]",
1657 "[THISMODULE] uses the atom names from the topology file. The atom names",
1658 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1659 "warnings when they do not match the atom names in the topology.",
1660 "Note that the atom names are irrelevant for the simulation as",
1661 "only the atom types are used for generating interaction parameters.[PAR]",
1663 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1664 "etc. The preprocessor supports the following keywords::",
1666 " #ifdef VARIABLE",
1667 " #ifndef VARIABLE",
1668 " #else",
1669 " #endif",
1670 " #define VARIABLE",
1671 " #undef VARIABLE",
1672 " #include \"filename\"",
1673 " #include <filename>",
1675 "The functioning of these statements in your topology may be modulated by",
1676 "using the following two flags in your [REF].mdp[ref] file::",
1678 " define = -DVARIABLE1 -DVARIABLE2",
1679 " include = -I/home/john/doe",
1681 "For further information a C-programming textbook may help you out.",
1682 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1683 "topology file written out so that you can verify its contents.[PAR]",
1685 "When using position restraints, a file with restraint coordinates",
1686 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1687 "for [TT]-c[tt]). For free energy calculations, separate reference",
1688 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1689 "otherwise they will be equal to those of the A topology.[PAR]",
1691 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1692 "The last frame with coordinates and velocities will be read,",
1693 "unless the [TT]-time[tt] option is used. Only if this information",
1694 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1695 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1696 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1697 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1698 "variables.[PAR]",
1700 "[THISMODULE] can be used to restart simulations (preserving",
1701 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1702 "However, for simply changing the number of run steps to extend",
1703 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1704 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1705 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1706 "like output frequency, then supplying the checkpoint file to",
1707 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1708 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1709 "the ensemble (if possible) still requires passing the checkpoint",
1710 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1712 "By default, all bonded interactions which have constant energy due to",
1713 "virtual site constructions will be removed. If this constant energy is",
1714 "not zero, this will result in a shift in the total energy. All bonded",
1715 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1716 "all constraints for distances which will be constant anyway because",
1717 "of virtual site constructions will be removed. If any constraints remain",
1718 "which involve virtual sites, a fatal error will result.[PAR]",
1720 "To verify your run input file, please take note of all warnings",
1721 "on the screen, and correct where necessary. Do also look at the contents",
1722 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1723 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1724 "with the [TT]-debug[tt] option which will give you more information",
1725 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1726 "can see the contents of the run input file with the [gmx-dump]",
1727 "program. [gmx-check] can be used to compare the contents of two",
1728 "run input files.[PAR]",
1730 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1731 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1732 "harmless, but usually they are not. The user is advised to carefully",
1733 "interpret the output messages before attempting to bypass them with",
1734 "this option."
1736 t_gromppopts *opts;
1737 std::vector<MoleculeInformation> mi;
1738 std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1739 int nvsite, comb;
1740 real fudgeQQ;
1741 double reppow;
1742 const char *mdparin;
1743 bool bNeedVel, bGenVel;
1744 gmx_bool have_atomnumber;
1745 gmx_output_env_t *oenv;
1746 gmx_bool bVerbose = FALSE;
1747 warninp *wi;
1748 char warn_buf[STRLEN];
1750 t_filenm fnm[] = {
1751 { efMDP, nullptr, nullptr, ffREAD },
1752 { efMDP, "-po", "mdout", ffWRITE },
1753 { efSTX, "-c", nullptr, ffREAD },
1754 { efSTX, "-r", "restraint", ffOPTRD },
1755 { efSTX, "-rb", "restraint", ffOPTRD },
1756 { efNDX, nullptr, nullptr, ffOPTRD },
1757 { efTOP, nullptr, nullptr, ffREAD },
1758 { efTOP, "-pp", "processed", ffOPTWR },
1759 { efTPR, "-o", nullptr, ffWRITE },
1760 { efTRN, "-t", nullptr, ffOPTRD },
1761 { efEDR, "-e", nullptr, ffOPTRD },
1762 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1763 { efGRO, "-imd", "imdgroup", ffOPTWR },
1764 { efTRN, "-ref", "rotref", ffOPTRW }
1766 #define NFILE asize(fnm)
1768 /* Command line options */
1769 gmx_bool bRenum = TRUE;
1770 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1771 int i, maxwarn = 0;
1772 real fr_time = -1;
1773 t_pargs pa[] = {
1774 { "-v", FALSE, etBOOL, {&bVerbose},
1775 "Be loud and noisy" },
1776 { "-time", FALSE, etREAL, {&fr_time},
1777 "Take frame at or first after this time." },
1778 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1779 "Remove constant bonded interactions with virtual sites" },
1780 { "-maxwarn", FALSE, etINT, {&maxwarn},
1781 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1782 { "-zero", FALSE, etBOOL, {&bZero},
1783 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1784 { "-renum", FALSE, etBOOL, {&bRenum},
1785 "Renumber atomtypes and minimize number of atomtypes" }
1788 /* Parse the command line */
1789 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1790 asize(desc), desc, 0, nullptr, &oenv))
1792 return 0;
1795 /* Initiate some variables */
1796 gmx::MDModules mdModules;
1797 t_inputrec irInstance;
1798 t_inputrec *ir = &irInstance;
1799 snew(opts, 1);
1800 snew(opts->include, STRLEN);
1801 snew(opts->define, STRLEN);
1803 wi = init_warning(TRUE, maxwarn);
1805 /* PARAMETER file processing */
1806 mdparin = opt2fn("-f", NFILE, fnm);
1807 set_warning_line(wi, mdparin, -1);
1810 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1812 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1814 if (bVerbose)
1816 fprintf(stderr, "checking input for internal consistency...\n");
1818 check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
1820 if (ir->ld_seed == -1)
1822 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1823 fprintf(stderr, "Setting the LD random seed to %" PRId64 "\n", ir->ld_seed);
1826 if (ir->expandedvals->lmc_seed == -1)
1828 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1829 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1832 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1833 bGenVel = (bNeedVel && opts->bGenVel);
1834 if (bGenVel && ir->bContinuation)
1836 sprintf(warn_buf,
1837 "Generating velocities is inconsistent with attempting "
1838 "to continue a previous run. Choose only one of "
1839 "gen-vel = yes and continuation = yes.");
1840 warning_error(wi, warn_buf);
1843 std::array<InteractionsOfType, F_NRE> interactions;
1844 gmx_mtop_t sys;
1845 PreprocessingAtomTypes atypes;
1846 if (debug)
1848 pr_symtab(debug, 0, "Just opened", &sys.symtab);
1851 const char *fn = ftp2fn(efTOP, NFILE, fnm);
1852 if (!gmx_fexist(fn))
1854 gmx_fatal(FARGS, "%s does not exist", fn);
1857 t_state state;
1858 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1859 opts, ir, bZero, bGenVel, bVerbose, &state,
1860 &atypes, &sys, &mi, &intermolecular_interactions,
1861 interactions, &comb, &reppow, &fudgeQQ,
1862 opts->bMorse,
1863 wi);
1865 if (debug)
1867 pr_symtab(debug, 0, "After new_status", &sys.symtab);
1870 nvsite = 0;
1871 /* set parameters for virtual site construction (not for vsiten) */
1872 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1874 nvsite +=
1875 set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions);
1877 /* now throw away all obsolete bonds, angles and dihedrals: */
1878 /* note: constraints are ALWAYS removed */
1879 if (nvsite)
1881 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1883 clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds);
1887 if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == econtSHAKE))
1889 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1891 sprintf(warn_buf, "Can not do %s with %s, use %s",
1892 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1893 warning_error(wi, warn_buf);
1895 if (ir->bPeriodicMols)
1897 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1898 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1899 warning_error(wi, warn_buf);
1903 if (EI_SD (ir->eI) && ir->etc != etcNO)
1905 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1908 /* If we are doing QM/MM, check that we got the atom numbers */
1909 have_atomnumber = TRUE;
1910 for (gmx::index i = 0; i < gmx::ssize(atypes); i++)
1912 have_atomnumber = have_atomnumber && (atypes.atomNumberFromAtomType(i) >= 0);
1914 if (!have_atomnumber && ir->bQMMM)
1916 warning_error(wi,
1917 "\n"
1918 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1919 "field you are using does not contain atom numbers fields. This is an\n"
1920 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1921 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1922 "an integer just before the mass column in ffXXXnb.itp.\n"
1923 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1926 /* Check for errors in the input now, since they might cause problems
1927 * during processing further down.
1929 check_warning_error(wi, FARGS);
1931 if (nint_ftype(&sys, mi, F_POSRES) > 0 ||
1932 nint_ftype(&sys, mi, F_FBPOSRES) > 0)
1934 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1936 sprintf(warn_buf, "You are combining position restraints with %s pressure coupling, which can lead to instabilities. If you really want to combine position restraints with pressure coupling, we suggest to use %s pressure coupling instead.",
1937 EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
1938 warning_note(wi, warn_buf);
1941 const char *fn = opt2fn("-r", NFILE, fnm);
1942 const char *fnB;
1944 if (!gmx_fexist(fn))
1946 gmx_fatal(FARGS,
1947 "Cannot find position restraint file %s (option -r).\n"
1948 "From GROMACS-2018, you need to specify the position restraint "
1949 "coordinate files explicitly to avoid mistakes, although you can "
1950 "still use the same file as you specify for the -c option.", fn);
1953 if (opt2bSet("-rb", NFILE, fnm))
1955 fnB = opt2fn("-rb", NFILE, fnm);
1956 if (!gmx_fexist(fnB))
1958 gmx_fatal(FARGS,
1959 "Cannot find B-state position restraint file %s (option -rb).\n"
1960 "From GROMACS-2018, you need to specify the position restraint "
1961 "coordinate files explicitly to avoid mistakes, although you can "
1962 "still use the same file as you specify for the -c option.", fn);
1965 else
1967 fnB = fn;
1970 if (bVerbose)
1972 fprintf(stderr, "Reading position restraint coords from %s", fn);
1973 if (strcmp(fn, fnB) == 0)
1975 fprintf(stderr, "\n");
1977 else
1979 fprintf(stderr, " and %s\n", fnB);
1982 gen_posres(&sys, mi, fn, fnB,
1983 ir->refcoord_scaling, ir->ePBC,
1984 ir->posres_com, ir->posres_comB,
1985 wi);
1988 /* If we are using CMAP, setup the pre-interpolation grid */
1989 if (interactions[F_CMAP].ncmap() > 0)
1991 init_cmap_grid(&sys.ffparams.cmap_grid, interactions[F_CMAP].cmapAngles, interactions[F_CMAP].cmakeGridSpacing);
1992 setup_cmap(interactions[F_CMAP].cmakeGridSpacing, interactions[F_CMAP].cmapAngles, interactions[F_CMAP].cmap, &sys.ffparams.cmap_grid);
1995 set_wall_atomtype(&atypes, opts, ir, wi);
1996 if (bRenum)
1998 atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
2001 if (ir->implicit_solvent)
2003 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
2006 /* PELA: Copy the atomtype data to the topology atomtype list */
2007 atypes.copyTot_atomtypes(&(sys.atomtypes));
2009 if (debug)
2011 pr_symtab(debug, 0, "After atype.renumberTypes", &sys.symtab);
2014 if (bVerbose)
2016 fprintf(stderr, "converting bonded parameters...\n");
2019 const int ntype = atypes.size();
2020 convertInteractionsOfType(ntype, interactions, mi, intermolecular_interactions.get(),
2021 comb, reppow, fudgeQQ, &sys);
2023 if (debug)
2025 pr_symtab(debug, 0, "After converInteractionsOfType", &sys.symtab);
2028 /* set ptype to VSite for virtual sites */
2029 for (gmx_moltype_t &moltype : sys.moltype)
2031 set_vsites_ptype(FALSE, &moltype);
2033 if (debug)
2035 pr_symtab(debug, 0, "After virtual sites", &sys.symtab);
2037 /* Check velocity for virtual sites and shells */
2038 if (bGenVel)
2040 check_vel(&sys, state.v.rvec_array());
2043 /* check for shells and inpurecs */
2044 check_shells_inputrec(&sys, ir, wi);
2046 /* check masses */
2047 check_mol(&sys, wi);
2049 checkForUnboundAtoms(&sys, bVerbose, wi);
2051 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2053 check_bonds_timestep(&sys, ir->delta_t, wi);
2056 checkDecoupledModeAccuracy(&sys, ir, wi);
2058 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2060 warning_note(wi, "Zero-step energy minimization will alter the coordinates before calculating the energy. If you just want the energy of a single point, try zero-step MD (with unconstrained_start = yes). To do multiple single-point energy evaluations of different configurations of the same topology, use mdrun -rerun.");
2063 check_warning_error(wi, FARGS);
2065 if (bVerbose)
2067 fprintf(stderr, "initialising group options...\n");
2069 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
2070 &sys, bVerbose, mdModules.notifier(), ir, wi);
2072 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2074 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2076 real buffer_temp;
2078 if (EI_MD(ir->eI) && ir->etc == etcNO)
2080 if (bGenVel)
2082 buffer_temp = opts->tempi;
2084 else
2086 buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2088 if (buffer_temp > 0)
2090 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
2091 warning_note(wi, warn_buf);
2093 else
2095 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
2096 gmx::roundToInt(verlet_buffer_ratio_NVE_T0*100));
2097 warning_note(wi, warn_buf);
2100 else
2102 buffer_temp = get_max_reference_temp(ir, wi);
2105 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2107 /* NVE with initial T=0: we add a fixed ratio to rlist.
2108 * Since we don't actually use verletbuf_tol, we set it to -1
2109 * so it can't be misused later.
2111 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2112 ir->verletbuf_tol = -1;
2114 else
2116 /* We warn for NVE simulations with a drift tolerance that
2117 * might result in a 1(.1)% drift over the total run-time.
2118 * Note that we can't warn when nsteps=0, since we don't
2119 * know how many steps the user intends to run.
2121 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
2122 ir->nsteps > 0)
2124 const real driftTolerance = 0.01;
2125 /* We use 2 DOF per atom = 2kT pot+kin energy,
2126 * to be on the safe side with constraints.
2128 const real totalEnergyDriftPerAtomPerPicosecond = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
2130 if (ir->verletbuf_tol > 1.1*driftTolerance*totalEnergyDriftPerAtomPerPicosecond)
2132 sprintf(warn_buf, "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an NVE simulation of length %g ps, which can give a final drift of %d%%. For conserving energy to %d%% when using constraints, you might need to set verlet-buffer-tolerance to %.1e.",
2133 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
2134 gmx::roundToInt(ir->verletbuf_tol/totalEnergyDriftPerAtomPerPicosecond*100),
2135 gmx::roundToInt(100*driftTolerance),
2136 driftTolerance*totalEnergyDriftPerAtomPerPicosecond);
2137 warning_note(wi, warn_buf);
2141 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi);
2146 /* Init the temperature coupling state */
2147 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2149 if (debug)
2151 pr_symtab(debug, 0, "After index", &sys.symtab);
2154 triple_check(mdparin, ir, &sys, wi);
2155 close_symtab(&sys.symtab);
2156 if (debug)
2158 pr_symtab(debug, 0, "After close", &sys.symtab);
2161 /* make exclusions between QM atoms and remove charges if needed */
2162 if (ir->bQMMM)
2164 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_ORIGINAL);
2165 if (ir->QMMMscheme != eQMMMschemeoniom)
2167 std::vector<int> qmmmAtoms = qmmmAtomIndices(*ir, sys);
2168 removeQmmmAtomCharges(&sys, qmmmAtoms);
2172 if (ir->eI == eiMimic)
2174 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_MIMIC);
2177 if (ftp2bSet(efTRN, NFILE, fnm))
2179 if (bVerbose)
2181 fprintf(stderr, "getting data from old trajectory ...\n");
2183 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
2184 bNeedVel, bGenVel, fr_time, ir, &state, &sys, oenv);
2187 if (ir->ePBC == epbcXY && ir->nwall != 2)
2189 clear_rvec(state.box[ZZ]);
2192 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2194 /* Calculate the optimal grid dimensions */
2195 matrix scaledBox;
2196 EwaldBoxZScaler boxScaler(*ir);
2197 boxScaler.scaleBox(state.box, scaledBox);
2199 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2201 /* Mark fourier_spacing as not used */
2202 ir->fourier_spacing = 0;
2204 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2206 set_warning_line(wi, mdparin, -1);
2207 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2209 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2210 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize,
2211 &(ir->nkx), &(ir->nky), &(ir->nkz));
2212 if (ir->nkx < minGridSize ||
2213 ir->nky < minGridSize ||
2214 ir->nkz < minGridSize)
2216 warning_error(wi, "The PME grid size should be >= 2*(pme-order - 1); either manually increase the grid size or decrease pme-order");
2220 /* MRS: eventually figure out better logic for initializing the fep
2221 values that makes declaring the lambda and declaring the state not
2222 potentially conflict if not handled correctly. */
2223 if (ir->efep != efepNO)
2225 state.fep_state = ir->fepvals->init_fep_state;
2226 for (i = 0; i < efptNR; i++)
2228 /* init_lambda trumps state definitions*/
2229 if (ir->fepvals->init_lambda >= 0)
2231 state.lambda[i] = ir->fepvals->init_lambda;
2233 else
2235 if (ir->fepvals->all_lambda[i] == nullptr)
2237 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2239 else
2241 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2247 pull_t *pull = nullptr;
2249 if (ir->bPull)
2251 pull = set_pull_init(ir, &sys, state.x.rvec_array(), state.box, state.lambda[efptMASS], wi);
2254 /* Modules that supply external potential for pull coordinates
2255 * should register those potentials here. finish_pull() will check
2256 * that providers have been registerd for all external potentials.
2259 if (ir->bDoAwh)
2261 tensor compressibility = { { 0 } };
2262 if (ir->epc != epcNO)
2264 copy_mat(ir->compress, compressibility);
2266 setStateDependentAwhParams(ir->awhParams, ir->pull, pull,
2267 state.box, ir->ePBC, compressibility,
2268 &ir->opts, wi);
2271 if (ir->bPull)
2273 finish_pull(pull);
2276 if (ir->bRot)
2278 set_reference_positions(ir->rot, state.x.rvec_array(), state.box,
2279 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2280 wi);
2283 /* reset_multinr(sys); */
2285 if (EEL_PME(ir->coulombtype))
2287 float ratio = pme_load_estimate(sys, *ir, state.box);
2288 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2289 /* With free energy we might need to do PME both for the A and B state
2290 * charges. This will double the cost, but the optimal performance will
2291 * then probably be at a slightly larger cut-off and grid spacing.
2293 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2294 (ir->efep != efepNO && ratio > 2.0/3.0))
2296 warning_note(wi,
2297 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2298 "and for highly parallel simulations between 0.25 and 0.33,\n"
2299 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2300 if (ir->efep != efepNO)
2302 warning_note(wi,
2303 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2309 char warn_buf[STRLEN];
2310 double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2311 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2312 if (cio > 2000)
2314 set_warning_line(wi, mdparin, -1);
2315 warning_note(wi, warn_buf);
2317 else
2319 printf("%s\n", warn_buf);
2323 // Add the md modules internal parameters that are not mdp options
2324 // e.g., atom indices
2327 gmx::KeyValueTreeBuilder internalParameterBuilder;
2328 mdModules.notifier().notifier_.notify(internalParameterBuilder.rootObject());
2329 ir->internalParameters = std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2332 if (bVerbose)
2334 fprintf(stderr, "writing run input file...\n");
2337 done_warning(wi, FARGS);
2338 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, &sys);
2340 /* Output IMD group, if bIMD is TRUE */
2341 gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, &sys, NFILE, fnm);
2343 sfree(opts->define);
2344 sfree(opts->include);
2345 sfree(opts);
2346 for (auto &mol : mi)
2348 // Some of the contents of molinfo have been stolen, so
2349 // fullCleanUp can't be called.
2350 mol.partialCleanUp();
2352 done_inputrec_strings();
2353 output_env_done(oenv);
2355 return 0;