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.
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/listoflists.h"
106 #include "gromacs/utility/mdmodulenotification.h"
107 #include "gromacs/utility/smalloc.h"
108 #include "gromacs/utility/snprintf.h"
110 /* TODO The implementation details should move to their own source file. */
111 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int> atoms,
112 gmx::ArrayRef<const real> params,
113 const std::string& name) :
114 atoms_(atoms.begin(), atoms.end()),
115 interactionTypeName_(name)
118 params.size() <= forceParam_.size(),
119 gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM)
121 auto forceParamIt = forceParam_.begin();
122 for (const auto param : params)
124 *forceParamIt++ = param;
126 while (forceParamIt != forceParam_.end())
128 *forceParamIt++ = NOTSET;
132 const int& InteractionOfType::ai() const
134 GMX_RELEASE_ASSERT(!atoms_.empty(), "Need to have at least one atom set");
138 const int& InteractionOfType::aj() const
140 GMX_RELEASE_ASSERT(atoms_.size() > 1, "Need to have at least two atoms set");
144 const int& InteractionOfType::ak() const
146 GMX_RELEASE_ASSERT(atoms_.size() > 2, "Need to have at least three atoms set");
150 const int& InteractionOfType::al() const
152 GMX_RELEASE_ASSERT(atoms_.size() > 3, "Need to have at least four atoms set");
156 const int& InteractionOfType::am() const
158 GMX_RELEASE_ASSERT(atoms_.size() > 4, "Need to have at least five atoms set");
162 const real& InteractionOfType::c0() const
164 return forceParam_[0];
167 const real& InteractionOfType::c1() const
169 return forceParam_[1];
172 const real& InteractionOfType::c2() const
174 return forceParam_[2];
177 const std::string& InteractionOfType::interactionTypeName() const
179 return interactionTypeName_;
182 void InteractionOfType::sortBondAtomIds()
186 std::swap(atoms_[0], atoms_[1]);
190 void InteractionOfType::sortAngleAtomIds()
194 std::swap(atoms_[0], atoms_[2]);
198 void InteractionOfType::sortDihedralAtomIds()
202 std::swap(atoms_[0], atoms_[3]);
203 std::swap(atoms_[1], atoms_[2]);
207 void InteractionOfType::sortAtomIds()
217 else if (isDihedral())
219 sortDihedralAtomIds();
223 GMX_THROW(gmx::InternalError(
224 "Can not sort parameters other than bonds, angles or dihedrals"));
228 void InteractionOfType::setForceParameter(int pos, real value)
230 GMX_RELEASE_ASSERT(pos < MAXFORCEPARAM,
231 "Can't set parameter beyond the maximum number of parameters");
232 forceParam_[pos] = value;
235 void MoleculeInformation::initMolInfo()
239 init_t_atoms(&atoms, 0, FALSE);
242 void MoleculeInformation::partialCleanUp()
247 void MoleculeInformation::fullCleanUp()
253 static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
256 /* For all the molecule types */
257 for (auto& mol : mols)
259 n += mol.interactions[ifunc].size();
260 mol.interactions[ifunc].interactionTypes.clear();
265 static int check_atom_names(const char* fn1, const char* fn2, gmx_mtop_t* mtop, const t_atoms* at)
267 int m, i, j, nmismatch;
269 #define MAXMISMATCH 20
271 if (mtop->natoms != at->nr)
273 gmx_incons("comparing atom names");
278 for (const gmx_molblock_t& molb : mtop->molblock)
280 tat = &mtop->moltype[molb.type].atoms;
281 for (m = 0; m < molb.nmol; m++)
283 for (j = 0; j < tat->nr; j++)
285 if (strcmp(*(tat->atomname[j]), *(at->atomname[i])) != 0)
287 if (nmismatch < MAXMISMATCH)
290 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
291 i + 1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
293 else if (nmismatch == MAXMISMATCH)
295 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
307 static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
309 /* This check is not intended to ensure accurate integration,
310 * rather it is to signal mistakes in the mdp settings.
311 * A common mistake is to forget to turn on constraints
312 * for MD after energy minimization with flexible bonds.
313 * This check can also detect too large time steps for flexible water
314 * models, but such errors will often be masked by the constraints
315 * mdp options, which turns flexible water into water with bond constraints,
316 * but without an angle constraint. Unfortunately such incorrect use
317 * of water models can not easily be detected without checking
318 * for specific model names.
320 * The stability limit of leap-frog or velocity verlet is 4.44 steps
321 * per oscillational period.
322 * But accurate bonds distributions are lost far before that limit.
323 * To allow relatively common schemes (although not common with Gromacs)
324 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
325 * we set the note limit to 10.
327 int min_steps_warn = 5;
328 int min_steps_note = 10;
330 int i, a1, a2, w_a1, w_a2, j;
331 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
332 bool bFound, bWater, bWarn;
333 char warn_buf[STRLEN];
335 /* Get the interaction parameters */
336 gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
338 twopi2 = gmx::square(2 * M_PI);
340 limit2 = gmx::square(min_steps_note * dt);
345 const gmx_moltype_t* w_moltype = nullptr;
346 for (const gmx_moltype_t& moltype : mtop->moltype)
348 const t_atom* atom = moltype.atoms.atom;
349 const InteractionLists& ilist = moltype.ilist;
350 const InteractionList& ilc = ilist[F_CONSTR];
351 const InteractionList& ils = ilist[F_SETTLE];
352 for (ftype = 0; ftype < F_NRE; ftype++)
354 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
359 const InteractionList& ilb = ilist[ftype];
360 for (i = 0; i < ilb.size(); i += 3)
362 fc = ip[ilb.iatoms[i]].harmonic.krA;
363 re = ip[ilb.iatoms[i]].harmonic.rA;
364 if (ftype == F_G96BONDS)
366 /* Convert squared sqaure fc to harmonic fc */
369 a1 = ilb.iatoms[i + 1];
370 a2 = ilb.iatoms[i + 2];
373 if (fc > 0 && m1 > 0 && m2 > 0)
375 period2 = twopi2 * m1 * m2 / ((m1 + m2) * fc);
379 period2 = GMX_FLOAT_MAX;
383 fprintf(debug, "fc %g m1 %g m2 %g period %g\n", fc, m1, m2, std::sqrt(period2));
385 if (period2 < limit2)
388 for (j = 0; j < ilc.size(); j += 3)
390 if ((ilc.iatoms[j + 1] == a1 && ilc.iatoms[j + 2] == a2)
391 || (ilc.iatoms[j + 1] == a2 && ilc.iatoms[j + 2] == a1))
396 for (j = 0; j < ils.size(); j += 4)
398 if ((a1 == ils.iatoms[j + 1] || a1 == ils.iatoms[j + 2] || a1 == ils.iatoms[j + 3])
399 && (a2 == ils.iatoms[j + 1] || a2 == ils.iatoms[j + 2]
400 || a2 == ils.iatoms[j + 3]))
405 if (!bFound && (w_moltype == nullptr || period2 < w_period2))
407 w_moltype = &moltype;
417 if (w_moltype != nullptr)
419 bWarn = (w_period2 < gmx::square(min_steps_warn * dt));
420 /* A check that would recognize most water models */
421 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' && w_moltype->atoms.nr <= 5);
423 "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated "
424 "oscillational period of %.1e ps, which is less than %d times the time step of "
427 *w_moltype->name, w_a1 + 1, *w_moltype->atoms.atomname[w_a1], w_a2 + 1,
428 *w_moltype->atoms.atomname[w_a2], std::sqrt(w_period2),
429 bWarn ? min_steps_warn : min_steps_note, dt,
430 bWater ? "Maybe you asked for fexible water."
431 : "Maybe you forgot to change the constraints mdp option.");
434 warning(wi, warn_buf);
438 warning_note(wi, warn_buf);
443 static void check_vel(gmx_mtop_t* mtop, rvec v[])
445 for (const AtomProxy atomP : AtomRange(*mtop))
447 const t_atom& local = atomP.atom();
448 int i = atomP.globalAtomNumber();
449 if (local.ptype == eptShell || local.ptype == eptBond || local.ptype == eptVSite)
456 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
459 char warn_buf[STRLEN];
461 for (const AtomProxy atomP : AtomRange(*mtop))
463 const t_atom& local = atomP.atom();
464 if (local.ptype == eptShell || local.ptype == eptBond)
469 if ((nshells > 0) && (ir->nstcalcenergy != 1))
471 set_warning_line(wi, "unknown", -1);
472 snprintf(warn_buf, STRLEN, "There are %d shells, changing nstcalcenergy from %d to 1",
473 nshells, ir->nstcalcenergy);
474 ir->nstcalcenergy = 1;
475 warning(wi, warn_buf);
479 /* TODO Decide whether this function can be consolidated with
480 * gmx_mtop_ftype_count */
481 static int nint_ftype(gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
484 for (const gmx_molblock_t& molb : mtop->molblock)
486 nint += molb.nmol * mi[molb.type].interactions[ftype].size();
492 /* This routine reorders the molecule type array
493 * in the order of use in the molblocks,
494 * unused molecule types are deleted.
496 static void renumber_moltypes(gmx_mtop_t* sys, std::vector<MoleculeInformation>* molinfo)
499 std::vector<int> order;
500 for (gmx_molblock_t& molblock : sys->molblock)
502 const auto found = std::find_if(order.begin(), order.end(), [&molblock](const auto& entry) {
503 return molblock.type == entry;
505 if (found == order.end())
507 /* This type did not occur yet, add it */
508 order.push_back(molblock.type);
509 molblock.type = order.size() - 1;
513 molblock.type = std::distance(order.begin(), found);
517 /* We still need to reorder the molinfo structs */
518 std::vector<MoleculeInformation> minew(order.size());
520 for (auto& mi : *molinfo)
522 const auto found = std::find(order.begin(), order.end(), index);
523 if (found != order.end())
525 int pos = std::distance(order.begin(), found);
530 // Need to manually clean up memory ....
539 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t* mtop)
541 mtop->moltype.resize(mi.size());
543 for (const auto& mol : mi)
545 gmx_moltype_t& molt = mtop->moltype[pos];
546 molt.name = mol.name;
547 molt.atoms = mol.atoms;
548 /* ilists are copied later */
549 molt.excls = mol.excls;
554 static void new_status(const char* topfile,
555 const char* topppfile,
563 PreprocessingAtomTypes* atypes,
565 std::vector<MoleculeInformation>* mi,
566 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
567 gmx::ArrayRef<InteractionsOfType> interactions,
574 std::vector<gmx_molblock_t> molblock;
576 bool ffParametrizedWithHBondConstraints;
578 char warn_buf[STRLEN];
580 /* TOPOLOGY processing */
581 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab), interactions,
582 comb, reppow, fudgeQQ, atypes, mi, intermolecular_interactions, ir,
583 &molblock, &ffParametrizedWithHBondConstraints, wi);
585 sys->molblock.clear();
588 for (const gmx_molblock_t& molb : molblock)
590 if (!sys->molblock.empty() && molb.type == sys->molblock.back().type)
592 /* Merge consecutive blocks with the same molecule type */
593 sys->molblock.back().nmol += molb.nmol;
595 else if (molb.nmol > 0)
597 /* Add a new molblock to the topology */
598 sys->molblock.push_back(molb);
600 sys->natoms += molb.nmol * (*mi)[sys->molblock.back().type].atoms.nr;
602 if (sys->molblock.empty())
604 gmx_fatal(FARGS, "No molecules were defined in the system");
607 renumber_moltypes(sys, mi);
611 convert_harmonics(*mi, atypes);
614 if (ir->eDisre == edrNone)
616 i = rm_interactions(F_DISRES, *mi);
619 set_warning_line(wi, "unknown", -1);
620 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
621 warning_note(wi, warn_buf);
626 i = rm_interactions(F_ORIRES, *mi);
629 set_warning_line(wi, "unknown", -1);
630 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
631 warning_note(wi, warn_buf);
635 /* Copy structures from msys to sys */
636 molinfo2mtop(*mi, sys);
638 gmx_mtop_finalize(sys);
640 /* Discourage using the, previously common, setup of applying constraints
641 * to all bonds with force fields that have been parametrized with H-bond
642 * constraints only. Do not print note with large timesteps or vsites.
644 if (opts->nshake == eshALLBONDS && ffParametrizedWithHBondConstraints && ir->delta_t < 0.0026
645 && gmx_mtop_ftype_count(sys, F_VSITE3FD) == 0)
647 set_warning_line(wi, "unknown", -1);
649 "You are using constraints on all bonds, whereas the forcefield "
650 "has been parametrized only with constraints involving hydrogen atoms. "
651 "We suggest using constraints = h-bonds instead, this will also improve "
655 /* COORDINATE file processing */
658 fprintf(stderr, "processing coordinates...\n");
665 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
666 state->natoms = conftop->atoms.nr;
667 if (state->natoms != sys->natoms)
670 "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);
680 state->flags |= (1 << estV);
682 state_change_natoms(state, state->natoms);
683 std::copy(x, x + state->natoms, state->x.data());
687 std::copy(v, v + state->natoms, state->v.data());
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);
700 "%d non-matching atom name%s\n"
701 "atom names from %s will be used\n"
702 "atom names from %s will be ignored\n",
703 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
707 /* Do more checks, mostly related to constraints */
710 fprintf(stderr, "double-checking input for internal consistency...\n");
713 bool bHasNormalConstraints =
714 0 < (nint_ftype(sys, *mi, F_CONSTR) + nint_ftype(sys, *mi, F_CONSTRNC));
715 bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
716 double_check(ir, state->box, bHasNormalConstraints, bHasAnyConstraints, wi);
723 snew(mass, state->natoms);
724 for (const AtomProxy atomP : AtomRange(*sys))
726 const t_atom& local = atomP.atom();
727 int i = atomP.globalAtomNumber();
731 if (opts->seed == -1)
733 opts->seed = static_cast<int>(gmx::makeRandomSeed());
734 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
736 state->flags |= (1 << estV);
737 maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array());
739 stop_cm(stdout, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
744 static void copy_state(const char* slog, t_trxframe* fr, bool bReadVel, t_state* state, double* use_time)
746 if (fr->not_ok & FRAME_NOT_OK)
748 gmx_fatal(FARGS, "Can not start from an incomplete frame");
752 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s", slog);
755 std::copy(fr->x, fr->x + state->natoms, state->x.data());
760 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
762 std::copy(fr->v, fr->v + state->natoms, state->v.data());
766 copy_mat(fr->box, state->box);
769 *use_time = fr->time;
772 static void cont_status(const char* slog,
780 const gmx_output_env_t* oenv)
781 /* If fr_time == -1 read the last frame available which is complete */
788 bReadVel = (bNeedVel && !bGenVel);
790 fprintf(stderr, "Reading Coordinates%s and Box size from old trajectory\n",
791 bReadVel ? ", Velocities" : "");
794 fprintf(stderr, "Will read whole trajectory\n");
798 fprintf(stderr, "Will read till time %g\n", fr_time);
805 "Velocities generated: "
806 "ignoring velocities in input trajectory\n");
808 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
812 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
818 "WARNING: Did not find a frame with velocities in file %s,\n"
819 " all velocities will be set to zero!\n\n",
821 for (auto& vi : makeArrayRef(state->v))
826 /* Search for a frame without velocities */
828 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
832 state->natoms = fr.natoms;
834 if (sys->natoms != state->natoms)
837 "Number of atoms in Topology "
838 "is not the same as in Trajectory");
840 copy_state(slog, &fr, bReadVel, state, &use_time);
842 /* Find the appropriate frame */
843 while ((fr_time == -1 || fr.time < fr_time) && read_next_frame(oenv, fp, &fr))
845 copy_state(slog, &fr, bReadVel, state, &use_time);
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,
881 int natoms, npbcdim = 0;
882 char warn_buf[STRLEN];
887 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
888 natoms = top->atoms.nr;
891 if (natoms != mtop->natoms)
894 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
895 "(%d). Will assume that the first %d atoms in the topology and %s match.",
896 fn, natoms, mtop->natoms, std::min(mtop->natoms, natoms), fn);
897 warning(wi, warn_buf);
900 npbcdim = ePBC2npbcdim(ePBC);
901 GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
903 if (rc_scaling != erscNO)
905 copy_mat(box, invbox);
906 for (int j = npbcdim; j < DIM; j++)
908 clear_rvec(invbox[j]);
911 gmx::invertBoxMatrix(invbox, invbox);
914 /* Copy the reference coordinates to mtop */
918 snew(hadAtom, natoms);
919 for (gmx_molblock_t& molb : mtop->molblock)
921 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
922 const InteractionsOfType* pr = &(molinfo[molb.type].interactions[F_POSRES]);
923 const InteractionsOfType* prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
924 if (pr->size() > 0 || prfb->size() > 0)
926 atom = mtop->moltype[molb.type].atoms.atom;
927 for (const auto& restraint : pr->interactionTypes)
929 int ai = restraint.ai();
933 "Position restraint atom index (%d) in moltype '%s' is larger than "
934 "number of atoms in %s (%d).\n",
935 ai + 1, *molinfo[molb.type].name, fn, natoms);
938 if (rc_scaling == erscCOM)
940 /* Determine the center of mass of the posres reference coordinates */
941 for (int j = 0; j < npbcdim; j++)
943 sum[j] += atom[ai].m * x[a + ai][j];
945 totmass += atom[ai].m;
948 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
949 for (const auto& restraint : prfb->interactionTypes)
951 int ai = restraint.ai();
955 "Position restraint atom index (%d) in moltype '%s' is larger than "
956 "number of atoms in %s (%d).\n",
957 ai + 1, *molinfo[molb.type].name, fn, natoms);
959 if (rc_scaling == erscCOM && !hadAtom[ai])
961 /* Determine the center of mass of the posres reference coordinates */
962 for (int j = 0; j < npbcdim; j++)
964 sum[j] += atom[ai].m * x[a + ai][j];
966 totmass += atom[ai].m;
971 molb.posres_xA.resize(nat_molb);
972 for (int i = 0; i < nat_molb; i++)
974 copy_rvec(x[a + i], molb.posres_xA[i]);
979 molb.posres_xB.resize(nat_molb);
980 for (int i = 0; i < nat_molb; i++)
982 copy_rvec(x[a + i], molb.posres_xB[i]);
988 if (rc_scaling == erscCOM)
992 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
994 for (int j = 0; j < npbcdim; j++)
996 com[j] = sum[j] / totmass;
999 "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n",
1000 com[XX], com[YY], com[ZZ]);
1003 if (rc_scaling != erscNO)
1005 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
1007 for (gmx_molblock_t& molb : mtop->molblock)
1009 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
1010 if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1012 std::vector<gmx::RVec>& xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1013 for (int i = 0; i < nat_molb; i++)
1015 for (int j = 0; j < npbcdim; j++)
1017 if (rc_scaling == erscALL)
1019 /* Convert from Cartesian to crystal coordinates */
1020 xp[i][j] *= invbox[j][j];
1021 for (int k = j + 1; k < npbcdim; k++)
1023 xp[i][j] += invbox[k][j] * xp[i][k];
1026 else if (rc_scaling == erscCOM)
1028 /* Subtract the center of mass */
1036 if (rc_scaling == erscCOM)
1038 /* Convert the COM from Cartesian to crystal coordinates */
1039 for (int j = 0; j < npbcdim; j++)
1041 com[j] *= invbox[j][j];
1042 for (int k = j + 1; k < npbcdim; k++)
1044 com[j] += invbox[k][j] * com[k];
1055 static void gen_posres(gmx_mtop_t* mtop,
1056 gmx::ArrayRef<const MoleculeInformation> mi,
1065 read_posres(mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1066 /* It is safer to simply read the b-state posres rather than trying
1067 * to be smart and copy the positions.
1069 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1072 static void set_wall_atomtype(PreprocessingAtomTypes* at, t_gromppopts* opts, t_inputrec* ir, warninp* wi)
1075 char warn_buf[STRLEN];
1079 fprintf(stderr, "Searching the wall atom type(s)\n");
1081 for (i = 0; i < ir->nwall; i++)
1083 ir->wall_atomtype[i] = at->atomTypeFromName(opts->wall_atomtype[i]);
1084 if (ir->wall_atomtype[i] == NOTSET)
1086 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1087 warning_error(wi, warn_buf);
1092 static int nrdf_internal(const t_atoms* atoms)
1097 for (i = 0; i < atoms->nr; i++)
1099 /* Vsite ptype might not be set here yet, so also check the mass */
1100 if ((atoms->atom[i].ptype == eptAtom || atoms->atom[i].ptype == eptNucleus) && atoms->atom[i].m > 0)
1107 case 0: nrdf = 0; break;
1108 case 1: nrdf = 0; break;
1109 case 2: nrdf = 1; break;
1110 default: nrdf = nmass * 3 - 6; break;
1116 static void spline1d(double dx, const double* y, int n, double* u, double* y2)
1124 for (i = 1; i < n - 1; i++)
1126 p = 0.5 * y2[i - 1] + 2.0;
1128 q = (y[i + 1] - 2.0 * y[i] + y[i - 1]) / dx;
1129 u[i] = (3.0 * q / dx - 0.5 * u[i - 1]) / p;
1134 for (i = n - 2; i >= 0; i--)
1136 y2[i] = y2[i] * y2[i + 1] + u[i];
1142 interpolate1d(double xmin, double dx, const double* ya, const double* y2a, double x, double* y, double* y1)
1147 ix = static_cast<int>((x - xmin) / dx);
1149 a = (xmin + (ix + 1) * dx - x) / dx;
1150 b = (x - xmin - ix * dx) / dx;
1152 *y = a * ya[ix] + b * ya[ix + 1]
1153 + ((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]
1155 + (3.0 * b * b - 1.0) / 6.0 * dx * y2a[ix + 1];
1159 static void setup_cmap(int grid_spacing, int nc, gmx::ArrayRef<const real> grid, gmx_cmap_t* cmap_grid)
1161 int i, j, k, ii, jj, kk, idx;
1163 double dx, xmin, v, v1, v2, v12;
1166 std::vector<double> tmp_u(2 * grid_spacing, 0.0);
1167 std::vector<double> tmp_u2(2 * grid_spacing, 0.0);
1168 std::vector<double> tmp_yy(2 * grid_spacing, 0.0);
1169 std::vector<double> tmp_y1(2 * grid_spacing, 0.0);
1170 std::vector<double> tmp_t2(2 * grid_spacing * 2 * grid_spacing, 0.0);
1171 std::vector<double> tmp_grid(2 * grid_spacing * 2 * grid_spacing, 0.0);
1173 dx = 360.0 / grid_spacing;
1174 xmin = -180.0 - dx * grid_spacing / 2;
1176 for (kk = 0; kk < nc; kk++)
1178 /* Compute an offset depending on which cmap we are using
1179 * Offset will be the map number multiplied with the
1180 * grid_spacing * grid_spacing * 2
1182 offset = kk * grid_spacing * grid_spacing * 2;
1184 for (i = 0; i < 2 * grid_spacing; i++)
1186 ii = (i + grid_spacing - grid_spacing / 2) % grid_spacing;
1188 for (j = 0; j < 2 * grid_spacing; j++)
1190 jj = (j + grid_spacing - grid_spacing / 2) % grid_spacing;
1191 tmp_grid[i * grid_spacing * 2 + j] = grid[offset + ii * grid_spacing + jj];
1195 for (i = 0; i < 2 * grid_spacing; i++)
1197 spline1d(dx, &(tmp_grid[2 * grid_spacing * i]), 2 * grid_spacing, tmp_u.data(),
1198 &(tmp_t2[2 * grid_spacing * i]));
1201 for (i = grid_spacing / 2; i < grid_spacing + grid_spacing / 2; i++)
1203 ii = i - grid_spacing / 2;
1204 phi = ii * dx - 180.0;
1206 for (j = grid_spacing / 2; j < grid_spacing + grid_spacing / 2; j++)
1208 jj = j - grid_spacing / 2;
1209 psi = jj * dx - 180.0;
1211 for (k = 0; k < 2 * grid_spacing; k++)
1213 interpolate1d(xmin, dx, &(tmp_grid[2 * grid_spacing * k]),
1214 &(tmp_t2[2 * grid_spacing * k]), psi, &tmp_yy[k], &tmp_y1[k]);
1217 spline1d(dx, tmp_yy.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1218 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1219 spline1d(dx, tmp_y1.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1220 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1222 idx = ii * grid_spacing + jj;
1223 cmap_grid->cmapdata[kk].cmap[idx * 4] = grid[offset + ii * grid_spacing + jj];
1224 cmap_grid->cmapdata[kk].cmap[idx * 4 + 1] = v1;
1225 cmap_grid->cmapdata[kk].cmap[idx * 4 + 2] = v2;
1226 cmap_grid->cmapdata[kk].cmap[idx * 4 + 3] = v12;
1232 static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
1236 cmap_grid->grid_spacing = grid_spacing;
1237 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
1239 cmap_grid->cmapdata.resize(ngrid);
1241 for (i = 0; i < ngrid; i++)
1243 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
1248 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
1250 int count, count_mol;
1254 for (const gmx_molblock_t& molb : mtop->molblock)
1257 gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1259 for (int i = 0; i < F_NRE; i++)
1263 count_mol += 3 * interactions[i].size();
1265 else if (interaction_function[i].flags & IF_CONSTRAINT)
1267 count_mol += interactions[i].size();
1271 if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1274 "Molecule type '%s' has %d constraints.\n"
1275 "For stability and efficiency there should not be more constraints than "
1276 "internal number of degrees of freedom: %d.\n",
1277 *mi[molb.type].name, count_mol, nrdf_internal(&mi[molb.type].atoms));
1280 count += molb.nmol * count_mol;
1286 static real calc_temp(const gmx_mtop_t* mtop, const t_inputrec* ir, rvec* v)
1289 for (const AtomProxy atomP : AtomRange(*mtop))
1291 const t_atom& local = atomP.atom();
1292 int i = atomP.globalAtomNumber();
1293 sum_mv2 += local.m * norm2(v[i]);
1297 for (int g = 0; g < ir->opts.ngtc; g++)
1299 nrdf += ir->opts.nrdf[g];
1302 return sum_mv2 / (nrdf * BOLTZ);
1305 static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
1313 for (i = 0; i < ir->opts.ngtc; i++)
1315 if (ir->opts.tau_t[i] < 0)
1321 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1330 "Some temperature coupling groups do not use temperature coupling. We will assume "
1331 "their temperature is not more than %.3f K. If their temperature is higher, the "
1332 "energy error and the Verlet buffer might be underestimated.",
1340 /* Checks if there are unbound atoms in moleculetype molt.
1341 * Prints a note for each unbound atoms and a warning if any is present.
1343 static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi)
1345 const t_atoms* atoms = &molt->atoms;
1349 /* Only one atom, there can't be unbound atoms */
1353 std::vector<int> count(atoms->nr, 0);
1355 for (int ftype = 0; ftype < F_NRE; ftype++)
1357 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS)
1358 || (interaction_function[ftype].flags & IF_CONSTRAINT) || ftype == F_SETTLE)
1360 const InteractionList& il = molt->ilist[ftype];
1361 const int nral = NRAL(ftype);
1363 for (int i = 0; i < il.size(); i += 1 + nral)
1365 for (int j = 0; j < nral; j++)
1367 count[il.iatoms[i + 1 + j]]++;
1373 int numDanglingAtoms = 0;
1374 for (int a = 0; a < atoms->nr; a++)
1376 if (atoms->atom[a].ptype != eptVSite && count[a] == 0)
1381 "\nAtom %d '%s' in moleculetype '%s' is not bound by a potential or "
1382 "constraint to any other atom in the same moleculetype.\n",
1383 a + 1, *atoms->atomname[a], *molt->name);
1389 if (numDanglingAtoms > 0)
1393 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
1394 "other atom in the same moleculetype. Although technically this might not cause "
1395 "issues in a simulation, this often means that the user forgot to add a "
1396 "bond/potential/constraint or put multiple molecules in the same moleculetype "
1397 "definition by mistake. Run with -v to get information for each atom.",
1398 *molt->name, numDanglingAtoms);
1399 warning_note(wi, buf);
1403 /* Checks all moleculetypes for unbound atoms */
1404 static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi)
1406 for (const gmx_moltype_t& molt : mtop->moltype)
1408 checkForUnboundAtoms(&molt, bVerbose, wi);
1412 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1414 * The specific decoupled modes this routine check for are angle modes
1415 * where the two bonds are constrained and the atoms a both ends are only
1416 * involved in a single constraint; the mass of the two atoms needs to
1417 * differ by more than \p massFactorThreshold.
1419 static bool haveDecoupledModeInMol(const gmx_moltype_t& molt,
1420 gmx::ArrayRef<const t_iparams> iparams,
1421 real massFactorThreshold)
1423 if (molt.ilist[F_CONSTR].size() == 0 && molt.ilist[F_CONSTRNC].size() == 0)
1428 const t_atom* atom = molt.atoms.atom;
1430 const auto atomToConstraints =
1431 gmx::make_at2con(molt, iparams, gmx::FlexibleConstraintTreatment::Exclude);
1433 bool haveDecoupledMode = false;
1434 for (int ftype = 0; ftype < F_NRE; ftype++)
1436 if (interaction_function[ftype].flags & IF_ATYPE)
1438 const int nral = NRAL(ftype);
1439 const InteractionList& il = molt.ilist[ftype];
1440 for (int i = 0; i < il.size(); i += 1 + nral)
1442 /* Here we check for the mass difference between the atoms
1443 * at both ends of the angle, that the atoms at the ends have
1444 * 1 contraint and the atom in the middle at least 3; we check
1445 * that the 3 atoms are linked by constraints below.
1446 * We check for at least three constraints for the middle atom,
1447 * since with only the two bonds in the angle, we have 3-atom
1448 * molecule, which has much more thermal exhange in this single
1449 * angle mode than molecules with more atoms.
1450 * Note that this check also catches molecules where more atoms
1451 * are connected to one or more atoms in the angle, but by
1452 * bond potentials instead of angles. But such cases will not
1453 * occur in "normal" molecules and it doesn't hurt running
1454 * those with higher accuracy settings as well.
1456 int a0 = il.iatoms[1 + i];
1457 int a1 = il.iatoms[1 + i + 1];
1458 int a2 = il.iatoms[1 + i + 2];
1459 if ((atom[a0].m > atom[a2].m * massFactorThreshold || atom[a2].m > atom[a0].m * massFactorThreshold)
1460 && atomToConstraints[a0].ssize() == 1 && atomToConstraints[a2].ssize() == 1
1461 && atomToConstraints[a1].ssize() >= 3)
1463 int constraint0 = atomToConstraints[a0][0];
1464 int constraint2 = atomToConstraints[a2][0];
1466 bool foundAtom0 = false;
1467 bool foundAtom2 = false;
1468 for (const int constraint : atomToConstraints[a1])
1470 if (constraint == constraint0)
1474 if (constraint == constraint2)
1479 if (foundAtom0 && foundAtom2)
1481 haveDecoupledMode = true;
1488 return haveDecoupledMode;
1491 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1493 * When decoupled modes are present and the accuray in insufficient,
1494 * this routine issues a warning if the accuracy is insufficient.
1496 static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec* ir, warninp* wi)
1498 /* We only have issues with decoupled modes with normal MD.
1499 * With stochastic dynamics equipartitioning is enforced strongly.
1506 /* When atoms of very different mass are involved in an angle potential
1507 * and both bonds in the angle are constrained, the dynamic modes in such
1508 * angles have very different periods and significant energy exchange
1509 * takes several nanoseconds. Thus even a small amount of error in
1510 * different algorithms can lead to violation of equipartitioning.
1511 * The parameters below are mainly based on an all-atom chloroform model
1512 * with all bonds constrained. Equipartitioning is off by more than 1%
1513 * (need to run 5-10 ns) when the difference in mass between H and Cl
1514 * is more than a factor 13 and the accuracy is less than the thresholds
1515 * given below. This has been verified on some other molecules.
1517 * Note that the buffer and shake parameters have unit length and
1518 * energy/time, respectively, so they will "only" work correctly
1519 * for atomistic force fields using MD units.
1521 const real massFactorThreshold = 13.0;
1522 const real bufferToleranceThreshold = 1e-4;
1523 const int lincsIterationThreshold = 2;
1524 const int lincsOrderThreshold = 4;
1525 const real shakeToleranceThreshold = 0.005 * ir->delta_t;
1527 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold
1528 && ir->nProjOrder >= lincsOrderThreshold);
1529 bool shakeWithSufficientTolerance =
1530 (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1 * shakeToleranceThreshold);
1531 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol <= 1.1 * bufferToleranceThreshold
1532 && (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1537 bool haveDecoupledMode = false;
1538 for (const gmx_moltype_t& molt : mtop->moltype)
1540 if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams, massFactorThreshold))
1542 haveDecoupledMode = true;
1546 if (haveDecoupledMode)
1548 std::string message = gmx::formatString(
1549 "There are atoms at both ends of an angle, connected by constraints "
1550 "and with masses that differ by more than a factor of %g. This means "
1551 "that there are likely dynamic modes that are only very weakly coupled.",
1552 massFactorThreshold);
1553 if (ir->cutoff_scheme == ecutsVERLET)
1555 message += gmx::formatString(
1556 " To ensure good equipartitioning, you need to either not use "
1557 "constraints on all bonds (but, if possible, only on bonds involving "
1558 "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1559 "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1560 ">= %d or SHAKE tolerance <= %g",
1561 ei_names[eiSD1], bufferToleranceThreshold, lincsIterationThreshold,
1562 lincsOrderThreshold, shakeToleranceThreshold);
1566 message += gmx::formatString(
1567 " To ensure good equipartitioning, we suggest to switch to the %s "
1568 "cutoff-scheme, since that allows for better control over the Verlet "
1569 "buffer size and thus over the energy drift.",
1570 ecutscheme_names[ecutsVERLET]);
1572 warning(wi, message);
1576 static void set_verlet_buffer(const gmx_mtop_t* mtop, t_inputrec* ir, real buffer_temp, matrix box, warninp* wi)
1578 char warn_buf[STRLEN];
1580 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol,
1583 /* Calculate the buffer size for simple atom vs atoms list */
1584 VerletbufListSetup listSetup1x1;
1585 listSetup1x1.cluster_size_i = 1;
1586 listSetup1x1.cluster_size_j = 1;
1587 const real rlist_1x1 = calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1588 buffer_temp, listSetup1x1);
1590 /* Set the pair-list buffer size in ir */
1591 VerletbufListSetup listSetup4x4 = verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1592 ir->rlist = calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1,
1593 buffer_temp, listSetup4x4);
1595 const int n_nonlin_vsite = countNonlinearVsites(*mtop);
1596 if (n_nonlin_vsite > 0)
1599 "There are %d non-linear virtual site constructions. Their contribution to the "
1600 "energy error is approximated. In most cases this does not affect the error "
1603 warning_note(wi, warn_buf);
1606 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n", 1, 1,
1607 rlist_1x1, rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
1609 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1610 listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j, ir->rlist,
1611 ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
1613 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1615 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1618 "The pair-list cut-off (%g nm) is longer than half the shortest box vector or "
1619 "longer than the smallest box diagonal element (%g nm). Increase the box size or "
1620 "decrease nstlist or increase verlet-buffer-tolerance.",
1621 ir->rlist, std::sqrt(max_cutoff2(ir->ePBC, box)));
1625 int gmx_grompp(int argc, char* argv[])
1627 const char* desc[] = {
1628 "[THISMODULE] (the gromacs preprocessor)",
1629 "reads a molecular topology file, checks the validity of the",
1630 "file, expands the topology from a molecular description to an atomic",
1631 "description. The topology file contains information about",
1632 "molecule types and the number of molecules, the preprocessor",
1633 "copies each molecule as needed. ",
1634 "There is no limitation on the number of molecule types. ",
1635 "Bonds and bond-angles can be converted into constraints, separately",
1636 "for hydrogens and heavy atoms.",
1637 "Then a coordinate file is read and velocities can be generated",
1638 "from a Maxwellian distribution if requested.",
1639 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1640 "(eg. number of MD steps, time step, cut-off), and others such as",
1641 "NEMD parameters, which are corrected so that the net acceleration",
1643 "Eventually a binary file is produced that can serve as the sole input",
1644 "file for the MD program.[PAR]",
1646 "[THISMODULE] uses the atom names from the topology file. The atom names",
1647 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1648 "warnings when they do not match the atom names in the topology.",
1649 "Note that the atom names are irrelevant for the simulation as",
1650 "only the atom types are used for generating interaction parameters.[PAR]",
1652 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1653 "etc. The preprocessor supports the following keywords::",
1656 " #ifndef VARIABLE",
1659 " #define VARIABLE",
1661 " #include \"filename\"",
1662 " #include <filename>",
1664 "The functioning of these statements in your topology may be modulated by",
1665 "using the following two flags in your [REF].mdp[ref] file::",
1667 " define = -DVARIABLE1 -DVARIABLE2",
1668 " include = -I/home/john/doe",
1670 "For further information a C-programming textbook may help you out.",
1671 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1672 "topology file written out so that you can verify its contents.[PAR]",
1674 "When using position restraints, a file with restraint coordinates",
1675 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1676 "for [TT]-c[tt]). For free energy calculations, separate reference",
1677 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1678 "otherwise they will be equal to those of the A topology.[PAR]",
1680 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1681 "The last frame with coordinates and velocities will be read,",
1682 "unless the [TT]-time[tt] option is used. Only if this information",
1683 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1684 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1685 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1686 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1689 "[THISMODULE] can be used to restart simulations (preserving",
1690 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1691 "However, for simply changing the number of run steps to extend",
1692 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1693 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1694 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1695 "like output frequency, then supplying the checkpoint file to",
1696 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1697 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1698 "the ensemble (if possible) still requires passing the checkpoint",
1699 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1701 "By default, all bonded interactions which have constant energy due to",
1702 "virtual site constructions will be removed. If this constant energy is",
1703 "not zero, this will result in a shift in the total energy. All bonded",
1704 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1705 "all constraints for distances which will be constant anyway because",
1706 "of virtual site constructions will be removed. If any constraints remain",
1707 "which involve virtual sites, a fatal error will result.[PAR]",
1709 "To verify your run input file, please take note of all warnings",
1710 "on the screen, and correct where necessary. Do also look at the contents",
1711 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1712 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1713 "with the [TT]-debug[tt] option which will give you more information",
1714 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1715 "can see the contents of the run input file with the [gmx-dump]",
1716 "program. [gmx-check] can be used to compare the contents of two",
1717 "run input files.[PAR]",
1719 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1720 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1721 "harmless, but usually they are not. The user is advised to carefully",
1722 "interpret the output messages before attempting to bypass them with",
1726 std::vector<MoleculeInformation> mi;
1727 std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1731 const char* mdparin;
1732 bool bNeedVel, bGenVel;
1733 gmx_bool have_atomnumber;
1734 gmx_output_env_t* oenv;
1735 gmx_bool bVerbose = FALSE;
1737 char warn_buf[STRLEN];
1739 t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
1740 { efMDP, "-po", "mdout", ffWRITE },
1741 { efSTX, "-c", nullptr, ffREAD },
1742 { efSTX, "-r", "restraint", ffOPTRD },
1743 { efSTX, "-rb", "restraint", ffOPTRD },
1744 { efNDX, nullptr, nullptr, ffOPTRD },
1745 { efTOP, nullptr, nullptr, ffREAD },
1746 { efTOP, "-pp", "processed", ffOPTWR },
1747 { efTPR, "-o", nullptr, ffWRITE },
1748 { efTRN, "-t", nullptr, ffOPTRD },
1749 { efEDR, "-e", nullptr, ffOPTRD },
1750 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1751 { efGRO, "-imd", "imdgroup", ffOPTWR },
1752 { efTRN, "-ref", "rotref", ffOPTRW } };
1753 #define NFILE asize(fnm)
1755 /* Command line options */
1756 gmx_bool bRenum = TRUE;
1757 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1761 { "-v", FALSE, etBOOL, { &bVerbose }, "Be loud and noisy" },
1762 { "-time", FALSE, etREAL, { &fr_time }, "Take frame at or first after this time." },
1767 "Remove constant bonded interactions with virtual sites" },
1772 "Number of allowed warnings during input processing. Not for normal use and may "
1773 "generate unstable systems" },
1778 "Set parameters for bonded interactions without defaults to zero instead of "
1779 "generating an error" },
1784 "Renumber atomtypes and minimize number of atomtypes" }
1787 /* Parse the command line */
1788 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
1793 /* Initiate some variables */
1794 gmx::MDModules mdModules;
1795 t_inputrec irInstance;
1796 t_inputrec* ir = &irInstance;
1798 snew(opts->include, STRLEN);
1799 snew(opts->define, STRLEN);
1801 wi = init_warning(TRUE, maxwarn);
1803 /* PARAMETER file processing */
1804 mdparin = opt2fn("-f", NFILE, fnm);
1805 set_warning_line(wi, mdparin, -1);
1808 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1810 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1814 fprintf(stderr, "checking input for internal consistency...\n");
1816 check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
1818 if (ir->ld_seed == -1)
1820 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1821 fprintf(stderr, "Setting the LD random seed to %" PRId64 "\n", ir->ld_seed);
1824 if (ir->expandedvals->lmc_seed == -1)
1826 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1827 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1830 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1831 bGenVel = (bNeedVel && opts->bGenVel);
1832 if (bGenVel && ir->bContinuation)
1835 "Generating velocities is inconsistent with attempting "
1836 "to continue a previous run. Choose only one of "
1837 "gen-vel = yes and continuation = yes.");
1838 warning_error(wi, warn_buf);
1841 std::array<InteractionsOfType, F_NRE> interactions;
1843 PreprocessingAtomTypes atypes;
1846 pr_symtab(debug, 0, "Just opened", &sys.symtab);
1849 const char* fn = ftp2fn(efTOP, NFILE, fnm);
1850 if (!gmx_fexist(fn))
1852 gmx_fatal(FARGS, "%s does not exist", fn);
1856 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm), opts, ir, bZero,
1857 bGenVel, bVerbose, &state, &atypes, &sys, &mi, &intermolecular_interactions,
1858 interactions, &comb, &reppow, &fudgeQQ, opts->bMorse, wi);
1862 pr_symtab(debug, 0, "After new_status", &sys.symtab);
1866 /* set parameters for virtual site construction (not for vsiten) */
1867 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1869 nvsite += set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions);
1871 /* now throw away all obsolete bonds, angles and dihedrals: */
1872 /* note: constraints are ALWAYS removed */
1875 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
1877 clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds);
1881 if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == econtSHAKE))
1883 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1885 sprintf(warn_buf, "Can not do %s with %s, use %s", EI(ir->eI),
1886 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1887 warning_error(wi, warn_buf);
1889 if (ir->bPeriodicMols)
1891 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1892 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1893 warning_error(wi, warn_buf);
1897 if (EI_SD(ir->eI) && ir->etc != etcNO)
1899 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1902 /* If we are doing QM/MM, check that we got the atom numbers */
1903 have_atomnumber = TRUE;
1904 for (gmx::index i = 0; i < gmx::ssize(atypes); i++)
1906 have_atomnumber = have_atomnumber && (atypes.atomNumberFromAtomType(i) >= 0);
1908 if (!have_atomnumber && ir->bQMMM)
1913 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1914 "field you are using does not contain atom numbers fields. This is an\n"
1915 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1916 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1917 "an integer just before the mass column in ffXXXnb.itp.\n"
1918 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1921 /* Check for errors in the input now, since they might cause problems
1922 * during processing further down.
1924 check_warning_error(wi, FARGS);
1926 if (nint_ftype(&sys, mi, F_POSRES) > 0 || nint_ftype(&sys, mi, F_FBPOSRES) > 0)
1928 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1931 "You are combining position restraints with %s pressure coupling, which can "
1932 "lead to instabilities. If you really want to combine position restraints with "
1933 "pressure coupling, we suggest to use %s pressure coupling instead.",
1934 EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
1935 warning_note(wi, warn_buf);
1938 const char* fn = opt2fn("-r", NFILE, fnm);
1941 if (!gmx_fexist(fn))
1944 "Cannot find position restraint file %s (option -r).\n"
1945 "From GROMACS-2018, you need to specify the position restraint "
1946 "coordinate files explicitly to avoid mistakes, although you can "
1947 "still use the same file as you specify for the -c option.",
1951 if (opt2bSet("-rb", NFILE, fnm))
1953 fnB = opt2fn("-rb", NFILE, fnm);
1954 if (!gmx_fexist(fnB))
1957 "Cannot find B-state position restraint file %s (option -rb).\n"
1958 "From GROMACS-2018, you need to specify the position restraint "
1959 "coordinate files explicitly to avoid mistakes, although you can "
1960 "still use the same file as you specify for the -c option.",
1971 fprintf(stderr, "Reading position restraint coords from %s", fn);
1972 if (strcmp(fn, fnB) == 0)
1974 fprintf(stderr, "\n");
1978 fprintf(stderr, " and %s\n", fnB);
1981 gen_posres(&sys, mi, fn, fnB, ir->refcoord_scaling, ir->ePBC, ir->posres_com, ir->posres_comB, wi);
1984 /* If we are using CMAP, setup the pre-interpolation grid */
1985 if (interactions[F_CMAP].ncmap() > 0)
1987 init_cmap_grid(&sys.ffparams.cmap_grid, interactions[F_CMAP].cmapAngles,
1988 interactions[F_CMAP].cmakeGridSpacing);
1989 setup_cmap(interactions[F_CMAP].cmakeGridSpacing, interactions[F_CMAP].cmapAngles,
1990 interactions[F_CMAP].cmap, &sys.ffparams.cmap_grid);
1993 set_wall_atomtype(&atypes, opts, ir, wi);
1996 atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
1999 if (ir->implicit_solvent)
2001 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
2004 /* PELA: Copy the atomtype data to the topology atomtype list */
2005 atypes.copyTot_atomtypes(&(sys.atomtypes));
2009 pr_symtab(debug, 0, "After atype.renumberTypes", &sys.symtab);
2014 fprintf(stderr, "converting bonded parameters...\n");
2017 const int ntype = atypes.size();
2018 convertInteractionsOfType(ntype, interactions, mi, intermolecular_interactions.get(), comb,
2019 reppow, fudgeQQ, &sys);
2023 pr_symtab(debug, 0, "After converInteractionsOfType", &sys.symtab);
2026 /* set ptype to VSite for virtual sites */
2027 for (gmx_moltype_t& moltype : sys.moltype)
2029 set_vsites_ptype(FALSE, &moltype);
2033 pr_symtab(debug, 0, "After virtual sites", &sys.symtab);
2035 /* Check velocity for virtual sites and shells */
2038 check_vel(&sys, state.v.rvec_array());
2041 /* check for shells and inpurecs */
2042 check_shells_inputrec(&sys, ir, wi);
2045 check_mol(&sys, wi);
2047 checkForUnboundAtoms(&sys, bVerbose, wi);
2049 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2051 check_bonds_timestep(&sys, ir->delta_t, wi);
2054 checkDecoupledModeAccuracy(&sys, ir, wi);
2056 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2060 "Zero-step energy minimization will alter the coordinates before calculating the "
2061 "energy. If you just want the energy of a single point, try zero-step MD (with "
2062 "unconstrained_start = yes). To do multiple single-point energy evaluations of "
2063 "different configurations of the same topology, use mdrun -rerun.");
2066 check_warning_error(wi, FARGS);
2070 fprintf(stderr, "initialising group options...\n");
2072 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifier(), ir, wi);
2074 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2076 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2080 if (EI_MD(ir->eI) && ir->etc == etcNO)
2084 buffer_temp = opts->tempi;
2088 buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2090 if (buffer_temp > 0)
2093 "NVE simulation: will use the initial temperature of %.3f K for "
2094 "determining the Verlet buffer size",
2096 warning_note(wi, warn_buf);
2101 "NVE simulation with an initial temperature of zero: will use a Verlet "
2102 "buffer of %d%%. Check your energy drift!",
2103 gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
2104 warning_note(wi, warn_buf);
2109 buffer_temp = get_max_reference_temp(ir, wi);
2112 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2114 /* NVE with initial T=0: we add a fixed ratio to rlist.
2115 * Since we don't actually use verletbuf_tol, we set it to -1
2116 * so it can't be misused later.
2118 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2119 ir->verletbuf_tol = -1;
2123 /* We warn for NVE simulations with a drift tolerance that
2124 * might result in a 1(.1)% drift over the total run-time.
2125 * Note that we can't warn when nsteps=0, since we don't
2126 * know how many steps the user intends to run.
2128 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 && ir->nsteps > 0)
2130 const real driftTolerance = 0.01;
2131 /* We use 2 DOF per atom = 2kT pot+kin energy,
2132 * to be on the safe side with constraints.
2134 const real totalEnergyDriftPerAtomPerPicosecond =
2135 2 * BOLTZ * buffer_temp / (ir->nsteps * ir->delta_t);
2137 if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
2140 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
2141 "NVE simulation of length %g ps, which can give a final drift of "
2142 "%d%%. For conserving energy to %d%% when using constraints, you "
2143 "might need to set verlet-buffer-tolerance to %.1e.",
2144 ir->verletbuf_tol, ir->nsteps * ir->delta_t,
2145 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
2146 gmx::roundToInt(100 * driftTolerance),
2147 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
2148 warning_note(wi, warn_buf);
2152 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi);
2157 /* Init the temperature coupling state */
2158 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2162 pr_symtab(debug, 0, "After index", &sys.symtab);
2165 triple_check(mdparin, ir, &sys, wi);
2166 close_symtab(&sys.symtab);
2169 pr_symtab(debug, 0, "After close", &sys.symtab);
2172 /* make exclusions between QM atoms and remove charges if needed */
2175 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_ORIGINAL);
2176 if (ir->QMMMscheme != eQMMMschemeoniom)
2178 std::vector<int> qmmmAtoms = qmmmAtomIndices(*ir, sys);
2179 removeQmmmAtomCharges(&sys, qmmmAtoms);
2183 if (ir->eI == eiMimic)
2185 generate_qmexcl(&sys, ir, wi, GmxQmmmMode::GMX_QMMM_MIMIC);
2188 if (ftp2bSet(efTRN, NFILE, fnm))
2192 fprintf(stderr, "getting data from old trajectory ...\n");
2194 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm), bNeedVel, bGenVel,
2195 fr_time, ir, &state, &sys, oenv);
2198 if (ir->ePBC == epbcXY && ir->nwall != 2)
2200 clear_rvec(state.box[ZZ]);
2203 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2205 /* Calculate the optimal grid dimensions */
2207 EwaldBoxZScaler boxScaler(*ir);
2208 boxScaler.scaleBox(state.box, scaledBox);
2210 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2212 /* Mark fourier_spacing as not used */
2213 ir->fourier_spacing = 0;
2215 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2217 set_warning_line(wi, mdparin, -1);
2219 wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2221 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2222 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize, &(ir->nkx), &(ir->nky),
2224 if (ir->nkx < minGridSize || ir->nky < minGridSize || ir->nkz < minGridSize)
2227 "The PME grid size should be >= 2*(pme-order - 1); either manually "
2228 "increase the grid size or decrease pme-order");
2232 /* MRS: eventually figure out better logic for initializing the fep
2233 values that makes declaring the lambda and declaring the state not
2234 potentially conflict if not handled correctly. */
2235 if (ir->efep != efepNO)
2237 state.fep_state = ir->fepvals->init_fep_state;
2238 for (i = 0; i < efptNR; i++)
2240 /* init_lambda trumps state definitions*/
2241 if (ir->fepvals->init_lambda >= 0)
2243 state.lambda[i] = ir->fepvals->init_lambda;
2247 if (ir->fepvals->all_lambda[i] == nullptr)
2249 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2253 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2259 pull_t* pull = nullptr;
2263 pull = set_pull_init(ir, &sys, state.x.rvec_array(), state.box, state.lambda[efptMASS], wi);
2266 /* Modules that supply external potential for pull coordinates
2267 * should register those potentials here. finish_pull() will check
2268 * that providers have been registerd for all external potentials.
2273 tensor compressibility = { { 0 } };
2274 if (ir->epc != epcNO)
2276 copy_mat(ir->compress, compressibility);
2278 setStateDependentAwhParams(ir->awhParams, ir->pull, pull, state.box, ir->ePBC,
2279 compressibility, &ir->opts, wi);
2289 set_reference_positions(ir->rot, state.x.rvec_array(), state.box,
2290 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm), wi);
2293 /* reset_multinr(sys); */
2295 if (EEL_PME(ir->coulombtype))
2297 float ratio = pme_load_estimate(sys, *ir, state.box);
2298 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2299 /* With free energy we might need to do PME both for the A and B state
2300 * charges. This will double the cost, but the optimal performance will
2301 * then probably be at a slightly larger cut-off and grid spacing.
2303 if ((ir->efep == efepNO && ratio > 1.0 / 2.0) || (ir->efep != efepNO && ratio > 2.0 / 3.0))
2307 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2308 "and for highly parallel simulations between 0.25 and 0.33,\n"
2309 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2310 if (ir->efep != efepNO)
2313 "For free energy simulations, the optimal load limit increases from "
2320 char warn_buf[STRLEN];
2321 double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2322 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2325 set_warning_line(wi, mdparin, -1);
2326 warning_note(wi, warn_buf);
2330 printf("%s\n", warn_buf);
2334 // Add the md modules internal parameters that are not mdp options
2335 // e.g., atom indices
2338 gmx::KeyValueTreeBuilder internalParameterBuilder;
2339 mdModules.notifier().notifier_.notify(internalParameterBuilder.rootObject());
2340 ir->internalParameters =
2341 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2346 fprintf(stderr, "writing run input file...\n");
2349 done_warning(wi, FARGS);
2350 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, &sys);
2352 /* Output IMD group, if bIMD is TRUE */
2353 gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, &sys, NFILE, fnm);
2355 sfree(opts->define);
2356 sfree(opts->include);
2358 for (auto& mol : mi)
2360 // Some of the contents of molinfo have been stolen, so
2361 // fullCleanUp can't be called.
2362 mol.partialCleanUp();
2364 done_inputrec_strings();
2365 output_env_done(oenv);