Make PBC type enumeration into PbcType enum class
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.cpp
blob28562ffeb4904b5d305e5decd32a893e363fb8c0
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 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.
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.
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.
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.
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.
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.
38 #include "gmxpre.h"
40 #include "grompp.h"
42 #include <cerrno>
43 #include <climits>
44 #include <cmath>
45 #include <cstring>
47 #include <algorithm>
48 #include <memory>
49 #include <vector>
51 #include <sys/types.h>
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/futil.h"
104 #include "gromacs/utility/gmxassert.h"
105 #include "gromacs/utility/keyvaluetreebuilder.h"
106 #include "gromacs/utility/listoflists.h"
107 #include "gromacs/utility/mdmodulenotification.h"
108 #include "gromacs/utility/smalloc.h"
109 #include "gromacs/utility/snprintf.h"
111 /* TODO The implementation details should move to their own source file. */
112 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int> atoms,
113 gmx::ArrayRef<const real> params,
114 const std::string& name) :
115 atoms_(atoms.begin(), atoms.end()),
116 interactionTypeName_(name)
118 GMX_RELEASE_ASSERT(
119 params.size() <= forceParam_.size(),
120 gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM)
121 .c_str());
122 auto forceParamIt = forceParam_.begin();
123 for (const auto param : params)
125 *forceParamIt++ = param;
127 while (forceParamIt != forceParam_.end())
129 *forceParamIt++ = NOTSET;
133 const int& InteractionOfType::ai() const
135 GMX_RELEASE_ASSERT(!atoms_.empty(), "Need to have at least one atom set");
136 return atoms_[0];
139 const int& InteractionOfType::aj() const
141 GMX_RELEASE_ASSERT(atoms_.size() > 1, "Need to have at least two atoms set");
142 return atoms_[1];
145 const int& InteractionOfType::ak() const
147 GMX_RELEASE_ASSERT(atoms_.size() > 2, "Need to have at least three atoms set");
148 return atoms_[2];
151 const int& InteractionOfType::al() const
153 GMX_RELEASE_ASSERT(atoms_.size() > 3, "Need to have at least four atoms set");
154 return atoms_[3];
157 const int& InteractionOfType::am() const
159 GMX_RELEASE_ASSERT(atoms_.size() > 4, "Need to have at least five atoms set");
160 return atoms_[4];
163 const real& InteractionOfType::c0() const
165 return forceParam_[0];
168 const real& InteractionOfType::c1() const
170 return forceParam_[1];
173 const real& InteractionOfType::c2() const
175 return forceParam_[2];
178 const std::string& InteractionOfType::interactionTypeName() const
180 return interactionTypeName_;
183 void InteractionOfType::sortBondAtomIds()
185 if (aj() < ai())
187 std::swap(atoms_[0], atoms_[1]);
191 void InteractionOfType::sortAngleAtomIds()
193 if (ak() < ai())
195 std::swap(atoms_[0], atoms_[2]);
199 void InteractionOfType::sortDihedralAtomIds()
201 if (al() < ai())
203 std::swap(atoms_[0], atoms_[3]);
204 std::swap(atoms_[1], atoms_[2]);
208 void InteractionOfType::sortAtomIds()
210 if (isBond())
212 sortBondAtomIds();
214 else if (isAngle())
216 sortAngleAtomIds();
218 else if (isDihedral())
220 sortDihedralAtomIds();
222 else
224 GMX_THROW(gmx::InternalError(
225 "Can not sort parameters other than bonds, angles or dihedrals"));
229 void InteractionOfType::setForceParameter(int pos, real value)
231 GMX_RELEASE_ASSERT(pos < MAXFORCEPARAM,
232 "Can't set parameter beyond the maximum number of parameters");
233 forceParam_[pos] = value;
236 void MoleculeInformation::initMolInfo()
238 init_block(&mols);
239 excls.clear();
240 init_t_atoms(&atoms, 0, FALSE);
243 void MoleculeInformation::partialCleanUp()
245 done_block(&mols);
248 void MoleculeInformation::fullCleanUp()
250 done_atom(&atoms);
251 done_block(&mols);
254 static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
256 int n = 0;
257 /* For all the molecule types */
258 for (auto& mol : mols)
260 n += mol.interactions[ifunc].size();
261 mol.interactions[ifunc].interactionTypes.clear();
263 return n;
266 static int check_atom_names(const char* fn1, const char* fn2, gmx_mtop_t* mtop, const t_atoms* at)
268 int m, i, j, nmismatch;
269 t_atoms* tat;
270 #define MAXMISMATCH 20
272 if (mtop->natoms != at->nr)
274 gmx_incons("comparing atom names");
277 nmismatch = 0;
278 i = 0;
279 for (const gmx_molblock_t& molb : mtop->molblock)
281 tat = &mtop->moltype[molb.type].atoms;
282 for (m = 0; m < molb.nmol; m++)
284 for (j = 0; j < tat->nr; j++)
286 if (strcmp(*(tat->atomname[j]), *(at->atomname[i])) != 0)
288 if (nmismatch < MAXMISMATCH)
290 fprintf(stderr,
291 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
292 i + 1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
294 else if (nmismatch == MAXMISMATCH)
296 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
298 nmismatch++;
300 i++;
305 return nmismatch;
308 static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
310 /* This check is not intended to ensure accurate integration,
311 * rather it is to signal mistakes in the mdp settings.
312 * A common mistake is to forget to turn on constraints
313 * for MD after energy minimization with flexible bonds.
314 * This check can also detect too large time steps for flexible water
315 * models, but such errors will often be masked by the constraints
316 * mdp options, which turns flexible water into water with bond constraints,
317 * but without an angle constraint. Unfortunately such incorrect use
318 * of water models can not easily be detected without checking
319 * for specific model names.
321 * The stability limit of leap-frog or velocity verlet is 4.44 steps
322 * per oscillational period.
323 * But accurate bonds distributions are lost far before that limit.
324 * To allow relatively common schemes (although not common with Gromacs)
325 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
326 * we set the note limit to 10.
328 int min_steps_warn = 5;
329 int min_steps_note = 10;
330 int ftype;
331 int i, a1, a2, w_a1, w_a2, j;
332 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
333 bool bFound, bWater, bWarn;
334 char warn_buf[STRLEN];
336 /* Get the interaction parameters */
337 gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
339 twopi2 = gmx::square(2 * M_PI);
341 limit2 = gmx::square(min_steps_note * dt);
343 w_a1 = w_a2 = -1;
344 w_period2 = -1.0;
346 const gmx_moltype_t* w_moltype = nullptr;
347 for (const gmx_moltype_t& moltype : mtop->moltype)
349 const t_atom* atom = moltype.atoms.atom;
350 const InteractionLists& ilist = moltype.ilist;
351 const InteractionList& ilc = ilist[F_CONSTR];
352 const InteractionList& ils = ilist[F_SETTLE];
353 for (ftype = 0; ftype < F_NRE; ftype++)
355 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
357 continue;
360 const InteractionList& ilb = ilist[ftype];
361 for (i = 0; i < ilb.size(); i += 3)
363 fc = ip[ilb.iatoms[i]].harmonic.krA;
364 re = ip[ilb.iatoms[i]].harmonic.rA;
365 if (ftype == F_G96BONDS)
367 /* Convert squared sqaure fc to harmonic fc */
368 fc = 2 * fc * re;
370 a1 = ilb.iatoms[i + 1];
371 a2 = ilb.iatoms[i + 2];
372 m1 = atom[a1].m;
373 m2 = atom[a2].m;
374 if (fc > 0 && m1 > 0 && m2 > 0)
376 period2 = twopi2 * m1 * m2 / ((m1 + m2) * fc);
378 else
380 period2 = GMX_FLOAT_MAX;
382 if (debug)
384 fprintf(debug, "fc %g m1 %g m2 %g period %g\n", fc, m1, m2, std::sqrt(period2));
386 if (period2 < limit2)
388 bFound = false;
389 for (j = 0; j < ilc.size(); j += 3)
391 if ((ilc.iatoms[j + 1] == a1 && ilc.iatoms[j + 2] == a2)
392 || (ilc.iatoms[j + 1] == a2 && ilc.iatoms[j + 2] == a1))
394 bFound = true;
397 for (j = 0; j < ils.size(); j += 4)
399 if ((a1 == ils.iatoms[j + 1] || a1 == ils.iatoms[j + 2] || a1 == ils.iatoms[j + 3])
400 && (a2 == ils.iatoms[j + 1] || a2 == ils.iatoms[j + 2]
401 || a2 == ils.iatoms[j + 3]))
403 bFound = true;
406 if (!bFound && (w_moltype == nullptr || period2 < w_period2))
408 w_moltype = &moltype;
409 w_a1 = a1;
410 w_a2 = a2;
411 w_period2 = period2;
418 if (w_moltype != nullptr)
420 bWarn = (w_period2 < gmx::square(min_steps_warn * dt));
421 /* A check that would recognize most water models */
422 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' && w_moltype->atoms.nr <= 5);
423 sprintf(warn_buf,
424 "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated "
425 "oscillational period of %.1e ps, which is less than %d times the time step of "
426 "%.1e ps.\n"
427 "%s",
428 *w_moltype->name, w_a1 + 1, *w_moltype->atoms.atomname[w_a1], w_a2 + 1,
429 *w_moltype->atoms.atomname[w_a2], std::sqrt(w_period2),
430 bWarn ? min_steps_warn : min_steps_note, dt,
431 bWater ? "Maybe you asked for fexible water."
432 : "Maybe you forgot to change the constraints mdp option.");
433 if (bWarn)
435 warning(wi, warn_buf);
437 else
439 warning_note(wi, warn_buf);
444 static void check_vel(gmx_mtop_t* mtop, rvec v[])
446 for (const AtomProxy atomP : AtomRange(*mtop))
448 const t_atom& local = atomP.atom();
449 int i = atomP.globalAtomNumber();
450 if (local.ptype == eptShell || local.ptype == eptBond || local.ptype == eptVSite)
452 clear_rvec(v[i]);
457 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, 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 || local.ptype == eptBond)
467 nshells++;
470 if ((nshells > 0) && (ir->nstcalcenergy != 1))
472 set_warning_line(wi, "unknown", -1);
473 snprintf(warn_buf, STRLEN, "There are %d shells, changing nstcalcenergy from %d to 1",
474 nshells, ir->nstcalcenergy);
475 ir->nstcalcenergy = 1;
476 warning(wi, warn_buf);
480 /* TODO Decide whether this function can be consolidated with
481 * gmx_mtop_ftype_count */
482 static int nint_ftype(gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
484 int nint = 0;
485 for (const gmx_molblock_t& molb : mtop->molblock)
487 nint += molb.nmol * mi[molb.type].interactions[ftype].size();
490 return nint;
493 /* This routine reorders the molecule type array
494 * in the order of use in the molblocks,
495 * unused molecule types are deleted.
497 static void renumber_moltypes(gmx_mtop_t* sys, std::vector<MoleculeInformation>* molinfo)
500 std::vector<int> order;
501 for (gmx_molblock_t& molblock : sys->molblock)
503 const auto found = std::find_if(order.begin(), order.end(), [&molblock](const auto& entry) {
504 return molblock.type == entry;
506 if (found == order.end())
508 /* This type did not occur yet, add it */
509 order.push_back(molblock.type);
510 molblock.type = order.size() - 1;
512 else
514 molblock.type = std::distance(order.begin(), found);
518 /* We still need to reorder the molinfo structs */
519 std::vector<MoleculeInformation> minew(order.size());
520 int index = 0;
521 for (auto& mi : *molinfo)
523 const auto found = std::find(order.begin(), order.end(), index);
524 if (found != order.end())
526 int pos = std::distance(order.begin(), found);
527 minew[pos] = mi;
529 else
531 // Need to manually clean up memory ....
532 mi.fullCleanUp();
534 index++;
537 *molinfo = minew;
540 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t* mtop)
542 mtop->moltype.resize(mi.size());
543 int pos = 0;
544 for (const auto& mol : mi)
546 gmx_moltype_t& molt = mtop->moltype[pos];
547 molt.name = mol.name;
548 molt.atoms = mol.atoms;
549 /* ilists are copied later */
550 molt.excls = mol.excls;
551 pos++;
555 static void new_status(const char* topfile,
556 const char* topppfile,
557 const char* confin,
558 t_gromppopts* opts,
559 t_inputrec* ir,
560 gmx_bool bZero,
561 bool bGenVel,
562 bool bVerbose,
563 t_state* state,
564 PreprocessingAtomTypes* atypes,
565 gmx_mtop_t* sys,
566 std::vector<MoleculeInformation>* mi,
567 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
568 gmx::ArrayRef<InteractionsOfType> interactions,
569 int* comb,
570 double* reppow,
571 real* fudgeQQ,
572 gmx_bool bMorse,
573 warninp* wi)
575 std::vector<gmx_molblock_t> molblock;
576 int i, nmismatch;
577 bool ffParametrizedWithHBondConstraints;
578 char buf[STRLEN];
579 char warn_buf[STRLEN];
581 /* TOPOLOGY processing */
582 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab), interactions,
583 comb, reppow, fudgeQQ, atypes, mi, intermolecular_interactions, ir,
584 &molblock, &ffParametrizedWithHBondConstraints, wi);
586 sys->molblock.clear();
588 sys->natoms = 0;
589 for (const gmx_molblock_t& molb : molblock)
591 if (!sys->molblock.empty() && 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 && ffParametrizedWithHBondConstraints && ir->delta_t < 0.0026
646 && gmx_mtop_ftype_count(sys, F_VSITE3FD) == 0)
648 set_warning_line(wi, "unknown", -1);
649 warning_note(wi,
650 "You are using constraints on all bonds, whereas the forcefield "
651 "has been parametrized only with constraints involving hydrogen atoms. "
652 "We suggest using constraints = h-bonds instead, this will also improve "
653 "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,
671 "number of coordinates in coordinate file (%s, %d)\n"
672 " does not match topology (%s, %d)",
673 confin, state->natoms, topfile, sys->natoms);
675 /* It would be nice to get rid of the copies below, but we don't know
676 * a priori if the number of atoms in confin matches what we expect.
678 state->flags |= (1 << estX);
679 if (v != nullptr)
681 state->flags |= (1 << estV);
683 state_change_natoms(state, state->natoms);
684 std::copy(x, x + state->natoms, state->x.data());
685 sfree(x);
686 if (v != nullptr)
688 std::copy(v, v + state->natoms, state->v.data());
689 sfree(v);
691 /* This call fixes the box shape for runs with pressure scaling */
692 set_box_rel(ir, state);
694 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
695 done_top(conftop);
696 sfree(conftop);
698 if (nmismatch)
700 sprintf(buf,
701 "%d non-matching atom name%s\n"
702 "atom names from %s will be used\n"
703 "atom names from %s will be ignored\n",
704 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
705 warning(wi, buf);
708 /* Do more checks, mostly related to constraints */
709 if (bVerbose)
711 fprintf(stderr, "double-checking input for internal consistency...\n");
714 bool bHasNormalConstraints =
715 0 < (nint_ftype(sys, *mi, F_CONSTR) + nint_ftype(sys, *mi, F_CONSTRNC));
716 bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
717 double_check(ir, state->box, bHasNormalConstraints, bHasAnyConstraints, wi);
720 if (bGenVel)
722 real* mass;
724 snew(mass, state->natoms);
725 for (const AtomProxy atomP : AtomRange(*sys))
727 const t_atom& local = atomP.atom();
728 int i = atomP.globalAtomNumber();
729 mass[i] = local.m;
732 if (opts->seed == -1)
734 opts->seed = static_cast<int>(gmx::makeRandomSeed());
735 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
737 state->flags |= (1 << estV);
738 maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array());
740 stop_cm(stdout, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
741 sfree(mass);
745 static void copy_state(const char* slog, t_trxframe* fr, bool bReadVel, t_state* state, double* use_time)
747 if (fr->not_ok & FRAME_NOT_OK)
749 gmx_fatal(FARGS, "Can not start from an incomplete frame");
751 if (!fr->bX)
753 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s", slog);
756 std::copy(fr->x, fr->x + state->natoms, state->x.data());
757 if (bReadVel)
759 if (!fr->bV)
761 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
763 std::copy(fr->v, fr->v + state->natoms, state->v.data());
765 if (fr->bBox)
767 copy_mat(fr->box, state->box);
770 *use_time = fr->time;
773 static void cont_status(const char* slog,
774 const char* ener,
775 bool bNeedVel,
776 bool bGenVel,
777 real fr_time,
778 t_inputrec* ir,
779 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, "Reading Coordinates%s and Box size from old trajectory\n",
792 bReadVel ? ", Velocities" : "");
793 if (fr_time == -1)
795 fprintf(stderr, "Will read whole trajectory\n");
797 else
799 fprintf(stderr, "Will read till time %g\n", fr_time);
801 if (!bReadVel)
803 if (bGenVel)
805 fprintf(stderr,
806 "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",
821 slog);
822 for (auto& vi : makeArrayRef(state->v))
824 vi = { 0, 0, 0 };
826 close_trx(fp);
827 /* Search for a frame without velocities */
828 bReadVel = false;
829 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
833 state->natoms = fr.natoms;
835 if (sys->natoms != state->natoms)
837 gmx_fatal(FARGS,
838 "Number of atoms in Topology "
839 "is not the same as in Trajectory");
841 copy_state(slog, &fr, bReadVel, state, &use_time);
843 /* Find the appropriate frame */
844 while ((fr_time == -1 || fr.time < fr_time) && read_next_frame(oenv, fp, &fr))
846 copy_state(slog, &fr, bReadVel, state, &use_time);
849 close_trx(fp);
851 /* Set the relative box lengths for preserving the box shape.
852 * Note that this call can lead to differences in the last bit
853 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
855 set_box_rel(ir, state);
857 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
858 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
860 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
862 get_enx_state(ener, use_time, sys->groups, ir, state);
863 preserve_box_shape(ir, state->box_rel, state->boxv);
867 static void read_posres(gmx_mtop_t* mtop,
868 gmx::ArrayRef<const MoleculeInformation> molinfo,
869 gmx_bool bTopB,
870 const char* fn,
871 int rc_scaling,
872 PbcType pbcType,
873 rvec com,
874 warninp* wi)
876 gmx_bool* hadAtom;
877 rvec * x, *v;
878 dvec sum;
879 double totmass;
880 t_topology* top;
881 matrix box, invbox;
882 int natoms, npbcdim = 0;
883 char warn_buf[STRLEN];
884 int a, nat_molb;
885 t_atom* atom;
887 snew(top, 1);
888 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
889 natoms = top->atoms.nr;
890 done_top(top);
891 sfree(top);
892 if (natoms != mtop->natoms)
894 sprintf(warn_buf,
895 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
896 "(%d). Will assume that the first %d atoms in the topology and %s match.",
897 fn, natoms, mtop->natoms, std::min(mtop->natoms, natoms), fn);
898 warning(wi, warn_buf);
901 npbcdim = numPbcDimensions(pbcType);
902 GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
903 clear_rvec(com);
904 if (rc_scaling != erscNO)
906 copy_mat(box, invbox);
907 for (int j = npbcdim; j < DIM; j++)
909 clear_rvec(invbox[j]);
910 invbox[j][j] = 1;
912 gmx::invertBoxMatrix(invbox, invbox);
915 /* Copy the reference coordinates to mtop */
916 clear_dvec(sum);
917 totmass = 0;
918 a = 0;
919 snew(hadAtom, natoms);
920 for (gmx_molblock_t& molb : mtop->molblock)
922 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
923 const InteractionsOfType* pr = &(molinfo[molb.type].interactions[F_POSRES]);
924 const InteractionsOfType* prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
925 if (pr->size() > 0 || prfb->size() > 0)
927 atom = mtop->moltype[molb.type].atoms.atom;
928 for (const auto& restraint : pr->interactionTypes)
930 int ai = restraint.ai();
931 if (ai >= natoms)
933 gmx_fatal(FARGS,
934 "Position restraint atom index (%d) in moltype '%s' is larger than "
935 "number of atoms in %s (%d).\n",
936 ai + 1, *molinfo[molb.type].name, fn, natoms);
938 hadAtom[ai] = TRUE;
939 if (rc_scaling == erscCOM)
941 /* Determine the center of mass of the posres reference coordinates */
942 for (int j = 0; j < npbcdim; j++)
944 sum[j] += atom[ai].m * x[a + ai][j];
946 totmass += atom[ai].m;
949 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
950 for (const auto& restraint : prfb->interactionTypes)
952 int ai = restraint.ai();
953 if (ai >= natoms)
955 gmx_fatal(FARGS,
956 "Position restraint atom index (%d) in moltype '%s' is larger than "
957 "number of atoms in %s (%d).\n",
958 ai + 1, *molinfo[molb.type].name, fn, natoms);
960 if (rc_scaling == erscCOM && !hadAtom[ai])
962 /* Determine the center of mass of the posres reference coordinates */
963 for (int j = 0; j < npbcdim; j++)
965 sum[j] += atom[ai].m * x[a + ai][j];
967 totmass += atom[ai].m;
970 if (!bTopB)
972 molb.posres_xA.resize(nat_molb);
973 for (int i = 0; i < nat_molb; i++)
975 copy_rvec(x[a + i], molb.posres_xA[i]);
978 else
980 molb.posres_xB.resize(nat_molb);
981 for (int i = 0; i < nat_molb; i++)
983 copy_rvec(x[a + i], molb.posres_xB[i]);
987 a += nat_molb;
989 if (rc_scaling == erscCOM)
991 if (totmass == 0)
993 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
995 for (int j = 0; j < npbcdim; j++)
997 com[j] = sum[j] / totmass;
999 fprintf(stderr,
1000 "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n",
1001 com[XX], com[YY], com[ZZ]);
1004 if (rc_scaling != erscNO)
1006 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
1008 for (gmx_molblock_t& molb : mtop->molblock)
1010 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
1011 if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1013 std::vector<gmx::RVec>& xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1014 for (int i = 0; i < nat_molb; i++)
1016 for (int j = 0; j < npbcdim; j++)
1018 if (rc_scaling == erscALL)
1020 /* Convert from Cartesian to crystal coordinates */
1021 xp[i][j] *= invbox[j][j];
1022 for (int k = j + 1; k < npbcdim; k++)
1024 xp[i][j] += invbox[k][j] * xp[i][k];
1027 else if (rc_scaling == erscCOM)
1029 /* Subtract the center of mass */
1030 xp[i][j] -= com[j];
1037 if (rc_scaling == erscCOM)
1039 /* Convert the COM from Cartesian to crystal coordinates */
1040 for (int j = 0; j < npbcdim; j++)
1042 com[j] *= invbox[j][j];
1043 for (int k = j + 1; k < npbcdim; k++)
1045 com[j] += invbox[k][j] * com[k];
1051 sfree(x);
1052 sfree(v);
1053 sfree(hadAtom);
1056 static void gen_posres(gmx_mtop_t* mtop,
1057 gmx::ArrayRef<const MoleculeInformation> mi,
1058 const char* fnA,
1059 const char* fnB,
1060 int rc_scaling,
1061 PbcType pbcType,
1062 rvec com,
1063 rvec comB,
1064 warninp* wi)
1066 read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi);
1067 /* It is safer to simply read the b-state posres rather than trying
1068 * to be smart and copy the positions.
1070 read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi);
1073 static void set_wall_atomtype(PreprocessingAtomTypes* at, t_gromppopts* opts, t_inputrec* ir, warninp* wi)
1075 int i;
1076 char warn_buf[STRLEN];
1078 if (ir->nwall > 0)
1080 fprintf(stderr, "Searching the wall atom type(s)\n");
1082 for (i = 0; i < ir->nwall; i++)
1084 ir->wall_atomtype[i] = at->atomTypeFromName(opts->wall_atomtype[i]);
1085 if (ir->wall_atomtype[i] == NOTSET)
1087 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1088 warning_error(wi, warn_buf);
1093 static int nrdf_internal(const t_atoms* atoms)
1095 int i, nmass, nrdf;
1097 nmass = 0;
1098 for (i = 0; i < atoms->nr; i++)
1100 /* Vsite ptype might not be set here yet, so also check the mass */
1101 if ((atoms->atom[i].ptype == eptAtom || atoms->atom[i].ptype == eptNucleus) && atoms->atom[i].m > 0)
1103 nmass++;
1106 switch (nmass)
1108 case 0: nrdf = 0; break;
1109 case 1: nrdf = 0; break;
1110 case 2: nrdf = 1; break;
1111 default: nrdf = nmass * 3 - 6; break;
1114 return nrdf;
1117 static void spline1d(double dx, const double* y, int n, double* u, double* y2)
1119 int i;
1120 double p, q;
1122 y2[0] = 0.0;
1123 u[0] = 0.0;
1125 for (i = 1; i < n - 1; i++)
1127 p = 0.5 * y2[i - 1] + 2.0;
1128 y2[i] = -0.5 / p;
1129 q = (y[i + 1] - 2.0 * y[i] + y[i - 1]) / dx;
1130 u[i] = (3.0 * q / dx - 0.5 * u[i - 1]) / p;
1133 y2[n - 1] = 0.0;
1135 for (i = n - 2; i >= 0; i--)
1137 y2[i] = y2[i] * y2[i + 1] + u[i];
1142 static void
1143 interpolate1d(double xmin, double dx, const double* ya, const double* y2a, double x, double* y, 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]
1154 + ((a * a * a - a) * y2a[ix] + (b * b * b - b) * y2a[ix + 1]) * (dx * dx) / 6.0;
1155 *y1 = (ya[ix + 1] - ya[ix]) / dx - (3.0 * a * a - 1.0) / 6.0 * dx * y2a[ix]
1156 + (3.0 * b * b - 1.0) / 6.0 * dx * y2a[ix + 1];
1160 static void setup_cmap(int grid_spacing, int nc, gmx::ArrayRef<const real> grid, gmx_cmap_t* cmap_grid)
1162 int i, j, k, ii, jj, kk, idx;
1163 int offset;
1164 double dx, xmin, v, v1, v2, v12;
1165 double phi, psi;
1167 std::vector<double> tmp_u(2 * grid_spacing, 0.0);
1168 std::vector<double> tmp_u2(2 * grid_spacing, 0.0);
1169 std::vector<double> tmp_yy(2 * grid_spacing, 0.0);
1170 std::vector<double> tmp_y1(2 * grid_spacing, 0.0);
1171 std::vector<double> tmp_t2(2 * grid_spacing * 2 * grid_spacing, 0.0);
1172 std::vector<double> tmp_grid(2 * grid_spacing * 2 * grid_spacing, 0.0);
1174 dx = 360.0 / grid_spacing;
1175 xmin = -180.0 - dx * grid_spacing / 2;
1177 for (kk = 0; kk < nc; kk++)
1179 /* Compute an offset depending on which cmap we are using
1180 * Offset will be the map number multiplied with the
1181 * grid_spacing * grid_spacing * 2
1183 offset = kk * grid_spacing * grid_spacing * 2;
1185 for (i = 0; i < 2 * grid_spacing; i++)
1187 ii = (i + grid_spacing - grid_spacing / 2) % grid_spacing;
1189 for (j = 0; j < 2 * grid_spacing; j++)
1191 jj = (j + grid_spacing - grid_spacing / 2) % grid_spacing;
1192 tmp_grid[i * grid_spacing * 2 + j] = grid[offset + ii * grid_spacing + jj];
1196 for (i = 0; i < 2 * grid_spacing; i++)
1198 spline1d(dx, &(tmp_grid[2 * grid_spacing * i]), 2 * grid_spacing, tmp_u.data(),
1199 &(tmp_t2[2 * grid_spacing * i]));
1202 for (i = grid_spacing / 2; i < grid_spacing + grid_spacing / 2; i++)
1204 ii = i - grid_spacing / 2;
1205 phi = ii * dx - 180.0;
1207 for (j = grid_spacing / 2; j < grid_spacing + grid_spacing / 2; j++)
1209 jj = j - grid_spacing / 2;
1210 psi = jj * dx - 180.0;
1212 for (k = 0; k < 2 * grid_spacing; k++)
1214 interpolate1d(xmin, dx, &(tmp_grid[2 * grid_spacing * k]),
1215 &(tmp_t2[2 * grid_spacing * k]), psi, &tmp_yy[k], &tmp_y1[k]);
1218 spline1d(dx, tmp_yy.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1219 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1220 spline1d(dx, tmp_y1.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1221 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1223 idx = ii * grid_spacing + jj;
1224 cmap_grid->cmapdata[kk].cmap[idx * 4] = grid[offset + ii * grid_spacing + jj];
1225 cmap_grid->cmapdata[kk].cmap[idx * 4 + 1] = v1;
1226 cmap_grid->cmapdata[kk].cmap[idx * 4 + 2] = v2;
1227 cmap_grid->cmapdata[kk].cmap[idx * 4 + 3] = v12;
1233 static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
1235 int i, nelem;
1237 cmap_grid->grid_spacing = grid_spacing;
1238 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
1240 cmap_grid->cmapdata.resize(ngrid);
1242 for (i = 0; i < ngrid; i++)
1244 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
1249 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
1251 int count, count_mol;
1252 char buf[STRLEN];
1254 count = 0;
1255 for (const gmx_molblock_t& molb : mtop->molblock)
1257 count_mol = 0;
1258 gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1260 for (int i = 0; i < F_NRE; i++)
1262 if (i == F_SETTLE)
1264 count_mol += 3 * interactions[i].size();
1266 else if (interaction_function[i].flags & IF_CONSTRAINT)
1268 count_mol += interactions[i].size();
1272 if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1274 sprintf(buf,
1275 "Molecule type '%s' has %d constraints.\n"
1276 "For stability and efficiency there should not be more constraints than "
1277 "internal number of degrees of freedom: %d.\n",
1278 *mi[molb.type].name, count_mol, nrdf_internal(&mi[molb.type].atoms));
1279 warning(wi, buf);
1281 count += molb.nmol * count_mol;
1284 return count;
1287 static real calc_temp(const gmx_mtop_t* mtop, const t_inputrec* ir, rvec* v)
1289 double sum_mv2 = 0;
1290 for (const AtomProxy atomP : AtomRange(*mtop))
1292 const t_atom& local = atomP.atom();
1293 int i = atomP.globalAtomNumber();
1294 sum_mv2 += local.m * norm2(v[i]);
1297 double nrdf = 0;
1298 for (int g = 0; g < ir->opts.ngtc; g++)
1300 nrdf += ir->opts.nrdf[g];
1303 return sum_mv2 / (nrdf * BOLTZ);
1306 static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
1308 real ref_t;
1309 int i;
1310 bool bNoCoupl;
1312 ref_t = 0;
1313 bNoCoupl = false;
1314 for (i = 0; i < ir->opts.ngtc; i++)
1316 if (ir->opts.tau_t[i] < 0)
1318 bNoCoupl = true;
1320 else
1322 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1326 if (bNoCoupl)
1328 char buf[STRLEN];
1330 sprintf(buf,
1331 "Some temperature coupling groups do not use temperature coupling. We will assume "
1332 "their temperature is not more than %.3f K. If their temperature is higher, the "
1333 "energy error and the Verlet buffer might be underestimated.",
1334 ref_t);
1335 warning(wi, buf);
1338 return ref_t;
1341 /* Checks if there are unbound atoms in moleculetype molt.
1342 * Prints a note for each unbound atoms and a warning if any is present.
1344 static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi)
1346 const t_atoms* atoms = &molt->atoms;
1348 if (atoms->nr == 1)
1350 /* Only one atom, there can't be unbound atoms */
1351 return;
1354 std::vector<int> count(atoms->nr, 0);
1356 for (int ftype = 0; ftype < F_NRE; ftype++)
1358 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS)
1359 || (interaction_function[ftype].flags & IF_CONSTRAINT) || ftype == F_SETTLE)
1361 const InteractionList& il = molt->ilist[ftype];
1362 const int nral = NRAL(ftype);
1364 for (int i = 0; i < il.size(); i += 1 + nral)
1366 for (int j = 0; j < nral; j++)
1368 count[il.iatoms[i + 1 + j]]++;
1374 int numDanglingAtoms = 0;
1375 for (int a = 0; a < atoms->nr; a++)
1377 if (atoms->atom[a].ptype != eptVSite && count[a] == 0)
1379 if (bVerbose)
1381 fprintf(stderr,
1382 "\nAtom %d '%s' in moleculetype '%s' is not bound by a potential or "
1383 "constraint to any other atom in the same moleculetype.\n",
1384 a + 1, *atoms->atomname[a], *molt->name);
1386 numDanglingAtoms++;
1390 if (numDanglingAtoms > 0)
1392 char buf[STRLEN];
1393 sprintf(buf,
1394 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
1395 "other atom in the same moleculetype. Although technically this might not cause "
1396 "issues in a simulation, this often means that the user forgot to add a "
1397 "bond/potential/constraint or put multiple molecules in the same moleculetype "
1398 "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, gmx_bool bVerbose, warninp* wi)
1407 for (const gmx_moltype_t& molt : mtop->moltype)
1409 checkForUnboundAtoms(&molt, bVerbose, wi);
1413 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1415 * The specific decoupled modes this routine check for are angle modes
1416 * where the two bonds are constrained and the atoms a both ends are only
1417 * involved in a single constraint; the mass of the two atoms needs to
1418 * differ by more than \p massFactorThreshold.
1420 static bool haveDecoupledModeInMol(const gmx_moltype_t& molt,
1421 gmx::ArrayRef<const t_iparams> iparams,
1422 real massFactorThreshold)
1424 if (molt.ilist[F_CONSTR].size() == 0 && molt.ilist[F_CONSTRNC].size() == 0)
1426 return false;
1429 const t_atom* atom = molt.atoms.atom;
1431 const auto atomToConstraints =
1432 gmx::make_at2con(molt, iparams, gmx::FlexibleConstraintTreatment::Exclude);
1434 bool haveDecoupledMode = false;
1435 for (int ftype = 0; ftype < F_NRE; ftype++)
1437 if (interaction_function[ftype].flags & IF_ATYPE)
1439 const int nral = NRAL(ftype);
1440 const InteractionList& il = molt.ilist[ftype];
1441 for (int i = 0; i < il.size(); i += 1 + nral)
1443 /* Here we check for the mass difference between the atoms
1444 * at both ends of the angle, that the atoms at the ends have
1445 * 1 contraint and the atom in the middle at least 3; we check
1446 * that the 3 atoms are linked by constraints below.
1447 * We check for at least three constraints for the middle atom,
1448 * since with only the two bonds in the angle, we have 3-atom
1449 * molecule, which has much more thermal exhange in this single
1450 * angle mode than molecules with more atoms.
1451 * Note that this check also catches molecules where more atoms
1452 * are connected to one or more atoms in the angle, but by
1453 * bond potentials instead of angles. But such cases will not
1454 * occur in "normal" molecules and it doesn't hurt running
1455 * those with higher accuracy settings as well.
1457 int a0 = il.iatoms[1 + i];
1458 int a1 = il.iatoms[1 + i + 1];
1459 int a2 = il.iatoms[1 + i + 2];
1460 if ((atom[a0].m > atom[a2].m * massFactorThreshold || atom[a2].m > atom[a0].m * massFactorThreshold)
1461 && atomToConstraints[a0].ssize() == 1 && atomToConstraints[a2].ssize() == 1
1462 && atomToConstraints[a1].ssize() >= 3)
1464 int constraint0 = atomToConstraints[a0][0];
1465 int constraint2 = atomToConstraints[a2][0];
1467 bool foundAtom0 = false;
1468 bool foundAtom2 = false;
1469 for (const int constraint : atomToConstraints[a1])
1471 if (constraint == constraint0)
1473 foundAtom0 = true;
1475 if (constraint == constraint2)
1477 foundAtom2 = true;
1480 if (foundAtom0 && foundAtom2)
1482 haveDecoupledMode = true;
1489 return haveDecoupledMode;
1492 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1494 * When decoupled modes are present and the accuray in insufficient,
1495 * this routine issues a warning if the accuracy is insufficient.
1497 static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec* ir, warninp* wi)
1499 /* We only have issues with decoupled modes with normal MD.
1500 * With stochastic dynamics equipartitioning is enforced strongly.
1502 if (!EI_MD(ir->eI))
1504 return;
1507 /* When atoms of very different mass are involved in an angle potential
1508 * and both bonds in the angle are constrained, the dynamic modes in such
1509 * angles have very different periods and significant energy exchange
1510 * takes several nanoseconds. Thus even a small amount of error in
1511 * different algorithms can lead to violation of equipartitioning.
1512 * The parameters below are mainly based on an all-atom chloroform model
1513 * with all bonds constrained. Equipartitioning is off by more than 1%
1514 * (need to run 5-10 ns) when the difference in mass between H and Cl
1515 * is more than a factor 13 and the accuracy is less than the thresholds
1516 * given below. This has been verified on some other molecules.
1518 * Note that the buffer and shake parameters have unit length and
1519 * energy/time, respectively, so they will "only" work correctly
1520 * for atomistic force fields using MD units.
1522 const real massFactorThreshold = 13.0;
1523 const real bufferToleranceThreshold = 1e-4;
1524 const int lincsIterationThreshold = 2;
1525 const int lincsOrderThreshold = 4;
1526 const real shakeToleranceThreshold = 0.005 * ir->delta_t;
1528 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold
1529 && ir->nProjOrder >= lincsOrderThreshold);
1530 bool shakeWithSufficientTolerance =
1531 (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1 * shakeToleranceThreshold);
1532 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol <= 1.1 * bufferToleranceThreshold
1533 && (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1535 return;
1538 bool haveDecoupledMode = false;
1539 for (const gmx_moltype_t& molt : mtop->moltype)
1541 if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams, massFactorThreshold))
1543 haveDecoupledMode = true;
1547 if (haveDecoupledMode)
1549 std::string message = gmx::formatString(
1550 "There are atoms at both ends of an angle, connected by constraints "
1551 "and with masses that differ by more than a factor of %g. This means "
1552 "that there are likely dynamic modes that are only very weakly coupled.",
1553 massFactorThreshold);
1554 if (ir->cutoff_scheme == ecutsVERLET)
1556 message += gmx::formatString(
1557 " To ensure good equipartitioning, you need to either not use "
1558 "constraints on all bonds (but, if possible, only on bonds involving "
1559 "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1560 "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1561 ">= %d or SHAKE tolerance <= %g",
1562 ei_names[eiSD1], bufferToleranceThreshold, lincsIterationThreshold,
1563 lincsOrderThreshold, shakeToleranceThreshold);
1565 else
1567 message += gmx::formatString(
1568 " To ensure good equipartitioning, we suggest to switch to the %s "
1569 "cutoff-scheme, since that allows for better control over the Verlet "
1570 "buffer size and thus over the energy drift.",
1571 ecutscheme_names[ecutsVERLET]);
1573 warning(wi, message);
1577 static void set_verlet_buffer(const gmx_mtop_t* mtop, t_inputrec* ir, real buffer_temp, matrix box, warninp* wi)
1579 char warn_buf[STRLEN];
1581 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol,
1582 buffer_temp);
1584 /* Calculate the buffer size for simple atom vs atoms list */
1585 VerletbufListSetup listSetup1x1;
1586 listSetup1x1.cluster_size_i = 1;
1587 listSetup1x1.cluster_size_j = 1;
1588 const real rlist_1x1 = calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1589 buffer_temp, listSetup1x1);
1591 /* Set the pair-list buffer size in ir */
1592 VerletbufListSetup listSetup4x4 = verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1593 ir->rlist = calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1594 buffer_temp, listSetup4x4);
1596 const int n_nonlin_vsite = countNonlinearVsites(*mtop);
1597 if (n_nonlin_vsite > 0)
1599 sprintf(warn_buf,
1600 "There are %d non-linear virtual site constructions. Their contribution to the "
1601 "energy error is approximated. In most cases this does not affect the error "
1602 "significantly.",
1603 n_nonlin_vsite);
1604 warning_note(wi, warn_buf);
1607 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n", 1, 1,
1608 rlist_1x1, rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
1610 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1611 listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j, ir->rlist,
1612 ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
1614 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1616 if (gmx::square(ir->rlist) >= max_cutoff2(ir->pbcType, box))
1618 gmx_fatal(FARGS,
1619 "The pair-list cut-off (%g nm) is longer than half the shortest box vector or "
1620 "longer than the smallest box diagonal element (%g nm). Increase the box size or "
1621 "decrease nstlist or increase verlet-buffer-tolerance.",
1622 ir->rlist, std::sqrt(max_cutoff2(ir->pbcType, box)));
1626 int gmx_grompp(int argc, char* argv[])
1628 const char* desc[] = {
1629 "[THISMODULE] (the gromacs preprocessor)",
1630 "reads a molecular topology file, checks the validity of the",
1631 "file, expands the topology from a molecular description to an atomic",
1632 "description. The topology file contains information about",
1633 "molecule types and the number of molecules, the preprocessor",
1634 "copies each molecule as needed. ",
1635 "There is no limitation on the number of molecule types. ",
1636 "Bonds and bond-angles can be converted into constraints, separately",
1637 "for hydrogens and heavy atoms.",
1638 "Then a coordinate file is read and velocities can be generated",
1639 "from a Maxwellian distribution if requested.",
1640 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1641 "(eg. number of MD steps, time step, cut-off), and others such as",
1642 "NEMD parameters, which are corrected so that the net acceleration",
1643 "is zero.",
1644 "Eventually a binary file is produced that can serve as the sole input",
1645 "file for the MD program.[PAR]",
1647 "[THISMODULE] uses the atom names from the topology file. The atom names",
1648 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1649 "warnings when they do not match the atom names in the topology.",
1650 "Note that the atom names are irrelevant for the simulation as",
1651 "only the atom types are used for generating interaction parameters.[PAR]",
1653 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1654 "etc. The preprocessor supports the following keywords::",
1656 " #ifdef VARIABLE",
1657 " #ifndef VARIABLE",
1658 " #else",
1659 " #endif",
1660 " #define VARIABLE",
1661 " #undef VARIABLE",
1662 " #include \"filename\"",
1663 " #include <filename>",
1665 "The functioning of these statements in your topology may be modulated by",
1666 "using the following two flags in your [REF].mdp[ref] file::",
1668 " define = -DVARIABLE1 -DVARIABLE2",
1669 " include = -I/home/john/doe",
1671 "For further information a C-programming textbook may help you out.",
1672 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1673 "topology file written out so that you can verify its contents.[PAR]",
1675 "When using position restraints, a file with restraint coordinates",
1676 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1677 "for [TT]-c[tt]). For free energy calculations, separate reference",
1678 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1679 "otherwise they will be equal to those of the A topology.[PAR]",
1681 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1682 "The last frame with coordinates and velocities will be read,",
1683 "unless the [TT]-time[tt] option is used. Only if this information",
1684 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1685 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1686 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1687 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1688 "variables.[PAR]",
1690 "[THISMODULE] can be used to restart simulations (preserving",
1691 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1692 "However, for simply changing the number of run steps to extend",
1693 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1694 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1695 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1696 "like output frequency, then supplying the checkpoint file to",
1697 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1698 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1699 "the ensemble (if possible) still requires passing the checkpoint",
1700 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1702 "By default, all bonded interactions which have constant energy due to",
1703 "virtual site constructions will be removed. If this constant energy is",
1704 "not zero, this will result in a shift in the total energy. All bonded",
1705 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1706 "all constraints for distances which will be constant anyway because",
1707 "of virtual site constructions will be removed. If any constraints remain",
1708 "which involve virtual sites, a fatal error will result.[PAR]",
1710 "To verify your run input file, please take note of all warnings",
1711 "on the screen, and correct where necessary. Do also look at the contents",
1712 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1713 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1714 "with the [TT]-debug[tt] option which will give you more information",
1715 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1716 "can see the contents of the run input file with the [gmx-dump]",
1717 "program. [gmx-check] can be used to compare the contents of two",
1718 "run input files.[PAR]",
1720 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1721 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1722 "harmless, but usually they are not. The user is advised to carefully",
1723 "interpret the output messages before attempting to bypass them with",
1724 "this option."
1726 t_gromppopts* opts;
1727 std::vector<MoleculeInformation> mi;
1728 std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1729 int nvsite, comb;
1730 real fudgeQQ;
1731 double reppow;
1732 const char* mdparin;
1733 bool bNeedVel, bGenVel;
1734 gmx_bool have_atomnumber;
1735 gmx_output_env_t* oenv;
1736 gmx_bool bVerbose = FALSE;
1737 warninp* wi;
1738 char warn_buf[STRLEN];
1740 t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
1741 { efMDP, "-po", "mdout", ffWRITE },
1742 { efSTX, "-c", nullptr, ffREAD },
1743 { efSTX, "-r", "restraint", ffOPTRD },
1744 { efSTX, "-rb", "restraint", ffOPTRD },
1745 { efNDX, nullptr, nullptr, ffOPTRD },
1746 { efTOP, nullptr, nullptr, ffREAD },
1747 { efTOP, "-pp", "processed", ffOPTWR },
1748 { efTPR, "-o", nullptr, ffWRITE },
1749 { efTRN, "-t", nullptr, ffOPTRD },
1750 { efEDR, "-e", nullptr, ffOPTRD },
1751 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1752 { efGRO, "-imd", "imdgroup", ffOPTWR },
1753 { efTRN, "-ref", "rotref", ffOPTRW } };
1754 #define NFILE asize(fnm)
1756 /* Command line options */
1757 gmx_bool bRenum = TRUE;
1758 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1759 int i, maxwarn = 0;
1760 real fr_time = -1;
1761 t_pargs pa[] = {
1762 { "-v", FALSE, etBOOL, { &bVerbose }, "Be loud and noisy" },
1763 { "-time", FALSE, etREAL, { &fr_time }, "Take frame at or first after this time." },
1764 { "-rmvsbds",
1765 FALSE,
1766 etBOOL,
1767 { &bRmVSBds },
1768 "Remove constant bonded interactions with virtual sites" },
1769 { "-maxwarn",
1770 FALSE,
1771 etINT,
1772 { &maxwarn },
1773 "Number of allowed warnings during input processing. Not for normal use and may "
1774 "generate unstable systems" },
1775 { "-zero",
1776 FALSE,
1777 etBOOL,
1778 { &bZero },
1779 "Set parameters for bonded interactions without defaults to zero instead of "
1780 "generating an error" },
1781 { "-renum",
1782 FALSE,
1783 etBOOL,
1784 { &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, asize(desc), desc, 0, nullptr, &oenv))
1791 return 0;
1794 /* Initiate some variables */
1795 gmx::MDModules mdModules;
1796 t_inputrec irInstance;
1797 t_inputrec* ir = &irInstance;
1798 snew(opts, 1);
1799 snew(opts->include, STRLEN);
1800 snew(opts->define, STRLEN);
1802 wi = init_warning(TRUE, maxwarn);
1804 /* PARAMETER file processing */
1805 mdparin = opt2fn("-f", NFILE, fnm);
1806 set_warning_line(wi, mdparin, -1);
1809 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1811 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1813 if (bVerbose)
1815 fprintf(stderr, "checking input for internal consistency...\n");
1817 check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
1819 if (ir->ld_seed == -1)
1821 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1822 fprintf(stderr, "Setting the LD random seed to %" PRId64 "\n", ir->ld_seed);
1825 if (ir->expandedvals->lmc_seed == -1)
1827 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1828 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1831 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1832 bGenVel = (bNeedVel && opts->bGenVel);
1833 if (bGenVel && ir->bContinuation)
1835 sprintf(warn_buf,
1836 "Generating velocities is inconsistent with attempting "
1837 "to continue a previous run. Choose only one of "
1838 "gen-vel = yes and continuation = yes.");
1839 warning_error(wi, warn_buf);
1842 std::array<InteractionsOfType, F_NRE> interactions;
1843 gmx_mtop_t sys;
1844 PreprocessingAtomTypes atypes;
1845 if (debug)
1847 pr_symtab(debug, 0, "Just opened", &sys.symtab);
1850 const char* fn = ftp2fn(efTOP, NFILE, fnm);
1851 if (!gmx_fexist(fn))
1853 gmx_fatal(FARGS, "%s does not exist", fn);
1856 t_state state;
1857 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm), opts, ir, bZero,
1858 bGenVel, bVerbose, &state, &atypes, &sys, &mi, &intermolecular_interactions,
1859 interactions, &comb, &reppow, &fudgeQQ, opts->bMorse, wi);
1861 if (debug)
1863 pr_symtab(debug, 0, "After new_status", &sys.symtab);
1866 nvsite = 0;
1867 /* set parameters for virtual site construction (not for vsiten) */
1868 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1870 nvsite += set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions);
1872 /* now throw away all obsolete bonds, angles and dihedrals: */
1873 /* note: constraints are ALWAYS removed */
1874 if (nvsite)
1876 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1878 clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds);
1882 if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == econtSHAKE))
1884 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1886 sprintf(warn_buf, "Can not do %s with %s, use %s", EI(ir->eI),
1887 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1888 warning_error(wi, warn_buf);
1890 if (ir->bPeriodicMols)
1892 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1893 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1894 warning_error(wi, warn_buf);
1898 if (EI_SD(ir->eI) && ir->etc != etcNO)
1900 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1903 /* If we are doing QM/MM, check that we got the atom numbers */
1904 have_atomnumber = TRUE;
1905 for (gmx::index i = 0; i < gmx::ssize(atypes); i++)
1907 have_atomnumber = have_atomnumber && (atypes.atomNumberFromAtomType(i) >= 0);
1909 if (!have_atomnumber && ir->bQMMM)
1911 warning_error(
1913 "\n"
1914 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1915 "field you are using does not contain atom numbers fields. This is an\n"
1916 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1917 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1918 "an integer just before the mass column in ffXXXnb.itp.\n"
1919 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1922 /* Check for errors in the input now, since they might cause problems
1923 * during processing further down.
1925 check_warning_error(wi, FARGS);
1927 if (nint_ftype(&sys, mi, F_POSRES) > 0 || nint_ftype(&sys, mi, F_FBPOSRES) > 0)
1929 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1931 sprintf(warn_buf,
1932 "You are combining position restraints with %s pressure coupling, which can "
1933 "lead to instabilities. If you really want to combine position restraints with "
1934 "pressure coupling, we suggest to use %s pressure coupling instead.",
1935 EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
1936 warning_note(wi, warn_buf);
1939 const char* fn = opt2fn("-r", NFILE, fnm);
1940 const char* fnB;
1942 if (!gmx_fexist(fn))
1944 gmx_fatal(FARGS,
1945 "Cannot find position restraint file %s (option -r).\n"
1946 "From GROMACS-2018, you need to specify the position restraint "
1947 "coordinate files explicitly to avoid mistakes, although you can "
1948 "still use the same file as you specify for the -c option.",
1949 fn);
1952 if (opt2bSet("-rb", NFILE, fnm))
1954 fnB = opt2fn("-rb", NFILE, fnm);
1955 if (!gmx_fexist(fnB))
1957 gmx_fatal(FARGS,
1958 "Cannot find B-state position restraint file %s (option -rb).\n"
1959 "From GROMACS-2018, you need to specify the position restraint "
1960 "coordinate files explicitly to avoid mistakes, although you can "
1961 "still use the same file as you specify for the -c option.",
1962 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, ir->refcoord_scaling, ir->pbcType, ir->posres_com,
1983 ir->posres_comB, wi);
1986 /* If we are using CMAP, setup the pre-interpolation grid */
1987 if (interactions[F_CMAP].ncmap() > 0)
1989 init_cmap_grid(&sys.ffparams.cmap_grid, interactions[F_CMAP].cmapAngles,
1990 interactions[F_CMAP].cmakeGridSpacing);
1991 setup_cmap(interactions[F_CMAP].cmakeGridSpacing, interactions[F_CMAP].cmapAngles,
1992 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(), comb,
2021 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(
2062 "Zero-step energy minimization will alter the coordinates before calculating the "
2063 "energy. If you just want the energy of a single point, try zero-step MD (with "
2064 "unconstrained_start = yes). To do multiple single-point energy evaluations of "
2065 "different configurations of the same topology, use mdrun -rerun.");
2068 check_warning_error(wi, FARGS);
2070 if (bVerbose)
2072 fprintf(stderr, "initialising group options...\n");
2074 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifier(), ir, wi);
2076 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2078 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2080 real buffer_temp;
2082 if (EI_MD(ir->eI) && ir->etc == etcNO)
2084 if (bGenVel)
2086 buffer_temp = opts->tempi;
2088 else
2090 buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2092 if (buffer_temp > 0)
2094 sprintf(warn_buf,
2095 "NVE simulation: will use the initial temperature of %.3f K for "
2096 "determining the Verlet buffer size",
2097 buffer_temp);
2098 warning_note(wi, warn_buf);
2100 else
2102 sprintf(warn_buf,
2103 "NVE simulation with an initial temperature of zero: will use a Verlet "
2104 "buffer of %d%%. Check your energy drift!",
2105 gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
2106 warning_note(wi, warn_buf);
2109 else
2111 buffer_temp = get_max_reference_temp(ir, wi);
2114 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2116 /* NVE with initial T=0: we add a fixed ratio to rlist.
2117 * Since we don't actually use verletbuf_tol, we set it to -1
2118 * so it can't be misused later.
2120 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2121 ir->verletbuf_tol = -1;
2123 else
2125 /* We warn for NVE simulations with a drift tolerance that
2126 * might result in a 1(.1)% drift over the total run-time.
2127 * Note that we can't warn when nsteps=0, since we don't
2128 * know how many steps the user intends to run.
2130 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 && ir->nsteps > 0)
2132 const real driftTolerance = 0.01;
2133 /* We use 2 DOF per atom = 2kT pot+kin energy,
2134 * to be on the safe side with constraints.
2136 const real totalEnergyDriftPerAtomPerPicosecond =
2137 2 * BOLTZ * buffer_temp / (ir->nsteps * ir->delta_t);
2139 if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
2141 sprintf(warn_buf,
2142 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
2143 "NVE simulation of length %g ps, which can give a final drift of "
2144 "%d%%. For conserving energy to %d%% when using constraints, you "
2145 "might need to set verlet-buffer-tolerance to %.1e.",
2146 ir->verletbuf_tol, ir->nsteps * ir->delta_t,
2147 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
2148 gmx::roundToInt(100 * driftTolerance),
2149 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
2150 warning_note(wi, warn_buf);
2154 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi);
2159 /* Init the temperature coupling state */
2160 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2162 if (debug)
2164 pr_symtab(debug, 0, "After index", &sys.symtab);
2167 triple_check(mdparin, ir, &sys, wi);
2168 close_symtab(&sys.symtab);
2169 if (debug)
2171 pr_symtab(debug, 0, "After close", &sys.symtab);
2174 /* make exclusions between QM atoms and remove charges if needed */
2175 if (ir->bQMMM)
2177 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_ORIGINAL);
2178 if (ir->QMMMscheme != eQMMMschemeoniom)
2180 std::vector<int> qmmmAtoms = qmmmAtomIndices(*ir, sys);
2181 removeQmmmAtomCharges(&sys, qmmmAtoms);
2185 if (ir->eI == eiMimic)
2187 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_MIMIC);
2190 if (ftp2bSet(efTRN, NFILE, fnm))
2192 if (bVerbose)
2194 fprintf(stderr, "getting data from old trajectory ...\n");
2196 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm), bNeedVel, bGenVel,
2197 fr_time, ir, &state, &sys, oenv);
2200 if (ir->pbcType == PbcType::XY && ir->nwall != 2)
2202 clear_rvec(state.box[ZZ]);
2205 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2207 /* Calculate the optimal grid dimensions */
2208 matrix scaledBox;
2209 EwaldBoxZScaler boxScaler(*ir);
2210 boxScaler.scaleBox(state.box, scaledBox);
2212 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2214 /* Mark fourier_spacing as not used */
2215 ir->fourier_spacing = 0;
2217 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2219 set_warning_line(wi, mdparin, -1);
2220 warning_error(
2221 wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2223 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2224 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize, &(ir->nkx), &(ir->nky),
2225 &(ir->nkz));
2226 if (ir->nkx < minGridSize || ir->nky < minGridSize || ir->nkz < minGridSize)
2228 warning_error(wi,
2229 "The PME grid size should be >= 2*(pme-order - 1); either manually "
2230 "increase the grid size or decrease pme-order");
2234 /* MRS: eventually figure out better logic for initializing the fep
2235 values that makes declaring the lambda and declaring the state not
2236 potentially conflict if not handled correctly. */
2237 if (ir->efep != efepNO)
2239 state.fep_state = ir->fepvals->init_fep_state;
2240 for (i = 0; i < efptNR; i++)
2242 /* init_lambda trumps state definitions*/
2243 if (ir->fepvals->init_lambda >= 0)
2245 state.lambda[i] = ir->fepvals->init_lambda;
2247 else
2249 if (ir->fepvals->all_lambda[i] == nullptr)
2251 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2253 else
2255 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2261 pull_t* pull = nullptr;
2263 if (ir->bPull)
2265 pull = set_pull_init(ir, &sys, state.x.rvec_array(), state.box, state.lambda[efptMASS], wi);
2268 /* Modules that supply external potential for pull coordinates
2269 * should register those potentials here. finish_pull() will check
2270 * that providers have been registerd for all external potentials.
2273 if (ir->bDoAwh)
2275 tensor compressibility = { { 0 } };
2276 if (ir->epc != epcNO)
2278 copy_mat(ir->compress, compressibility);
2280 setStateDependentAwhParams(ir->awhParams, ir->pull, pull, state.box, ir->pbcType,
2281 compressibility, &ir->opts, wi);
2284 if (ir->bPull)
2286 finish_pull(pull);
2289 if (ir->bRot)
2291 set_reference_positions(ir->rot, state.x.rvec_array(), state.box,
2292 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm), wi);
2295 /* reset_multinr(sys); */
2297 if (EEL_PME(ir->coulombtype))
2299 float ratio = pme_load_estimate(sys, *ir, state.box);
2300 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2301 /* With free energy we might need to do PME both for the A and B state
2302 * charges. This will double the cost, but the optimal performance will
2303 * then probably be at a slightly larger cut-off and grid spacing.
2305 if ((ir->efep == efepNO && ratio > 1.0 / 2.0) || (ir->efep != efepNO && ratio > 2.0 / 3.0))
2307 warning_note(
2309 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2310 "and for highly parallel simulations between 0.25 and 0.33,\n"
2311 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2312 if (ir->efep != efepNO)
2314 warning_note(wi,
2315 "For free energy simulations, the optimal load limit increases from "
2316 "0.5 to 0.667\n");
2322 char warn_buf[STRLEN];
2323 double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2324 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2325 if (cio > 2000)
2327 set_warning_line(wi, mdparin, -1);
2328 warning_note(wi, warn_buf);
2330 else
2332 printf("%s\n", warn_buf);
2336 // Add the md modules internal parameters that are not mdp options
2337 // e.g., atom indices
2340 gmx::KeyValueTreeBuilder internalParameterBuilder;
2341 mdModules.notifier().notifier_.notify(internalParameterBuilder.rootObject());
2342 ir->internalParameters =
2343 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2346 if (bVerbose)
2348 fprintf(stderr, "writing run input file...\n");
2351 done_warning(wi, FARGS);
2352 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, &sys);
2354 /* Output IMD group, if bIMD is TRUE */
2355 gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, &sys, NFILE, fnm);
2357 sfree(opts->define);
2358 sfree(opts->include);
2359 sfree(opts);
2360 for (auto& mol : mi)
2362 // Some of the contents of molinfo have been stolen, so
2363 // fullCleanUp can't be called.
2364 mol.partialCleanUp();
2366 done_inputrec_strings();
2367 output_env_done(oenv);
2369 return 0;