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,2021, 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.
51 #include <sys/types.h>
53 #include "gromacs/applied_forces/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/md_support.h"
83 #include "gromacs/mdlib/perf_est.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdrun/mdmodules.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/mdtypes/md_enums.h"
88 #include "gromacs/mdtypes/nblist.h"
89 #include "gromacs/mdtypes/state.h"
90 #include "gromacs/pbcutil/boxutilities.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/pulling/pull.h"
93 #include "gromacs/random/seed.h"
94 #include "gromacs/topology/ifunc.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/symtab.h"
97 #include "gromacs/topology/topology.h"
98 #include "gromacs/trajectory/trajectoryframe.h"
99 #include "gromacs/utility/arraysize.h"
100 #include "gromacs/utility/cstringutil.h"
101 #include "gromacs/utility/exceptions.h"
102 #include "gromacs/utility/fatalerror.h"
103 #include "gromacs/utility/filestream.h"
104 #include "gromacs/utility/futil.h"
105 #include "gromacs/utility/gmxassert.h"
106 #include "gromacs/utility/keyvaluetreebuilder.h"
107 #include "gromacs/utility/listoflists.h"
108 #include "gromacs/utility/logger.h"
109 #include "gromacs/utility/loggerbuilder.h"
110 #include "gromacs/utility/mdmodulenotification.h"
111 #include "gromacs/utility/smalloc.h"
112 #include "gromacs/utility/snprintf.h"
114 /* TODO The implementation details should move to their own source file. */
115 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int> atoms,
116 gmx::ArrayRef<const real> params,
117 const std::string& name) :
118 atoms_(atoms.begin(), atoms.end()),
119 interactionTypeName_(name)
122 params.size() <= forceParam_.size(),
123 gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM)
125 auto forceParamIt = forceParam_.begin();
126 for (const auto param : params)
128 *forceParamIt++ = param;
130 while (forceParamIt != forceParam_.end())
132 *forceParamIt++ = NOTSET;
136 const int& InteractionOfType::ai() const
138 GMX_RELEASE_ASSERT(!atoms_.empty(), "Need to have at least one atom set");
142 const int& InteractionOfType::aj() const
144 GMX_RELEASE_ASSERT(atoms_.size() > 1, "Need to have at least two atoms set");
148 const int& InteractionOfType::ak() const
150 GMX_RELEASE_ASSERT(atoms_.size() > 2, "Need to have at least three atoms set");
154 const int& InteractionOfType::al() const
156 GMX_RELEASE_ASSERT(atoms_.size() > 3, "Need to have at least four atoms set");
160 const int& InteractionOfType::am() const
162 GMX_RELEASE_ASSERT(atoms_.size() > 4, "Need to have at least five atoms set");
166 const real& InteractionOfType::c0() const
168 return forceParam_[0];
171 const real& InteractionOfType::c1() const
173 return forceParam_[1];
176 const real& InteractionOfType::c2() const
178 return forceParam_[2];
181 const std::string& InteractionOfType::interactionTypeName() const
183 return interactionTypeName_;
186 void InteractionOfType::sortBondAtomIds()
190 std::swap(atoms_[0], atoms_[1]);
194 void InteractionOfType::sortAngleAtomIds()
198 std::swap(atoms_[0], atoms_[2]);
202 void InteractionOfType::sortDihedralAtomIds()
206 std::swap(atoms_[0], atoms_[3]);
207 std::swap(atoms_[1], atoms_[2]);
211 void InteractionOfType::sortAtomIds()
221 else if (isDihedral())
223 sortDihedralAtomIds();
227 GMX_THROW(gmx::InternalError(
228 "Can not sort parameters other than bonds, angles or dihedrals"));
232 void InteractionOfType::setForceParameter(int pos, real value)
234 GMX_RELEASE_ASSERT(pos < MAXFORCEPARAM,
235 "Can't set parameter beyond the maximum number of parameters");
236 forceParam_[pos] = value;
239 void MoleculeInformation::initMolInfo()
243 init_t_atoms(&atoms, 0, FALSE);
246 void MoleculeInformation::partialCleanUp()
251 void MoleculeInformation::fullCleanUp()
257 static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
260 /* For all the molecule types */
261 for (auto& mol : mols)
263 n += mol.interactions[ifunc].size();
264 mol.interactions[ifunc].interactionTypes.clear();
269 static int check_atom_names(const char* fn1,
273 const gmx::MDLogger& logger)
275 int m, i, j, nmismatch;
278 constexpr int c_maxNumberOfMismatches = 20;
280 if (mtop->natoms != at->nr)
282 gmx_incons("comparing atom names");
287 for (const gmx_molblock_t& molb : mtop->molblock)
289 tat = &mtop->moltype[molb.type].atoms;
290 for (m = 0; m < molb.nmol; m++)
292 for (j = 0; j < tat->nr; j++)
294 if (strcmp(*(tat->atomname[j]), *(at->atomname[i])) != 0)
296 if (nmismatch < c_maxNumberOfMismatches)
298 GMX_LOG(logger.warning)
300 .appendTextFormatted(
301 "atom name %d in %s and %s does not match (%s - %s)",
308 else if (nmismatch == c_maxNumberOfMismatches)
310 GMX_LOG(logger.warning)
312 .appendTextFormatted("(more than %d non-matching atom names)",
313 c_maxNumberOfMismatches);
325 static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
327 /* This check is not intended to ensure accurate integration,
328 * rather it is to signal mistakes in the mdp settings.
329 * A common mistake is to forget to turn on constraints
330 * for MD after energy minimization with flexible bonds.
331 * This check can also detect too large time steps for flexible water
332 * models, but such errors will often be masked by the constraints
333 * mdp options, which turns flexible water into water with bond constraints,
334 * but without an angle constraint. Unfortunately such incorrect use
335 * of water models can not easily be detected without checking
336 * for specific model names.
338 * The stability limit of leap-frog or velocity verlet is 4.44 steps
339 * per oscillational period.
340 * But accurate bonds distributions are lost far before that limit.
341 * To allow relatively common schemes (although not common with Gromacs)
342 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
343 * we set the note limit to 10.
345 int min_steps_warn = 5;
346 int min_steps_note = 10;
348 int i, a1, a2, w_a1, w_a2, j;
349 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
350 bool bFound, bWater, bWarn;
352 /* Get the interaction parameters */
353 gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
355 twopi2 = gmx::square(2 * M_PI);
357 limit2 = gmx::square(min_steps_note * dt);
362 const gmx_moltype_t* w_moltype = nullptr;
363 for (const gmx_moltype_t& moltype : mtop->moltype)
365 const t_atom* atom = moltype.atoms.atom;
366 const InteractionLists& ilist = moltype.ilist;
367 const InteractionList& ilc = ilist[F_CONSTR];
368 const InteractionList& ils = ilist[F_SETTLE];
369 for (ftype = 0; ftype < F_NRE; ftype++)
371 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
376 const InteractionList& ilb = ilist[ftype];
377 for (i = 0; i < ilb.size(); i += 3)
379 fc = ip[ilb.iatoms[i]].harmonic.krA;
380 re = ip[ilb.iatoms[i]].harmonic.rA;
381 if (ftype == F_G96BONDS)
383 /* Convert squared sqaure fc to harmonic fc */
386 a1 = ilb.iatoms[i + 1];
387 a2 = ilb.iatoms[i + 2];
390 if (fc > 0 && m1 > 0 && m2 > 0)
392 period2 = twopi2 * m1 * m2 / ((m1 + m2) * fc);
396 period2 = GMX_FLOAT_MAX;
400 fprintf(debug, "fc %g m1 %g m2 %g period %g\n", fc, m1, m2, std::sqrt(period2));
402 if (period2 < limit2)
405 for (j = 0; j < ilc.size(); j += 3)
407 if ((ilc.iatoms[j + 1] == a1 && ilc.iatoms[j + 2] == a2)
408 || (ilc.iatoms[j + 1] == a2 && ilc.iatoms[j + 2] == a1))
413 for (j = 0; j < ils.size(); j += 4)
415 if ((a1 == ils.iatoms[j + 1] || a1 == ils.iatoms[j + 2] || a1 == ils.iatoms[j + 3])
416 && (a2 == ils.iatoms[j + 1] || a2 == ils.iatoms[j + 2]
417 || a2 == ils.iatoms[j + 3]))
422 if (!bFound && (w_moltype == nullptr || period2 < w_period2))
424 w_moltype = &moltype;
434 if (w_moltype != nullptr)
436 bWarn = (w_period2 < gmx::square(min_steps_warn * dt));
437 /* A check that would recognize most water models */
438 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' && w_moltype->atoms.nr <= 5);
439 std::string warningMessage = gmx::formatString(
440 "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated "
441 "oscillational period of %.1e ps, which is less than %d times the time step of "
446 *w_moltype->atoms.atomname[w_a1],
448 *w_moltype->atoms.atomname[w_a2],
449 std::sqrt(w_period2),
450 bWarn ? min_steps_warn : min_steps_note,
452 bWater ? "Maybe you asked for fexible water."
453 : "Maybe you forgot to change the constraints mdp option.");
456 warning(wi, warningMessage.c_str());
460 warning_note(wi, warningMessage.c_str());
465 static void check_vel(gmx_mtop_t* mtop, rvec v[])
467 for (const AtomProxy atomP : AtomRange(*mtop))
469 const t_atom& local = atomP.atom();
470 int i = atomP.globalAtomNumber();
471 if (local.ptype == eptShell || local.ptype == eptBond || local.ptype == eptVSite)
478 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
482 for (const AtomProxy atomP : AtomRange(*mtop))
484 const t_atom& local = atomP.atom();
485 if (local.ptype == eptShell || local.ptype == eptBond)
490 if ((nshells > 0) && (ir->nstcalcenergy != 1))
492 set_warning_line(wi, "unknown", -1);
493 std::string warningMessage = gmx::formatString(
494 "There are %d shells, changing nstcalcenergy from %d to 1", nshells, ir->nstcalcenergy);
495 ir->nstcalcenergy = 1;
496 warning(wi, warningMessage.c_str());
500 /* TODO Decide whether this function can be consolidated with
501 * gmx_mtop_ftype_count */
502 static int nint_ftype(gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
505 for (const gmx_molblock_t& molb : mtop->molblock)
507 nint += molb.nmol * mi[molb.type].interactions[ftype].size();
513 /* This routine reorders the molecule type array
514 * in the order of use in the molblocks,
515 * unused molecule types are deleted.
517 static void renumber_moltypes(gmx_mtop_t* sys, std::vector<MoleculeInformation>* molinfo)
520 std::vector<int> order;
521 for (gmx_molblock_t& molblock : sys->molblock)
523 const auto found = std::find_if(order.begin(), order.end(), [&molblock](const auto& entry) {
524 return molblock.type == entry;
526 if (found == order.end())
528 /* This type did not occur yet, add it */
529 order.push_back(molblock.type);
530 molblock.type = order.size() - 1;
534 molblock.type = std::distance(order.begin(), found);
538 /* We still need to reorder the molinfo structs */
539 std::vector<MoleculeInformation> minew(order.size());
541 for (auto& mi : *molinfo)
543 const auto found = std::find(order.begin(), order.end(), index);
544 if (found != order.end())
546 int pos = std::distance(order.begin(), found);
551 // Need to manually clean up memory ....
560 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t* mtop)
562 mtop->moltype.resize(mi.size());
564 for (const auto& mol : mi)
566 gmx_moltype_t& molt = mtop->moltype[pos];
567 molt.name = mol.name;
568 molt.atoms = mol.atoms;
569 /* ilists are copied later */
570 molt.excls = mol.excls;
575 static void new_status(const char* topfile,
576 const char* topppfile,
584 PreprocessingAtomTypes* atypes,
586 std::vector<MoleculeInformation>* mi,
587 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
588 gmx::ArrayRef<InteractionsOfType> interactions,
594 const gmx::MDLogger& logger)
596 std::vector<gmx_molblock_t> molblock;
598 bool ffParametrizedWithHBondConstraints;
600 /* TOPOLOGY processing */
601 sys->name = do_top(bVerbose,
613 intermolecular_interactions,
616 &ffParametrizedWithHBondConstraints,
620 sys->molblock.clear();
623 for (const gmx_molblock_t& molb : molblock)
625 if (!sys->molblock.empty() && molb.type == sys->molblock.back().type)
627 /* Merge consecutive blocks with the same molecule type */
628 sys->molblock.back().nmol += molb.nmol;
630 else if (molb.nmol > 0)
632 /* Add a new molblock to the topology */
633 sys->molblock.push_back(molb);
635 sys->natoms += molb.nmol * (*mi)[sys->molblock.back().type].atoms.nr;
637 if (sys->molblock.empty())
639 gmx_fatal(FARGS, "No molecules were defined in the system");
642 renumber_moltypes(sys, mi);
646 convert_harmonics(*mi, atypes);
649 if (ir->eDisre == edrNone)
651 i = rm_interactions(F_DISRES, *mi);
654 set_warning_line(wi, "unknown", -1);
655 std::string warningMessage =
656 gmx::formatString("disre = no, removed %d distance restraints", i);
657 warning_note(wi, warningMessage.c_str());
662 i = rm_interactions(F_ORIRES, *mi);
665 set_warning_line(wi, "unknown", -1);
666 std::string warningMessage =
667 gmx::formatString("orire = no, removed %d orientation restraints", i);
668 warning_note(wi, warningMessage.c_str());
672 /* Copy structures from msys to sys */
673 molinfo2mtop(*mi, sys);
677 /* Discourage using the, previously common, setup of applying constraints
678 * to all bonds with force fields that have been parametrized with H-bond
679 * constraints only. Do not print note with large timesteps or vsites.
681 if (opts->nshake == eshALLBONDS && ffParametrizedWithHBondConstraints && ir->delta_t < 0.0026
682 && gmx_mtop_ftype_count(sys, F_VSITE3FD) == 0)
684 set_warning_line(wi, "unknown", -1);
686 "You are using constraints on all bonds, whereas the forcefield "
687 "has been parametrized only with constraints involving hydrogen atoms. "
688 "We suggest using constraints = h-bonds instead, this will also improve "
692 /* COORDINATE file processing */
695 GMX_LOG(logger.info).asParagraph().appendTextFormatted("processing coordinates...");
702 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
703 state->natoms = conftop->atoms.nr;
704 if (state->natoms != sys->natoms)
707 "number of coordinates in coordinate file (%s, %d)\n"
708 " does not match topology (%s, %d)",
714 /* It would be nice to get rid of the copies below, but we don't know
715 * a priori if the number of atoms in confin matches what we expect.
717 state->flags |= (1 << estX);
720 state->flags |= (1 << estV);
722 state_change_natoms(state, state->natoms);
723 std::copy(x, x + state->natoms, state->x.data());
727 std::copy(v, v + state->natoms, state->v.data());
730 /* This call fixes the box shape for runs with pressure scaling */
731 set_box_rel(ir, state);
733 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms, logger);
739 std::string warningMessage = gmx::formatString(
740 "%d non-matching atom name%s\n"
741 "atom names from %s will be used\n"
742 "atom names from %s will be ignored\n",
744 (nmismatch == 1) ? "" : "s",
747 warning(wi, warningMessage.c_str());
750 /* Do more checks, mostly related to constraints */
755 .appendTextFormatted("double-checking input for internal consistency...");
758 bool bHasNormalConstraints =
759 0 < (nint_ftype(sys, *mi, F_CONSTR) + nint_ftype(sys, *mi, F_CONSTRNC));
760 bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
761 double_check(ir, state->box, bHasNormalConstraints, bHasAnyConstraints, wi);
768 snew(mass, state->natoms);
769 for (const AtomProxy atomP : AtomRange(*sys))
771 const t_atom& local = atomP.atom();
772 int i = atomP.globalAtomNumber();
776 if (opts->seed == -1)
778 opts->seed = static_cast<int>(gmx::makeRandomSeed());
779 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Setting gen_seed to %d", opts->seed);
781 state->flags |= (1 << estV);
782 maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array(), logger);
784 stop_cm(logger, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
789 static void copy_state(const char* slog, t_trxframe* fr, bool bReadVel, t_state* state, double* use_time)
791 if (fr->not_ok & FRAME_NOT_OK)
793 gmx_fatal(FARGS, "Can not start from an incomplete frame");
797 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s", slog);
800 std::copy(fr->x, fr->x + state->natoms, state->x.data());
805 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
807 std::copy(fr->v, fr->v + state->natoms, state->v.data());
811 copy_mat(fr->box, state->box);
814 *use_time = fr->time;
817 static void cont_status(const char* slog,
825 const gmx_output_env_t* oenv,
826 const gmx::MDLogger& logger)
827 /* If fr_time == -1 read the last frame available which is complete */
834 bReadVel = (bNeedVel && !bGenVel);
838 .appendTextFormatted("Reading Coordinates%s and Box size from old trajectory",
839 bReadVel ? ", Velocities" : "");
842 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read whole trajectory");
846 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read till time %g", fr_time);
854 .appendTextFormatted(
855 "Velocities generated: "
856 "ignoring velocities in input trajectory");
858 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
862 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
866 GMX_LOG(logger.warning)
868 .appendTextFormatted(
869 "WARNING: Did not find a frame with velocities in file %s,\n"
870 " all velocities will be set to zero!",
872 for (auto& vi : makeArrayRef(state->v))
877 /* Search for a frame without velocities */
879 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
883 state->natoms = fr.natoms;
885 if (sys->natoms != state->natoms)
888 "Number of atoms in Topology "
889 "is not the same as in Trajectory");
891 copy_state(slog, &fr, bReadVel, state, &use_time);
893 /* Find the appropriate frame */
894 while ((fr_time == -1 || fr.time < fr_time) && read_next_frame(oenv, fp, &fr))
896 copy_state(slog, &fr, bReadVel, state, &use_time);
901 /* Set the relative box lengths for preserving the box shape.
902 * Note that this call can lead to differences in the last bit
903 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
905 set_box_rel(ir, state);
907 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Using frame at t = %g ps", use_time);
908 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Starting time for run is %g ps", ir->init_t);
910 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
912 get_enx_state(ener, use_time, sys->groups, ir, state);
913 preserve_box_shape(ir, state->box_rel, state->boxv);
917 static void read_posres(gmx_mtop_t* mtop,
918 gmx::ArrayRef<const MoleculeInformation> molinfo,
925 const gmx::MDLogger& logger)
933 int natoms, npbcdim = 0;
938 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
939 natoms = top->atoms.nr;
942 if (natoms != mtop->natoms)
944 std::string warningMessage = gmx::formatString(
945 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
946 "(%d). Will assume that the first %d atoms in the topology and %s match.",
950 std::min(mtop->natoms, natoms),
952 warning(wi, warningMessage.c_str());
955 npbcdim = numPbcDimensions(pbcType);
956 GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
958 if (rc_scaling != erscNO)
960 copy_mat(box, invbox);
961 for (int j = npbcdim; j < DIM; j++)
963 clear_rvec(invbox[j]);
966 gmx::invertBoxMatrix(invbox, invbox);
969 /* Copy the reference coordinates to mtop */
973 snew(hadAtom, natoms);
974 for (gmx_molblock_t& molb : mtop->molblock)
976 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
977 const InteractionsOfType* pr = &(molinfo[molb.type].interactions[F_POSRES]);
978 const InteractionsOfType* prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
979 if (pr->size() > 0 || prfb->size() > 0)
981 atom = mtop->moltype[molb.type].atoms.atom;
982 for (const auto& restraint : pr->interactionTypes)
984 int ai = restraint.ai();
988 "Position restraint atom index (%d) in moltype '%s' is larger than "
989 "number of atoms in %s (%d).\n",
991 *molinfo[molb.type].name,
996 if (rc_scaling == erscCOM)
998 /* Determine the center of mass of the posres reference coordinates */
999 for (int j = 0; j < npbcdim; j++)
1001 sum[j] += atom[ai].m * x[a + ai][j];
1003 totmass += atom[ai].m;
1006 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
1007 for (const auto& restraint : prfb->interactionTypes)
1009 int ai = restraint.ai();
1013 "Position restraint atom index (%d) in moltype '%s' is larger than "
1014 "number of atoms in %s (%d).\n",
1016 *molinfo[molb.type].name,
1020 if (rc_scaling == erscCOM && !hadAtom[ai])
1022 /* Determine the center of mass of the posres reference coordinates */
1023 for (int j = 0; j < npbcdim; j++)
1025 sum[j] += atom[ai].m * x[a + ai][j];
1027 totmass += atom[ai].m;
1032 molb.posres_xA.resize(nat_molb);
1033 for (int i = 0; i < nat_molb; i++)
1035 copy_rvec(x[a + i], molb.posres_xA[i]);
1040 molb.posres_xB.resize(nat_molb);
1041 for (int i = 0; i < nat_molb; i++)
1043 copy_rvec(x[a + i], molb.posres_xB[i]);
1049 if (rc_scaling == erscCOM)
1053 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
1055 for (int j = 0; j < npbcdim; j++)
1057 com[j] = sum[j] / totmass;
1059 GMX_LOG(logger.info)
1061 .appendTextFormatted(
1062 "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f",
1068 if (rc_scaling != erscNO)
1070 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
1072 for (gmx_molblock_t& molb : mtop->molblock)
1074 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
1075 if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1077 std::vector<gmx::RVec>& xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1078 for (int i = 0; i < nat_molb; i++)
1080 for (int j = 0; j < npbcdim; j++)
1082 if (rc_scaling == erscALL)
1084 /* Convert from Cartesian to crystal coordinates */
1085 xp[i][j] *= invbox[j][j];
1086 for (int k = j + 1; k < npbcdim; k++)
1088 xp[i][j] += invbox[k][j] * xp[i][k];
1091 else if (rc_scaling == erscCOM)
1093 /* Subtract the center of mass */
1101 if (rc_scaling == erscCOM)
1103 /* Convert the COM from Cartesian to crystal coordinates */
1104 for (int j = 0; j < npbcdim; j++)
1106 com[j] *= invbox[j][j];
1107 for (int k = j + 1; k < npbcdim; k++)
1109 com[j] += invbox[k][j] * com[k];
1120 static void gen_posres(gmx_mtop_t* mtop,
1121 gmx::ArrayRef<const MoleculeInformation> mi,
1129 const gmx::MDLogger& logger)
1131 read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi, logger);
1132 /* It is safer to simply read the b-state posres rather than trying
1133 * to be smart and copy the positions.
1135 read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi, logger);
1138 static void set_wall_atomtype(PreprocessingAtomTypes* at,
1142 const gmx::MDLogger& logger)
1148 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Searching the wall atom type(s)");
1150 for (i = 0; i < ir->nwall; i++)
1152 ir->wall_atomtype[i] = at->atomTypeFromName(opts->wall_atomtype[i]);
1153 if (ir->wall_atomtype[i] == NOTSET)
1155 std::string warningMessage = gmx::formatString(
1156 "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1157 warning_error(wi, warningMessage.c_str());
1162 static int nrdf_internal(const t_atoms* atoms)
1167 for (i = 0; i < atoms->nr; i++)
1169 /* Vsite ptype might not be set here yet, so also check the mass */
1170 if ((atoms->atom[i].ptype == eptAtom || atoms->atom[i].ptype == eptNucleus) && atoms->atom[i].m > 0)
1177 case 0: // Fall through intended
1178 case 1: nrdf = 0; break;
1179 case 2: nrdf = 1; break;
1180 default: nrdf = nmass * 3 - 6; break;
1186 static void spline1d(double dx, const double* y, int n, double* u, double* y2)
1194 for (i = 1; i < n - 1; i++)
1196 p = 0.5 * y2[i - 1] + 2.0;
1198 q = (y[i + 1] - 2.0 * y[i] + y[i - 1]) / dx;
1199 u[i] = (3.0 * q / dx - 0.5 * u[i - 1]) / p;
1204 for (i = n - 2; i >= 0; i--)
1206 y2[i] = y2[i] * y2[i + 1] + u[i];
1212 interpolate1d(double xmin, double dx, const double* ya, const double* y2a, double x, double* y, double* y1)
1217 ix = static_cast<int>((x - xmin) / dx);
1219 a = (xmin + (ix + 1) * dx - x) / dx;
1220 b = (x - xmin - ix * dx) / dx;
1222 *y = a * ya[ix] + b * ya[ix + 1]
1223 + ((a * a * a - a) * y2a[ix] + (b * b * b - b) * y2a[ix + 1]) * (dx * dx) / 6.0;
1224 *y1 = (ya[ix + 1] - ya[ix]) / dx - (3.0 * a * a - 1.0) / 6.0 * dx * y2a[ix]
1225 + (3.0 * b * b - 1.0) / 6.0 * dx * y2a[ix + 1];
1229 static void setup_cmap(int grid_spacing, int nc, gmx::ArrayRef<const real> grid, gmx_cmap_t* cmap_grid)
1231 int i, j, k, ii, jj, kk, idx;
1233 double dx, xmin, v, v1, v2, v12;
1236 std::vector<double> tmp_u(2 * grid_spacing, 0.0);
1237 std::vector<double> tmp_u2(2 * grid_spacing, 0.0);
1238 std::vector<double> tmp_yy(2 * grid_spacing, 0.0);
1239 std::vector<double> tmp_y1(2 * grid_spacing, 0.0);
1240 std::vector<double> tmp_t2(2 * grid_spacing * 2 * grid_spacing, 0.0);
1241 std::vector<double> tmp_grid(2 * grid_spacing * 2 * grid_spacing, 0.0);
1243 dx = 360.0 / grid_spacing;
1244 xmin = -180.0 - dx * grid_spacing / 2;
1246 for (kk = 0; kk < nc; kk++)
1248 /* Compute an offset depending on which cmap we are using
1249 * Offset will be the map number multiplied with the
1250 * grid_spacing * grid_spacing * 2
1252 offset = kk * grid_spacing * grid_spacing * 2;
1254 for (i = 0; i < 2 * grid_spacing; i++)
1256 ii = (i + grid_spacing - grid_spacing / 2) % grid_spacing;
1258 for (j = 0; j < 2 * grid_spacing; j++)
1260 jj = (j + grid_spacing - grid_spacing / 2) % grid_spacing;
1261 tmp_grid[i * grid_spacing * 2 + j] = grid[offset + ii * grid_spacing + jj];
1265 for (i = 0; i < 2 * grid_spacing; i++)
1268 &(tmp_grid[2 * grid_spacing * i]),
1271 &(tmp_t2[2 * grid_spacing * i]));
1274 for (i = grid_spacing / 2; i < grid_spacing + grid_spacing / 2; i++)
1276 ii = i - grid_spacing / 2;
1277 phi = ii * dx - 180.0;
1279 for (j = grid_spacing / 2; j < grid_spacing + grid_spacing / 2; j++)
1281 jj = j - grid_spacing / 2;
1282 psi = jj * dx - 180.0;
1284 for (k = 0; k < 2 * grid_spacing; k++)
1288 &(tmp_grid[2 * grid_spacing * k]),
1289 &(tmp_t2[2 * grid_spacing * k]),
1295 spline1d(dx, tmp_yy.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1296 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1297 spline1d(dx, tmp_y1.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1298 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1300 idx = ii * grid_spacing + jj;
1301 cmap_grid->cmapdata[kk].cmap[idx * 4] = grid[offset + ii * grid_spacing + jj];
1302 cmap_grid->cmapdata[kk].cmap[idx * 4 + 1] = v1;
1303 cmap_grid->cmapdata[kk].cmap[idx * 4 + 2] = v2;
1304 cmap_grid->cmapdata[kk].cmap[idx * 4 + 3] = v12;
1310 static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
1314 cmap_grid->grid_spacing = grid_spacing;
1315 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
1317 cmap_grid->cmapdata.resize(ngrid);
1319 for (i = 0; i < ngrid; i++)
1321 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
1326 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
1328 int count, count_mol;
1331 for (const gmx_molblock_t& molb : mtop->molblock)
1334 gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1336 for (int i = 0; i < F_NRE; i++)
1340 count_mol += 3 * interactions[i].size();
1342 else if (interaction_function[i].flags & IF_CONSTRAINT)
1344 count_mol += interactions[i].size();
1348 if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1350 std::string warningMessage = gmx::formatString(
1351 "Molecule type '%s' has %d constraints.\n"
1352 "For stability and efficiency there should not be more constraints than "
1353 "internal number of degrees of freedom: %d.\n",
1354 *mi[molb.type].name,
1356 nrdf_internal(&mi[molb.type].atoms));
1357 warning(wi, warningMessage.c_str());
1359 count += molb.nmol * count_mol;
1365 static real calc_temp(const gmx_mtop_t* mtop, const t_inputrec* ir, rvec* v)
1368 for (const AtomProxy atomP : AtomRange(*mtop))
1370 const t_atom& local = atomP.atom();
1371 int i = atomP.globalAtomNumber();
1372 sum_mv2 += local.m * norm2(v[i]);
1376 for (int g = 0; g < ir->opts.ngtc; g++)
1378 nrdf += ir->opts.nrdf[g];
1381 return sum_mv2 / (nrdf * BOLTZ);
1384 static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
1392 for (i = 0; i < ir->opts.ngtc; i++)
1394 if (ir->opts.tau_t[i] < 0)
1400 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1406 std::string warningMessage = gmx::formatString(
1407 "Some temperature coupling groups do not use temperature coupling. We will assume "
1408 "their temperature is not more than %.3f K. If their temperature is higher, the "
1409 "energy error and the Verlet buffer might be underestimated.",
1411 warning(wi, warningMessage.c_str());
1417 /* Checks if there are unbound atoms in moleculetype molt.
1418 * Prints a note for each unbound atoms and a warning if any is present.
1420 static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1422 const t_atoms* atoms = &molt->atoms;
1426 /* Only one atom, there can't be unbound atoms */
1430 std::vector<int> count(atoms->nr, 0);
1432 for (int ftype = 0; ftype < F_NRE; ftype++)
1434 if (((interaction_function[ftype].flags & IF_BOND) && NRAL(ftype) == 2 && ftype != F_CONNBONDS)
1435 || (interaction_function[ftype].flags & IF_CONSTRAINT) || ftype == F_SETTLE)
1437 const InteractionList& il = molt->ilist[ftype];
1438 const int nral = NRAL(ftype);
1440 for (int i = 0; i < il.size(); i += 1 + nral)
1442 for (int j = 0; j < nral; j++)
1444 count[il.iatoms[i + 1 + j]]++;
1450 int numDanglingAtoms = 0;
1451 for (int a = 0; a < atoms->nr; a++)
1453 if (atoms->atom[a].ptype != eptVSite && count[a] == 0)
1457 GMX_LOG(logger.warning)
1459 .appendTextFormatted(
1460 "Atom %d '%s' in moleculetype '%s' is not bound by a potential or "
1461 "constraint to any other atom in the same moleculetype.",
1463 *atoms->atomname[a],
1470 if (numDanglingAtoms > 0)
1472 std::string warningMessage = gmx::formatString(
1473 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
1474 "other atom in the same moleculetype. Although technically this might not cause "
1475 "issues in a simulation, this often means that the user forgot to add a "
1476 "bond/potential/constraint or put multiple molecules in the same moleculetype "
1477 "definition by mistake. Run with -v to get information for each atom.",
1480 warning_note(wi, warningMessage.c_str());
1484 /* Checks all moleculetypes for unbound atoms */
1485 static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1487 for (const gmx_moltype_t& molt : mtop->moltype)
1489 checkForUnboundAtoms(&molt, bVerbose, wi, logger);
1493 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1495 * The specific decoupled modes this routine check for are angle modes
1496 * where the two bonds are constrained and the atoms a both ends are only
1497 * involved in a single constraint; the mass of the two atoms needs to
1498 * differ by more than \p massFactorThreshold.
1500 static bool haveDecoupledModeInMol(const gmx_moltype_t& molt,
1501 gmx::ArrayRef<const t_iparams> iparams,
1502 real massFactorThreshold)
1504 if (molt.ilist[F_CONSTR].empty() && molt.ilist[F_CONSTRNC].empty())
1509 const t_atom* atom = molt.atoms.atom;
1511 const auto atomToConstraints =
1512 gmx::make_at2con(molt, iparams, gmx::FlexibleConstraintTreatment::Exclude);
1514 bool haveDecoupledMode = false;
1515 for (int ftype = 0; ftype < F_NRE; ftype++)
1517 if (interaction_function[ftype].flags & IF_ATYPE)
1519 const int nral = NRAL(ftype);
1520 const InteractionList& il = molt.ilist[ftype];
1521 for (int i = 0; i < il.size(); i += 1 + nral)
1523 /* Here we check for the mass difference between the atoms
1524 * at both ends of the angle, that the atoms at the ends have
1525 * 1 contraint and the atom in the middle at least 3; we check
1526 * that the 3 atoms are linked by constraints below.
1527 * We check for at least three constraints for the middle atom,
1528 * since with only the two bonds in the angle, we have 3-atom
1529 * molecule, which has much more thermal exhange in this single
1530 * angle mode than molecules with more atoms.
1531 * Note that this check also catches molecules where more atoms
1532 * are connected to one or more atoms in the angle, but by
1533 * bond potentials instead of angles. But such cases will not
1534 * occur in "normal" molecules and it doesn't hurt running
1535 * those with higher accuracy settings as well.
1537 int a0 = il.iatoms[1 + i];
1538 int a1 = il.iatoms[1 + i + 1];
1539 int a2 = il.iatoms[1 + i + 2];
1540 if ((atom[a0].m > atom[a2].m * massFactorThreshold || atom[a2].m > atom[a0].m * massFactorThreshold)
1541 && atomToConstraints[a0].ssize() == 1 && atomToConstraints[a2].ssize() == 1
1542 && atomToConstraints[a1].ssize() >= 3)
1544 int constraint0 = atomToConstraints[a0][0];
1545 int constraint2 = atomToConstraints[a2][0];
1547 bool foundAtom0 = false;
1548 bool foundAtom2 = false;
1549 for (const int constraint : atomToConstraints[a1])
1551 if (constraint == constraint0)
1555 if (constraint == constraint2)
1560 if (foundAtom0 && foundAtom2)
1562 haveDecoupledMode = true;
1569 return haveDecoupledMode;
1572 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1574 * When decoupled modes are present and the accuray in insufficient,
1575 * this routine issues a warning if the accuracy is insufficient.
1577 static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec* ir, warninp* wi)
1579 /* We only have issues with decoupled modes with normal MD.
1580 * With stochastic dynamics equipartitioning is enforced strongly.
1587 /* When atoms of very different mass are involved in an angle potential
1588 * and both bonds in the angle are constrained, the dynamic modes in such
1589 * angles have very different periods and significant energy exchange
1590 * takes several nanoseconds. Thus even a small amount of error in
1591 * different algorithms can lead to violation of equipartitioning.
1592 * The parameters below are mainly based on an all-atom chloroform model
1593 * with all bonds constrained. Equipartitioning is off by more than 1%
1594 * (need to run 5-10 ns) when the difference in mass between H and Cl
1595 * is more than a factor 13 and the accuracy is less than the thresholds
1596 * given below. This has been verified on some other molecules.
1598 * Note that the buffer and shake parameters have unit length and
1599 * energy/time, respectively, so they will "only" work correctly
1600 * for atomistic force fields using MD units.
1602 const real massFactorThreshold = 13.0;
1603 const real bufferToleranceThreshold = 1e-4;
1604 const int lincsIterationThreshold = 2;
1605 const int lincsOrderThreshold = 4;
1606 const real shakeToleranceThreshold = 0.005 * ir->delta_t;
1608 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold
1609 && ir->nProjOrder >= lincsOrderThreshold);
1610 bool shakeWithSufficientTolerance =
1611 (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1 * shakeToleranceThreshold);
1612 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol <= 1.1 * bufferToleranceThreshold
1613 && (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1618 bool haveDecoupledMode = false;
1619 for (const gmx_moltype_t& molt : mtop->moltype)
1621 if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams, massFactorThreshold))
1623 haveDecoupledMode = true;
1627 if (haveDecoupledMode)
1629 std::string message = gmx::formatString(
1630 "There are atoms at both ends of an angle, connected by constraints "
1631 "and with masses that differ by more than a factor of %g. This means "
1632 "that there are likely dynamic modes that are only very weakly coupled.",
1633 massFactorThreshold);
1634 if (ir->cutoff_scheme == ecutsVERLET)
1636 message += gmx::formatString(
1637 " To ensure good equipartitioning, you need to either not use "
1638 "constraints on all bonds (but, if possible, only on bonds involving "
1639 "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1640 "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1641 ">= %d or SHAKE tolerance <= %g",
1643 bufferToleranceThreshold,
1644 lincsIterationThreshold,
1645 lincsOrderThreshold,
1646 shakeToleranceThreshold);
1650 message += gmx::formatString(
1651 " To ensure good equipartitioning, we suggest to switch to the %s "
1652 "cutoff-scheme, since that allows for better control over the Verlet "
1653 "buffer size and thus over the energy drift.",
1654 ecutscheme_names[ecutsVERLET]);
1656 warning(wi, message);
1660 static void set_verlet_buffer(const gmx_mtop_t* mtop,
1665 const gmx::MDLogger& logger)
1667 GMX_LOG(logger.info)
1669 .appendTextFormatted(
1670 "Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K",
1674 /* Calculate the buffer size for simple atom vs atoms list */
1675 VerletbufListSetup listSetup1x1;
1676 listSetup1x1.cluster_size_i = 1;
1677 listSetup1x1.cluster_size_j = 1;
1678 const real rlist_1x1 = calcVerletBufferSize(
1679 *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup1x1);
1681 /* Set the pair-list buffer size in ir */
1682 VerletbufListSetup listSetup4x4 = verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1683 ir->rlist = calcVerletBufferSize(
1684 *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup4x4);
1686 const int n_nonlin_vsite = gmx::countNonlinearVsites(*mtop);
1687 if (n_nonlin_vsite > 0)
1689 std::string warningMessage = gmx::formatString(
1690 "There are %d non-linear virtual site constructions. Their contribution to the "
1691 "energy error is approximated. In most cases this does not affect the error "
1694 warning_note(wi, warningMessage);
1697 GMX_LOG(logger.info)
1699 .appendTextFormatted(
1700 "Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm",
1704 rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
1706 GMX_LOG(logger.info)
1708 .appendTextFormatted(
1709 "Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm",
1710 listSetup4x4.cluster_size_i,
1711 listSetup4x4.cluster_size_j,
1713 ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
1715 GMX_LOG(logger.info)
1717 .appendTextFormatted(
1718 "Note that mdrun will redetermine rlist based on the actual pair-list setup");
1720 if (gmx::square(ir->rlist) >= max_cutoff2(ir->pbcType, box))
1723 "The pair-list cut-off (%g nm) is longer than half the shortest box vector or "
1724 "longer than the smallest box diagonal element (%g nm). Increase the box size or "
1725 "decrease nstlist or increase verlet-buffer-tolerance.",
1727 std::sqrt(max_cutoff2(ir->pbcType, box)));
1731 int gmx_grompp(int argc, char* argv[])
1733 const char* desc[] = {
1734 "[THISMODULE] (the gromacs preprocessor)",
1735 "reads a molecular topology file, checks the validity of the",
1736 "file, expands the topology from a molecular description to an atomic",
1737 "description. The topology file contains information about",
1738 "molecule types and the number of molecules, the preprocessor",
1739 "copies each molecule as needed. ",
1740 "There is no limitation on the number of molecule types. ",
1741 "Bonds and bond-angles can be converted into constraints, separately",
1742 "for hydrogens and heavy atoms.",
1743 "Then a coordinate file is read and velocities can be generated",
1744 "from a Maxwellian distribution if requested.",
1745 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1746 "(eg. number of MD steps, time step, cut-off), and others such as",
1747 "NEMD parameters, which are corrected so that the net acceleration",
1749 "Eventually a binary file is produced that can serve as the sole input",
1750 "file for the MD program.[PAR]",
1752 "[THISMODULE] uses the atom names from the topology file. The atom names",
1753 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1754 "warnings when they do not match the atom names in the topology.",
1755 "Note that the atom names are irrelevant for the simulation as",
1756 "only the atom types are used for generating interaction parameters.[PAR]",
1758 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1759 "etc. The preprocessor supports the following keywords::",
1762 " #ifndef VARIABLE",
1765 " #define VARIABLE",
1767 " #include \"filename\"",
1768 " #include <filename>",
1770 "The functioning of these statements in your topology may be modulated by",
1771 "using the following two flags in your [REF].mdp[ref] file::",
1773 " define = -DVARIABLE1 -DVARIABLE2",
1774 " include = -I/home/john/doe",
1776 "For further information a C-programming textbook may help you out.",
1777 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1778 "topology file written out so that you can verify its contents.[PAR]",
1780 "When using position restraints, a file with restraint coordinates",
1781 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1782 "for [TT]-c[tt]). For free energy calculations, separate reference",
1783 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1784 "otherwise they will be equal to those of the A topology.[PAR]",
1786 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1787 "The last frame with coordinates and velocities will be read,",
1788 "unless the [TT]-time[tt] option is used. Only if this information",
1789 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1790 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1791 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1792 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1795 "[THISMODULE] can be used to restart simulations (preserving",
1796 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1797 "However, for simply changing the number of run steps to extend",
1798 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1799 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1800 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1801 "like output frequency, then supplying the checkpoint file to",
1802 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1803 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1804 "the ensemble (if possible) still requires passing the checkpoint",
1805 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1807 "By default, all bonded interactions which have constant energy due to",
1808 "virtual site constructions will be removed. If this constant energy is",
1809 "not zero, this will result in a shift in the total energy. All bonded",
1810 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1811 "all constraints for distances which will be constant anyway because",
1812 "of virtual site constructions will be removed. If any constraints remain",
1813 "which involve virtual sites, a fatal error will result.[PAR]",
1815 "To verify your run input file, please take note of all warnings",
1816 "on the screen, and correct where necessary. Do also look at the contents",
1817 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1818 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1819 "with the [TT]-debug[tt] option which will give you more information",
1820 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1821 "can see the contents of the run input file with the [gmx-dump]",
1822 "program. [gmx-check] can be used to compare the contents of two",
1823 "run input files.[PAR]",
1825 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1826 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1827 "harmless, but usually they are not. The user is advised to carefully",
1828 "interpret the output messages before attempting to bypass them with",
1831 std::vector<MoleculeInformation> mi;
1832 std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1836 const char* mdparin;
1837 bool bNeedVel, bGenVel;
1838 gmx_output_env_t* oenv;
1839 gmx_bool bVerbose = FALSE;
1842 t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
1843 { efMDP, "-po", "mdout", ffWRITE },
1844 { efSTX, "-c", nullptr, ffREAD },
1845 { efSTX, "-r", "restraint", ffOPTRD },
1846 { efSTX, "-rb", "restraint", ffOPTRD },
1847 { efNDX, nullptr, nullptr, ffOPTRD },
1848 { efTOP, nullptr, nullptr, ffREAD },
1849 { efTOP, "-pp", "processed", ffOPTWR },
1850 { efTPR, "-o", nullptr, ffWRITE },
1851 { efTRN, "-t", nullptr, ffOPTRD },
1852 { efEDR, "-e", nullptr, ffOPTRD },
1853 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1854 { efGRO, "-imd", "imdgroup", ffOPTWR },
1855 { efTRN, "-ref", "rotref", ffOPTRW | ffALLOW_MISSING } };
1856 #define NFILE asize(fnm)
1858 /* Command line options */
1859 gmx_bool bRenum = TRUE;
1860 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1864 { "-v", FALSE, etBOOL, { &bVerbose }, "Be loud and noisy" },
1865 { "-time", FALSE, etREAL, { &fr_time }, "Take frame at or first after this time." },
1870 "Remove constant bonded interactions with virtual sites" },
1875 "Number of allowed warnings during input processing. Not for normal use and may "
1876 "generate unstable systems" },
1881 "Set parameters for bonded interactions without defaults to zero instead of "
1882 "generating an error" },
1887 "Renumber atomtypes and minimize number of atomtypes" }
1890 /* Parse the command line */
1891 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
1896 /* Initiate some variables */
1897 gmx::MDModules mdModules;
1898 t_inputrec irInstance;
1899 t_inputrec* ir = &irInstance;
1900 t_gromppopts optsInstance;
1901 t_gromppopts* opts = &optsInstance;
1902 snew(opts->include, STRLEN);
1903 snew(opts->define, STRLEN);
1905 gmx::LoggerBuilder builder;
1906 builder.addTargetStream(gmx::MDLogger::LogLevel::Info, &gmx::TextOutputFile::standardOutput());
1907 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
1908 gmx::LoggerOwner logOwner(builder.build());
1909 const gmx::MDLogger logger(logOwner.logger());
1912 wi = init_warning(TRUE, maxwarn);
1914 /* PARAMETER file processing */
1915 mdparin = opt2fn("-f", NFILE, fnm);
1916 set_warning_line(wi, mdparin, -1);
1919 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1921 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1923 // Now that the MdModules have their options assigned from get_ir, subscribe
1924 // to eventual notifications during pre-processing their data
1925 mdModules.subscribeToPreProcessingNotifications();
1930 GMX_LOG(logger.info)
1932 .appendTextFormatted("checking input for internal consistency...");
1934 check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
1936 if (ir->ld_seed == -1)
1938 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1939 GMX_LOG(logger.info)
1941 .appendTextFormatted("Setting the LD random seed to %" PRId64 "", ir->ld_seed);
1944 if (ir->expandedvals->lmc_seed == -1)
1946 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1947 GMX_LOG(logger.info)
1949 .appendTextFormatted("Setting the lambda MC random seed to %d", ir->expandedvals->lmc_seed);
1952 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1953 bGenVel = (bNeedVel && opts->bGenVel);
1954 if (bGenVel && ir->bContinuation)
1956 std::string warningMessage = gmx::formatString(
1957 "Generating velocities is inconsistent with attempting "
1958 "to continue a previous run. Choose only one of "
1959 "gen-vel = yes and continuation = yes.");
1960 warning_error(wi, warningMessage);
1963 std::array<InteractionsOfType, F_NRE> interactions;
1965 PreprocessingAtomTypes atypes;
1968 pr_symtab(debug, 0, "Just opened", &sys.symtab);
1971 const char* fn = ftp2fn(efTOP, NFILE, fnm);
1972 if (!gmx_fexist(fn))
1974 gmx_fatal(FARGS, "%s does not exist", fn);
1979 opt2fn_null("-pp", NFILE, fnm),
1980 opt2fn("-c", NFILE, fnm),
1990 &intermolecular_interactions,
2001 pr_symtab(debug, 0, "After new_status", &sys.symtab);
2005 /* set parameters for virtual site construction (not for vsiten) */
2006 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2008 nvsite += set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions, logger);
2010 /* now throw away all obsolete bonds, angles and dihedrals: */
2011 /* note: constraints are ALWAYS removed */
2014 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2016 clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds, logger);
2020 if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == econtSHAKE))
2022 if (ir->eI == eiCG || ir->eI == eiLBFGS)
2024 std::string warningMessage = gmx::formatString("Can not do %s with %s, use %s",
2026 econstr_names[econtSHAKE],
2027 econstr_names[econtLINCS]);
2028 warning_error(wi, warningMessage);
2030 if (ir->bPeriodicMols)
2032 std::string warningMessage =
2033 gmx::formatString("Can not do periodic molecules with %s, use %s",
2034 econstr_names[econtSHAKE],
2035 econstr_names[econtLINCS]);
2036 warning_error(wi, warningMessage);
2040 if (EI_SD(ir->eI) && ir->etc != etcNO)
2042 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
2045 /* Check for errors in the input now, since they might cause problems
2046 * during processing further down.
2048 check_warning_error(wi, FARGS);
2050 if (nint_ftype(&sys, mi, F_POSRES) > 0 || nint_ftype(&sys, mi, F_FBPOSRES) > 0)
2052 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
2054 std::string warningMessage = gmx::formatString(
2055 "You are combining position restraints with %s pressure coupling, which can "
2056 "lead to instabilities. If you really want to combine position restraints with "
2057 "pressure coupling, we suggest to use %s pressure coupling instead.",
2058 EPCOUPLTYPE(ir->epc),
2059 EPCOUPLTYPE(epcBERENDSEN));
2060 warning_note(wi, warningMessage);
2063 const char* fn = opt2fn("-r", NFILE, fnm);
2066 if (!gmx_fexist(fn))
2069 "Cannot find position restraint file %s (option -r).\n"
2070 "From GROMACS-2018, you need to specify the position restraint "
2071 "coordinate files explicitly to avoid mistakes, although you can "
2072 "still use the same file as you specify for the -c option.",
2076 if (opt2bSet("-rb", NFILE, fnm))
2078 fnB = opt2fn("-rb", NFILE, fnm);
2079 if (!gmx_fexist(fnB))
2082 "Cannot find B-state position restraint file %s (option -rb).\n"
2083 "From GROMACS-2018, you need to specify the position restraint "
2084 "coordinate files explicitly to avoid mistakes, although you can "
2085 "still use the same file as you specify for the -c option.",
2096 std::string message = gmx::formatString("Reading position restraint coords from %s", fn);
2097 if (strcmp(fn, fnB) != 0)
2099 message += gmx::formatString(" and %s", fnB);
2101 GMX_LOG(logger.info).asParagraph().appendText(message);
2103 gen_posres(&sys, mi, fn, fnB, ir->refcoord_scaling, ir->pbcType, ir->posres_com, ir->posres_comB, wi, logger);
2106 /* If we are using CMAP, setup the pre-interpolation grid */
2107 if (interactions[F_CMAP].ncmap() > 0)
2109 init_cmap_grid(&sys.ffparams.cmap_grid,
2110 interactions[F_CMAP].cmapAngles,
2111 interactions[F_CMAP].cmakeGridSpacing);
2112 setup_cmap(interactions[F_CMAP].cmakeGridSpacing,
2113 interactions[F_CMAP].cmapAngles,
2114 interactions[F_CMAP].cmap,
2115 &sys.ffparams.cmap_grid);
2118 set_wall_atomtype(&atypes, opts, ir, wi, logger);
2121 atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
2124 if (ir->implicit_solvent)
2126 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
2129 /* PELA: Copy the atomtype data to the topology atomtype list */
2130 atypes.copyTot_atomtypes(&(sys.atomtypes));
2134 pr_symtab(debug, 0, "After atype.renumberTypes", &sys.symtab);
2139 GMX_LOG(logger.info).asParagraph().appendTextFormatted("converting bonded parameters...");
2142 const int ntype = atypes.size();
2143 convertInteractionsOfType(
2144 ntype, interactions, mi, intermolecular_interactions.get(), comb, reppow, fudgeQQ, &sys);
2148 pr_symtab(debug, 0, "After converInteractionsOfType", &sys.symtab);
2151 /* set ptype to VSite for virtual sites */
2152 for (gmx_moltype_t& moltype : sys.moltype)
2154 set_vsites_ptype(FALSE, &moltype, logger);
2158 pr_symtab(debug, 0, "After virtual sites", &sys.symtab);
2160 /* Check velocity for virtual sites and shells */
2163 check_vel(&sys, state.v.rvec_array());
2166 /* check for shells and inpurecs */
2167 check_shells_inputrec(&sys, ir, wi);
2170 check_mol(&sys, wi);
2172 checkForUnboundAtoms(&sys, bVerbose, wi, logger);
2174 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2176 check_bonds_timestep(&sys, ir->delta_t, wi);
2179 checkDecoupledModeAccuracy(&sys, ir, wi);
2181 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2185 "Zero-step energy minimization will alter the coordinates before calculating the "
2186 "energy. If you just want the energy of a single point, try zero-step MD (with "
2187 "unconstrained_start = yes). To do multiple single-point energy evaluations of "
2188 "different configurations of the same topology, use mdrun -rerun.");
2191 check_warning_error(wi, FARGS);
2195 GMX_LOG(logger.info).asParagraph().appendTextFormatted("initialising group options...");
2197 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifier(), ir, wi);
2199 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2201 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2205 if (EI_MD(ir->eI) && ir->etc == etcNO)
2209 buffer_temp = opts->tempi;
2213 buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2215 if (buffer_temp > 0)
2217 std::string warningMessage = gmx::formatString(
2218 "NVE simulation: will use the initial temperature of %.3f K for "
2219 "determining the Verlet buffer size",
2221 warning_note(wi, warningMessage);
2225 std::string warningMessage = gmx::formatString(
2226 "NVE simulation with an initial temperature of zero: will use a Verlet "
2227 "buffer of %d%%. Check your energy drift!",
2228 gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
2229 warning_note(wi, warningMessage);
2234 buffer_temp = get_max_reference_temp(ir, wi);
2237 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2239 /* NVE with initial T=0: we add a fixed ratio to rlist.
2240 * Since we don't actually use verletbuf_tol, we set it to -1
2241 * so it can't be misused later.
2243 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2244 ir->verletbuf_tol = -1;
2248 /* We warn for NVE simulations with a drift tolerance that
2249 * might result in a 1(.1)% drift over the total run-time.
2250 * Note that we can't warn when nsteps=0, since we don't
2251 * know how many steps the user intends to run.
2253 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 && ir->nsteps > 0)
2255 const real driftTolerance = 0.01;
2256 /* We use 2 DOF per atom = 2kT pot+kin energy,
2257 * to be on the safe side with constraints.
2259 const real totalEnergyDriftPerAtomPerPicosecond =
2260 2 * BOLTZ * buffer_temp / (ir->nsteps * ir->delta_t);
2262 if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
2264 std::string warningMessage = gmx::formatString(
2265 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
2266 "NVE simulation of length %g ps, which can give a final drift of "
2267 "%d%%. For conserving energy to %d%% when using constraints, you "
2268 "might need to set verlet-buffer-tolerance to %.1e.",
2270 ir->nsteps * ir->delta_t,
2271 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
2272 gmx::roundToInt(100 * driftTolerance),
2273 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
2274 warning_note(wi, warningMessage);
2278 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi, logger);
2283 /* Init the temperature coupling state */
2284 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2288 pr_symtab(debug, 0, "After index", &sys.symtab);
2291 triple_check(mdparin, ir, &sys, wi);
2292 close_symtab(&sys.symtab);
2295 pr_symtab(debug, 0, "After close", &sys.symtab);
2298 if (ir->eI == eiMimic)
2300 generate_qmexcl(&sys, ir, logger);
2303 if (ftp2bSet(efTRN, NFILE, fnm))
2307 GMX_LOG(logger.info)
2309 .appendTextFormatted("getting data from old trajectory ...");
2311 cont_status(ftp2fn(efTRN, NFILE, fnm),
2312 ftp2fn_null(efEDR, NFILE, fnm),
2323 if (ir->pbcType == PbcType::XY && ir->nwall != 2)
2325 clear_rvec(state.box[ZZ]);
2328 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2330 /* Calculate the optimal grid dimensions */
2332 EwaldBoxZScaler boxScaler(*ir);
2333 boxScaler.scaleBox(state.box, scaledBox);
2335 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2337 /* Mark fourier_spacing as not used */
2338 ir->fourier_spacing = 0;
2340 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2342 set_warning_line(wi, mdparin, -1);
2344 wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2346 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2347 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize, &(ir->nkx), &(ir->nky), &(ir->nkz));
2348 if (ir->nkx < minGridSize || ir->nky < minGridSize || ir->nkz < minGridSize)
2351 "The PME grid size should be >= 2*(pme-order - 1); either manually "
2352 "increase the grid size or decrease pme-order");
2356 /* MRS: eventually figure out better logic for initializing the fep
2357 values that makes declaring the lambda and declaring the state not
2358 potentially conflict if not handled correctly. */
2359 if (ir->efep != efepNO)
2361 state.fep_state = ir->fepvals->init_fep_state;
2362 for (i = 0; i < efptNR; i++)
2364 /* init_lambda trumps state definitions*/
2365 if (ir->fepvals->init_lambda >= 0)
2367 state.lambda[i] = ir->fepvals->init_lambda;
2371 if (ir->fepvals->all_lambda[i] == nullptr)
2373 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2377 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2383 pull_t* pull = nullptr;
2387 pull = set_pull_init(ir, &sys, state.x.rvec_array(), state.box, state.lambda[efptMASS], wi);
2390 /* Modules that supply external potential for pull coordinates
2391 * should register those potentials here. finish_pull() will check
2392 * that providers have been registerd for all external potentials.
2397 tensor compressibility = { { 0 } };
2398 if (ir->epc != epcNO)
2400 copy_mat(ir->compress, compressibility);
2402 setStateDependentAwhParams(
2410 ir->efep != efepNO ? ir->fepvals->all_lambda[efptFEP][ir->fepvals->init_fep_state] : 0,
2422 set_reference_positions(ir->rot,
2423 state.x.rvec_array(),
2425 opt2fn("-ref", NFILE, fnm),
2426 opt2bSet("-ref", NFILE, fnm),
2430 /* reset_multinr(sys); */
2432 if (EEL_PME(ir->coulombtype))
2434 float ratio = pme_load_estimate(sys, *ir, state.box);
2435 GMX_LOG(logger.info)
2437 .appendTextFormatted(
2438 "Estimate for the relative computational load of the PME mesh part: %.2f", ratio);
2439 /* With free energy we might need to do PME both for the A and B state
2440 * charges. This will double the cost, but the optimal performance will
2441 * then probably be at a slightly larger cut-off and grid spacing.
2443 if ((ir->efep == efepNO && ratio > 1.0 / 2.0) || (ir->efep != efepNO && ratio > 2.0 / 3.0))
2447 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2448 "and for highly parallel simulations between 0.25 and 0.33,\n"
2449 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2450 if (ir->efep != efepNO)
2453 "For free energy simulations, the optimal load limit increases from "
2460 double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2461 std::string warningMessage =
2462 gmx::formatString("This run will generate roughly %.0f Mb of data", cio);
2465 set_warning_line(wi, mdparin, -1);
2466 warning_note(wi, warningMessage);
2470 GMX_LOG(logger.info).asParagraph().appendText(warningMessage);
2474 // Add the md modules internal parameters that are not mdp options
2475 // e.g., atom indices
2478 gmx::KeyValueTreeBuilder internalParameterBuilder;
2479 mdModules.notifier().preProcessingNotifications_.notify(internalParameterBuilder.rootObject());
2480 ir->internalParameters =
2481 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2484 if (ir->comm_mode != ecmNO)
2486 const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
2487 if (ir->nstcomm % nstglobalcomm != 0)
2492 "COM removal frequency is set to (%d).\n"
2493 "Other settings require a global communication frequency of %d.\n"
2494 "Note that this will require additional global communication steps,\n"
2495 "which will reduce performance when using multiple ranks.\n"
2496 "Consider setting nstcomm to a multiple of %d.",
2505 GMX_LOG(logger.info).asParagraph().appendTextFormatted("writing run input file...");
2508 done_warning(wi, FARGS);
2509 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, &sys);
2511 /* Output IMD group, if bIMD is TRUE */
2512 gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, &sys, NFILE, fnm);
2514 sfree(opts->define);
2515 sfree(opts->wall_atomtype[0]);
2516 sfree(opts->wall_atomtype[1]);
2517 sfree(opts->include);
2518 for (auto& mol : mi)
2520 // Some of the contents of molinfo have been stolen, so
2521 // fullCleanUp can't be called.
2522 mol.partialCleanUp();
2524 done_inputrec_strings();
2525 output_env_done(oenv);