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.
52 #include <sys/types.h>
54 #include "gromacs/applied_forces/awh/read_params.h"
55 #include "gromacs/commandline/pargs.h"
56 #include "gromacs/ewald/ewald_utils.h"
57 #include "gromacs/ewald/pme.h"
58 #include "gromacs/fft/calcgrid.h"
59 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/enxio.h"
61 #include "gromacs/fileio/tpxio.h"
62 #include "gromacs/fileio/trxio.h"
63 #include "gromacs/fileio/warninp.h"
64 #include "gromacs/gmxpreprocess/add_par.h"
65 #include "gromacs/gmxpreprocess/convparm.h"
66 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
67 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
68 #include "gromacs/gmxpreprocess/grompp_impl.h"
69 #include "gromacs/gmxpreprocess/notset.h"
70 #include "gromacs/gmxpreprocess/readir.h"
71 #include "gromacs/gmxpreprocess/tomorse.h"
72 #include "gromacs/gmxpreprocess/topio.h"
73 #include "gromacs/gmxpreprocess/toputil.h"
74 #include "gromacs/gmxpreprocess/vsite_parm.h"
75 #include "gromacs/imd/imd.h"
76 #include "gromacs/math/functions.h"
77 #include "gromacs/math/invertmatrix.h"
78 #include "gromacs/math/units.h"
79 #include "gromacs/math/vec.h"
80 #include "gromacs/mdlib/calc_verletbuf.h"
81 #include "gromacs/mdlib/compute_io.h"
82 #include "gromacs/mdlib/constr.h"
83 #include "gromacs/mdlib/md_support.h"
84 #include "gromacs/mdlib/perf_est.h"
85 #include "gromacs/mdlib/vsite.h"
86 #include "gromacs/mdrun/mdmodules.h"
87 #include "gromacs/mdtypes/inputrec.h"
88 #include "gromacs/mdtypes/md_enums.h"
89 #include "gromacs/mdtypes/nblist.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/pbcutil/boxutilities.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/pulling/pull.h"
94 #include "gromacs/random/seed.h"
95 #include "gromacs/topology/ifunc.h"
96 #include "gromacs/topology/mtop_util.h"
97 #include "gromacs/topology/symtab.h"
98 #include "gromacs/topology/topology.h"
99 #include "gromacs/trajectory/trajectoryframe.h"
100 #include "gromacs/utility/arraysize.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/filestream.h"
105 #include "gromacs/utility/futil.h"
106 #include "gromacs/utility/gmxassert.h"
107 #include "gromacs/utility/keyvaluetreebuilder.h"
108 #include "gromacs/utility/listoflists.h"
109 #include "gromacs/utility/logger.h"
110 #include "gromacs/utility/loggerbuilder.h"
111 #include "gromacs/utility/mdmodulesnotifiers.h"
112 #include "gromacs/utility/smalloc.h"
113 #include "gromacs/utility/snprintf.h"
115 /* TODO The implementation details should move to their own source file. */
116 InteractionOfType::InteractionOfType(gmx::ArrayRef<const int> atoms,
117 gmx::ArrayRef<const real> params,
118 const std::string& name) :
119 atoms_(atoms.begin(), atoms.end()), interactionTypeName_(name)
122 params.size() <= forceParam_.size(),
123 gmx::formatString("Cannot have more parameters than the maximum number possible (%d)", MAXFORCEPARAM)
125 std::array<real, MAXFORCEPARAM>::iterator 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 == ParticleType::Shell || local.ptype == ParticleType::Bond
472 || local.ptype == ParticleType::VSite)
479 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
483 for (const AtomProxy atomP : AtomRange(*mtop))
485 const t_atom& local = atomP.atom();
486 if (local.ptype == ParticleType::Shell || local.ptype == ParticleType::Bond)
491 if ((nshells > 0) && (ir->nstcalcenergy != 1))
493 set_warning_line(wi, "unknown", -1);
494 std::string warningMessage = gmx::formatString(
495 "There are %d shells, changing nstcalcenergy from %d to 1", nshells, ir->nstcalcenergy);
496 ir->nstcalcenergy = 1;
497 warning(wi, warningMessage.c_str());
501 /* TODO Decide whether this function can be consolidated with
502 * gmx_mtop_ftype_count */
503 static int nint_ftype(gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, int ftype)
506 for (const gmx_molblock_t& molb : mtop->molblock)
508 nint += molb.nmol * mi[molb.type].interactions[ftype].size();
514 /* This routine reorders the molecule type array
515 * in the order of use in the molblocks,
516 * unused molecule types are deleted.
518 static void renumber_moltypes(gmx_mtop_t* sys, std::vector<MoleculeInformation>* molinfo)
521 std::vector<int> order;
522 for (gmx_molblock_t& molblock : sys->molblock)
524 const auto found = std::find_if(order.begin(), order.end(), [&molblock](const auto& entry) {
525 return molblock.type == entry;
527 if (found == order.end())
529 /* This type did not occur yet, add it */
530 order.push_back(molblock.type);
531 molblock.type = order.size() - 1;
535 molblock.type = std::distance(order.begin(), found);
539 /* We still need to reorder the molinfo structs */
540 std::vector<MoleculeInformation> minew(order.size());
542 for (auto& mi : *molinfo)
544 const auto found = std::find(order.begin(), order.end(), index);
545 if (found != order.end())
547 int pos = std::distance(order.begin(), found);
552 // Need to manually clean up memory ....
561 static void molinfo2mtop(gmx::ArrayRef<const MoleculeInformation> mi, gmx_mtop_t* mtop)
563 mtop->moltype.resize(mi.size());
565 for (const auto& mol : mi)
567 gmx_moltype_t& molt = mtop->moltype[pos];
568 molt.name = mol.name;
569 molt.atoms = mol.atoms;
570 /* ilists are copied later */
571 molt.excls = mol.excls;
576 static void new_status(const char* topfile,
577 const char* topppfile,
585 PreprocessingAtomTypes* atypes,
587 std::vector<MoleculeInformation>* mi,
588 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
589 gmx::ArrayRef<InteractionsOfType> interactions,
590 CombinationRule* comb,
595 const gmx::MDLogger& logger)
597 std::vector<gmx_molblock_t> molblock;
599 bool ffParametrizedWithHBondConstraints;
601 /* TOPOLOGY processing */
602 sys->name = do_top(bVerbose,
614 intermolecular_interactions,
617 &ffParametrizedWithHBondConstraints,
621 sys->molblock.clear();
624 for (const gmx_molblock_t& molb : molblock)
626 if (!sys->molblock.empty() && molb.type == sys->molblock.back().type)
628 /* Merge consecutive blocks with the same molecule type */
629 sys->molblock.back().nmol += molb.nmol;
631 else if (molb.nmol > 0)
633 /* Add a new molblock to the topology */
634 sys->molblock.push_back(molb);
636 sys->natoms += molb.nmol * (*mi)[sys->molblock.back().type].atoms.nr;
638 if (sys->molblock.empty())
640 gmx_fatal(FARGS, "No molecules were defined in the system");
643 renumber_moltypes(sys, mi);
647 convert_harmonics(*mi, atypes);
650 if (ir->eDisre == DistanceRestraintRefinement::None)
652 i = rm_interactions(F_DISRES, *mi);
655 set_warning_line(wi, "unknown", -1);
656 std::string warningMessage =
657 gmx::formatString("disre = no, removed %d distance restraints", i);
658 warning_note(wi, warningMessage.c_str());
663 i = rm_interactions(F_ORIRES, *mi);
666 set_warning_line(wi, "unknown", -1);
667 std::string warningMessage =
668 gmx::formatString("orire = no, removed %d orientation restraints", i);
669 warning_note(wi, warningMessage.c_str());
673 /* Copy structures from msys to sys */
674 molinfo2mtop(*mi, sys);
678 /* Discourage using the, previously common, setup of applying constraints
679 * to all bonds with force fields that have been parametrized with H-bond
680 * constraints only. Do not print note with large timesteps or vsites.
682 if (opts->nshake == eshALLBONDS && ffParametrizedWithHBondConstraints && ir->delta_t < 0.0026
683 && gmx_mtop_ftype_count(*sys, F_VSITE3FD) == 0)
685 set_warning_line(wi, "unknown", -1);
687 "You are using constraints on all bonds, whereas the forcefield "
688 "has been parametrized only with constraints involving hydrogen atoms. "
689 "We suggest using constraints = h-bonds instead, this will also improve "
693 /* COORDINATE file processing */
696 GMX_LOG(logger.info).asParagraph().appendTextFormatted("processing coordinates...");
703 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
704 state->natoms = conftop->atoms.nr;
705 if (state->natoms != sys->natoms)
708 "number of coordinates in coordinate file (%s, %d)\n"
709 " does not match topology (%s, %d)",
715 /* It would be nice to get rid of the copies below, but we don't know
716 * a priori if the number of atoms in confin matches what we expect.
718 state->flags |= enumValueToBitMask(StateEntry::X);
721 state->flags |= enumValueToBitMask(StateEntry::V);
723 state_change_natoms(state, state->natoms);
724 std::copy(x, x + state->natoms, state->x.data());
728 std::copy(v, v + state->natoms, state->v.data());
731 /* This call fixes the box shape for runs with pressure scaling */
732 set_box_rel(ir, state);
734 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms, logger);
740 std::string warningMessage = gmx::formatString(
741 "%d non-matching atom name%s\n"
742 "atom names from %s will be used\n"
743 "atom names from %s will be ignored\n",
745 (nmismatch == 1) ? "" : "s",
748 warning(wi, warningMessage.c_str());
751 /* Do more checks, mostly related to constraints */
756 .appendTextFormatted("double-checking input for internal consistency...");
759 bool bHasNormalConstraints =
760 0 < (nint_ftype(sys, *mi, F_CONSTR) + nint_ftype(sys, *mi, F_CONSTRNC));
761 bool bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, *mi, F_SETTLE);
762 double_check(ir, state->box, bHasNormalConstraints, bHasAnyConstraints, wi);
769 snew(mass, state->natoms);
770 for (const AtomProxy atomP : AtomRange(*sys))
772 const t_atom& local = atomP.atom();
773 int i = atomP.globalAtomNumber();
777 if (opts->seed == -1)
779 opts->seed = static_cast<int>(gmx::makeRandomSeed());
780 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Setting gen_seed to %d", opts->seed);
782 state->flags |= enumValueToBitMask(StateEntry::V);
783 maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array(), logger);
785 stop_cm(logger, state->natoms, mass, state->x.rvec_array(), state->v.rvec_array());
790 static void copy_state(const char* slog, t_trxframe* fr, bool bReadVel, t_state* state, double* use_time)
792 if (fr->not_ok & FRAME_NOT_OK)
794 gmx_fatal(FARGS, "Can not start from an incomplete frame");
798 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s", slog);
801 std::copy(fr->x, fr->x + state->natoms, state->x.data());
806 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
808 std::copy(fr->v, fr->v + state->natoms, state->v.data());
812 copy_mat(fr->box, state->box);
815 *use_time = fr->time;
818 static void cont_status(const char* slog,
826 const gmx_output_env_t* oenv,
827 const gmx::MDLogger& logger)
828 /* If fr_time == -1 read the last frame available which is complete */
835 bReadVel = (bNeedVel && !bGenVel);
839 .appendTextFormatted("Reading Coordinates%s and Box size from old trajectory",
840 bReadVel ? ", Velocities" : "");
843 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read whole trajectory");
847 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read till time %g", fr_time);
855 .appendTextFormatted(
856 "Velocities generated: "
857 "ignoring velocities in input trajectory");
859 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
863 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
867 GMX_LOG(logger.warning)
869 .appendTextFormatted(
870 "WARNING: Did not find a frame with velocities in file %s,\n"
871 " all velocities will be set to zero!",
873 for (auto& vi : makeArrayRef(state->v))
878 /* Search for a frame without velocities */
880 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
884 state->natoms = fr.natoms;
886 if (sys->natoms != state->natoms)
889 "Number of atoms in Topology "
890 "is not the same as in Trajectory");
892 copy_state(slog, &fr, bReadVel, state, &use_time);
894 /* Find the appropriate frame */
895 while ((fr_time == -1 || fr.time < fr_time) && read_next_frame(oenv, fp, &fr))
897 copy_state(slog, &fr, bReadVel, state, &use_time);
902 /* Set the relative box lengths for preserving the box shape.
903 * Note that this call can lead to differences in the last bit
904 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
906 set_box_rel(ir, state);
908 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Using frame at t = %g ps", use_time);
909 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Starting time for run is %g ps", ir->init_t);
911 if ((ir->epc != PressureCoupling::No || ir->etc == TemperatureCoupling::NoseHoover) && ener)
913 get_enx_state(ener, use_time, sys->groups, ir, state);
914 preserve_box_shape(ir, state->box_rel, state->boxv);
918 static void read_posres(gmx_mtop_t* mtop,
919 gmx::ArrayRef<const MoleculeInformation> molinfo,
922 RefCoordScaling rc_scaling,
926 const gmx::MDLogger& logger)
934 int natoms, npbcdim = 0;
939 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
940 natoms = top->atoms.nr;
943 if (natoms != mtop->natoms)
945 std::string warningMessage = gmx::formatString(
946 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
947 "(%d). Will assume that the first %d atoms in the topology and %s match.",
951 std::min(mtop->natoms, natoms),
953 warning(wi, warningMessage.c_str());
956 npbcdim = numPbcDimensions(pbcType);
957 GMX_RELEASE_ASSERT(npbcdim <= DIM, "Invalid npbcdim");
959 if (rc_scaling != RefCoordScaling::No)
961 copy_mat(box, invbox);
962 for (int j = npbcdim; j < DIM; j++)
964 clear_rvec(invbox[j]);
967 gmx::invertBoxMatrix(invbox, invbox);
970 /* Copy the reference coordinates to mtop */
974 snew(hadAtom, natoms);
975 for (gmx_molblock_t& molb : mtop->molblock)
977 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
978 const InteractionsOfType* pr = &(molinfo[molb.type].interactions[F_POSRES]);
979 const InteractionsOfType* prfb = &(molinfo[molb.type].interactions[F_FBPOSRES]);
980 if (pr->size() > 0 || prfb->size() > 0)
982 atom = mtop->moltype[molb.type].atoms.atom;
983 for (const auto& restraint : pr->interactionTypes)
985 int ai = restraint.ai();
989 "Position restraint atom index (%d) in moltype '%s' is larger than "
990 "number of atoms in %s (%d).\n",
992 *molinfo[molb.type].name,
997 if (rc_scaling == RefCoordScaling::Com)
999 /* Determine the center of mass of the posres reference coordinates */
1000 for (int j = 0; j < npbcdim; j++)
1002 sum[j] += atom[ai].m * x[a + ai][j];
1004 totmass += atom[ai].m;
1007 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
1008 for (const auto& restraint : prfb->interactionTypes)
1010 int ai = restraint.ai();
1014 "Position restraint atom index (%d) in moltype '%s' is larger than "
1015 "number of atoms in %s (%d).\n",
1017 *molinfo[molb.type].name,
1021 if (rc_scaling == RefCoordScaling::Com && !hadAtom[ai])
1023 /* Determine the center of mass of the posres reference coordinates */
1024 for (int j = 0; j < npbcdim; j++)
1026 sum[j] += atom[ai].m * x[a + ai][j];
1028 totmass += atom[ai].m;
1033 molb.posres_xA.resize(nat_molb);
1034 for (int i = 0; i < nat_molb; i++)
1036 copy_rvec(x[a + i], molb.posres_xA[i]);
1041 molb.posres_xB.resize(nat_molb);
1042 for (int i = 0; i < nat_molb; i++)
1044 copy_rvec(x[a + i], molb.posres_xB[i]);
1050 if (rc_scaling == RefCoordScaling::Com)
1054 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
1056 for (int j = 0; j < npbcdim; j++)
1058 com[j] = sum[j] / totmass;
1060 GMX_LOG(logger.info)
1062 .appendTextFormatted(
1063 "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f",
1069 if (rc_scaling != RefCoordScaling::No)
1071 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
1073 for (gmx_molblock_t& molb : mtop->molblock)
1075 nat_molb = molb.nmol * mtop->moltype[molb.type].atoms.nr;
1076 if (!molb.posres_xA.empty() || !molb.posres_xB.empty())
1078 std::vector<gmx::RVec>& xp = (!bTopB ? molb.posres_xA : molb.posres_xB);
1079 for (int i = 0; i < nat_molb; i++)
1081 for (int j = 0; j < npbcdim; j++)
1083 if (rc_scaling == RefCoordScaling::All)
1085 /* Convert from Cartesian to crystal coordinates */
1086 xp[i][j] *= invbox[j][j];
1087 for (int k = j + 1; k < npbcdim; k++)
1089 xp[i][j] += invbox[k][j] * xp[i][k];
1092 else if (rc_scaling == RefCoordScaling::Com)
1094 /* Subtract the center of mass */
1102 if (rc_scaling == RefCoordScaling::Com)
1104 /* Convert the COM from Cartesian to crystal coordinates */
1105 for (int j = 0; j < npbcdim; j++)
1107 com[j] *= invbox[j][j];
1108 for (int k = j + 1; k < npbcdim; k++)
1110 com[j] += invbox[k][j] * com[k];
1121 static void gen_posres(gmx_mtop_t* mtop,
1122 gmx::ArrayRef<const MoleculeInformation> mi,
1125 RefCoordScaling rc_scaling,
1130 const gmx::MDLogger& logger)
1132 read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi, logger);
1133 /* It is safer to simply read the b-state posres rather than trying
1134 * to be smart and copy the positions.
1136 read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi, logger);
1139 static void set_wall_atomtype(PreprocessingAtomTypes* at,
1143 const gmx::MDLogger& logger)
1149 GMX_LOG(logger.info).asParagraph().appendTextFormatted("Searching the wall atom type(s)");
1151 for (i = 0; i < ir->nwall; i++)
1153 auto atomType = at->atomTypeFromName(opts->wall_atomtype[i]);
1154 if (!atomType.has_value())
1156 std::string warningMessage = gmx::formatString(
1157 "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1158 warning_error(wi, warningMessage.c_str());
1162 ir->wall_atomtype[i] = *atomType;
1167 static int nrdf_internal(const t_atoms* atoms)
1172 for (i = 0; i < atoms->nr; i++)
1174 /* Vsite ptype might not be set here yet, so also check the mass */
1175 if ((atoms->atom[i].ptype == ParticleType::Atom || atoms->atom[i].ptype == ParticleType::Nucleus)
1176 && atoms->atom[i].m > 0)
1183 case 0: // Fall through intended
1184 case 1: nrdf = 0; break;
1185 case 2: nrdf = 1; break;
1186 default: nrdf = nmass * 3 - 6; break;
1192 static void spline1d(double dx, const double* y, int n, double* u, double* y2)
1200 for (i = 1; i < n - 1; i++)
1202 p = 0.5 * y2[i - 1] + 2.0;
1204 q = (y[i + 1] - 2.0 * y[i] + y[i - 1]) / dx;
1205 u[i] = (3.0 * q / dx - 0.5 * u[i - 1]) / p;
1210 for (i = n - 2; i >= 0; i--)
1212 y2[i] = y2[i] * y2[i + 1] + u[i];
1218 interpolate1d(double xmin, double dx, const double* ya, const double* y2a, double x, double* y, double* y1)
1223 ix = static_cast<int>((x - xmin) / dx);
1225 a = (xmin + (ix + 1) * dx - x) / dx;
1226 b = (x - xmin - ix * dx) / dx;
1228 *y = a * ya[ix] + b * ya[ix + 1]
1229 + ((a * a * a - a) * y2a[ix] + (b * b * b - b) * y2a[ix + 1]) * (dx * dx) / 6.0;
1230 *y1 = (ya[ix + 1] - ya[ix]) / dx - (3.0 * a * a - 1.0) / 6.0 * dx * y2a[ix]
1231 + (3.0 * b * b - 1.0) / 6.0 * dx * y2a[ix + 1];
1235 static void setup_cmap(int grid_spacing, int nc, gmx::ArrayRef<const real> grid, gmx_cmap_t* cmap_grid)
1237 int i, j, k, ii, jj, kk, idx;
1239 double dx, xmin, v, v1, v2, v12;
1242 std::vector<double> tmp_u(2 * grid_spacing, 0.0);
1243 std::vector<double> tmp_u2(2 * grid_spacing, 0.0);
1244 std::vector<double> tmp_yy(2 * grid_spacing, 0.0);
1245 std::vector<double> tmp_y1(2 * grid_spacing, 0.0);
1246 std::vector<double> tmp_t2(2 * grid_spacing * 2 * grid_spacing, 0.0);
1247 std::vector<double> tmp_grid(2 * grid_spacing * 2 * grid_spacing, 0.0);
1249 dx = 360.0 / grid_spacing;
1250 xmin = -180.0 - dx * grid_spacing / 2;
1252 for (kk = 0; kk < nc; kk++)
1254 /* Compute an offset depending on which cmap we are using
1255 * Offset will be the map number multiplied with the
1256 * grid_spacing * grid_spacing * 2
1258 offset = kk * grid_spacing * grid_spacing * 2;
1260 for (i = 0; i < 2 * grid_spacing; i++)
1262 ii = (i + grid_spacing - grid_spacing / 2) % grid_spacing;
1264 for (j = 0; j < 2 * grid_spacing; j++)
1266 jj = (j + grid_spacing - grid_spacing / 2) % grid_spacing;
1267 tmp_grid[i * grid_spacing * 2 + j] = grid[offset + ii * grid_spacing + jj];
1271 for (i = 0; i < 2 * grid_spacing; i++)
1274 &(tmp_grid[2 * grid_spacing * i]),
1277 &(tmp_t2[2 * grid_spacing * i]));
1280 for (i = grid_spacing / 2; i < grid_spacing + grid_spacing / 2; i++)
1282 ii = i - grid_spacing / 2;
1283 phi = ii * dx - 180.0;
1285 for (j = grid_spacing / 2; j < grid_spacing + grid_spacing / 2; j++)
1287 jj = j - grid_spacing / 2;
1288 psi = jj * dx - 180.0;
1290 for (k = 0; k < 2 * grid_spacing; k++)
1294 &(tmp_grid[2 * grid_spacing * k]),
1295 &(tmp_t2[2 * grid_spacing * k]),
1301 spline1d(dx, tmp_yy.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1302 interpolate1d(xmin, dx, tmp_yy.data(), tmp_u2.data(), phi, &v, &v1);
1303 spline1d(dx, tmp_y1.data(), 2 * grid_spacing, tmp_u.data(), tmp_u2.data());
1304 interpolate1d(xmin, dx, tmp_y1.data(), tmp_u2.data(), phi, &v2, &v12);
1306 idx = ii * grid_spacing + jj;
1307 cmap_grid->cmapdata[kk].cmap[idx * 4] = grid[offset + ii * grid_spacing + jj];
1308 cmap_grid->cmapdata[kk].cmap[idx * 4 + 1] = v1;
1309 cmap_grid->cmapdata[kk].cmap[idx * 4 + 2] = v2;
1310 cmap_grid->cmapdata[kk].cmap[idx * 4 + 3] = v12;
1316 static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
1320 cmap_grid->grid_spacing = grid_spacing;
1321 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
1323 cmap_grid->cmapdata.resize(ngrid);
1325 for (i = 0; i < ngrid; i++)
1327 cmap_grid->cmapdata[i].cmap.resize(4 * nelem);
1332 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
1334 int count, count_mol;
1337 for (const gmx_molblock_t& molb : mtop->molblock)
1340 gmx::ArrayRef<const InteractionsOfType> interactions = mi[molb.type].interactions;
1342 for (int i = 0; i < F_NRE; i++)
1346 count_mol += 3 * interactions[i].size();
1348 else if (interaction_function[i].flags & IF_CONSTRAINT)
1350 count_mol += interactions[i].size();
1354 if (count_mol > nrdf_internal(&mi[molb.type].atoms))
1356 std::string warningMessage = gmx::formatString(
1357 "Molecule type '%s' has %d constraints.\n"
1358 "For stability and efficiency there should not be more constraints than "
1359 "internal number of degrees of freedom: %d.\n",
1360 *mi[molb.type].name,
1362 nrdf_internal(&mi[molb.type].atoms));
1363 warning(wi, warningMessage.c_str());
1365 count += molb.nmol * count_mol;
1371 static real calc_temp(const gmx_mtop_t* mtop, const t_inputrec* ir, rvec* v)
1374 for (const AtomProxy atomP : AtomRange(*mtop))
1376 const t_atom& local = atomP.atom();
1377 int i = atomP.globalAtomNumber();
1378 sum_mv2 += local.m * norm2(v[i]);
1382 for (int g = 0; g < ir->opts.ngtc; g++)
1384 nrdf += ir->opts.nrdf[g];
1387 return sum_mv2 / (nrdf * gmx::c_boltz);
1390 static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
1398 for (i = 0; i < ir->opts.ngtc; i++)
1400 if (ir->opts.tau_t[i] < 0)
1406 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1412 std::string warningMessage = gmx::formatString(
1413 "Some temperature coupling groups do not use temperature coupling. We will assume "
1414 "their temperature is not more than %.3f K. If their temperature is higher, the "
1415 "energy error and the Verlet buffer might be underestimated.",
1417 warning(wi, warningMessage.c_str());
1423 /* Checks if there are unbound atoms in moleculetype molt.
1424 * Prints a note for each unbound atoms and a warning if any is present.
1426 static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1428 const t_atoms* atoms = &molt->atoms;
1432 /* Only one atom, there can't be unbound atoms */
1436 std::vector<int> count(atoms->nr, 0);
1438 for (int ftype = 0; ftype < F_NRE; ftype++)
1440 if (((interaction_function[ftype].flags & IF_BOND) && NRAL(ftype) == 2 && ftype != F_CONNBONDS)
1441 || (interaction_function[ftype].flags & IF_CONSTRAINT) || ftype == F_SETTLE)
1443 const InteractionList& il = molt->ilist[ftype];
1444 const int nral = NRAL(ftype);
1446 for (int i = 0; i < il.size(); i += 1 + nral)
1448 for (int j = 0; j < nral; j++)
1450 count[il.iatoms[i + 1 + j]]++;
1456 int numDanglingAtoms = 0;
1457 for (int a = 0; a < atoms->nr; a++)
1459 if (atoms->atom[a].ptype != ParticleType::VSite && count[a] == 0)
1463 GMX_LOG(logger.warning)
1465 .appendTextFormatted(
1466 "Atom %d '%s' in moleculetype '%s' is not bound by a potential or "
1467 "constraint to any other atom in the same moleculetype.",
1469 *atoms->atomname[a],
1476 if (numDanglingAtoms > 0)
1478 std::string warningMessage = gmx::formatString(
1479 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
1480 "other atom in the same moleculetype. Although technically this might not cause "
1481 "issues in a simulation, this often means that the user forgot to add a "
1482 "bond/potential/constraint or put multiple molecules in the same moleculetype "
1483 "definition by mistake. Run with -v to get information for each atom.",
1486 warning_note(wi, warningMessage.c_str());
1490 /* Checks all moleculetypes for unbound atoms */
1491 static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
1493 for (const gmx_moltype_t& molt : mtop->moltype)
1495 checkForUnboundAtoms(&molt, bVerbose, wi, logger);
1499 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1501 * The specific decoupled modes this routine check for are angle modes
1502 * where the two bonds are constrained and the atoms a both ends are only
1503 * involved in a single constraint; the mass of the two atoms needs to
1504 * differ by more than \p massFactorThreshold.
1506 static bool haveDecoupledModeInMol(const gmx_moltype_t& molt,
1507 gmx::ArrayRef<const t_iparams> iparams,
1508 real massFactorThreshold)
1510 if (molt.ilist[F_CONSTR].empty() && molt.ilist[F_CONSTRNC].empty())
1515 const t_atom* atom = molt.atoms.atom;
1517 const auto atomToConstraints =
1518 gmx::make_at2con(molt, iparams, gmx::FlexibleConstraintTreatment::Exclude);
1520 bool haveDecoupledMode = false;
1521 for (int ftype = 0; ftype < F_NRE; ftype++)
1523 if (interaction_function[ftype].flags & IF_ATYPE)
1525 const int nral = NRAL(ftype);
1526 const InteractionList& il = molt.ilist[ftype];
1527 for (int i = 0; i < il.size(); i += 1 + nral)
1529 /* Here we check for the mass difference between the atoms
1530 * at both ends of the angle, that the atoms at the ends have
1531 * 1 contraint and the atom in the middle at least 3; we check
1532 * that the 3 atoms are linked by constraints below.
1533 * We check for at least three constraints for the middle atom,
1534 * since with only the two bonds in the angle, we have 3-atom
1535 * molecule, which has much more thermal exhange in this single
1536 * angle mode than molecules with more atoms.
1537 * Note that this check also catches molecules where more atoms
1538 * are connected to one or more atoms in the angle, but by
1539 * bond potentials instead of angles. But such cases will not
1540 * occur in "normal" molecules and it doesn't hurt running
1541 * those with higher accuracy settings as well.
1543 int a0 = il.iatoms[1 + i];
1544 int a1 = il.iatoms[1 + i + 1];
1545 int a2 = il.iatoms[1 + i + 2];
1546 if ((atom[a0].m > atom[a2].m * massFactorThreshold || atom[a2].m > atom[a0].m * massFactorThreshold)
1547 && atomToConstraints[a0].ssize() == 1 && atomToConstraints[a2].ssize() == 1
1548 && atomToConstraints[a1].ssize() >= 3)
1550 int constraint0 = atomToConstraints[a0][0];
1551 int constraint2 = atomToConstraints[a2][0];
1553 bool foundAtom0 = false;
1554 bool foundAtom2 = false;
1555 for (const int constraint : atomToConstraints[a1])
1557 if (constraint == constraint0)
1561 if (constraint == constraint2)
1566 if (foundAtom0 && foundAtom2)
1568 haveDecoupledMode = true;
1575 return haveDecoupledMode;
1578 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1580 * When decoupled modes are present and the accuray in insufficient,
1581 * this routine issues a warning if the accuracy is insufficient.
1583 static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec* ir, warninp* wi)
1585 /* We only have issues with decoupled modes with normal MD.
1586 * With stochastic dynamics equipartitioning is enforced strongly.
1593 /* When atoms of very different mass are involved in an angle potential
1594 * and both bonds in the angle are constrained, the dynamic modes in such
1595 * angles have very different periods and significant energy exchange
1596 * takes several nanoseconds. Thus even a small amount of error in
1597 * different algorithms can lead to violation of equipartitioning.
1598 * The parameters below are mainly based on an all-atom chloroform model
1599 * with all bonds constrained. Equipartitioning is off by more than 1%
1600 * (need to run 5-10 ns) when the difference in mass between H and Cl
1601 * is more than a factor 13 and the accuracy is less than the thresholds
1602 * given below. This has been verified on some other molecules.
1604 * Note that the buffer and shake parameters have unit length and
1605 * energy/time, respectively, so they will "only" work correctly
1606 * for atomistic force fields using MD units.
1608 const real massFactorThreshold = 13.0;
1609 const real bufferToleranceThreshold = 1e-4;
1610 const int lincsIterationThreshold = 2;
1611 const int lincsOrderThreshold = 4;
1612 const real shakeToleranceThreshold = 0.005 * ir->delta_t;
1614 bool lincsWithSufficientTolerance =
1615 (ir->eConstrAlg == ConstraintAlgorithm::Lincs
1616 && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1617 bool shakeWithSufficientTolerance = (ir->eConstrAlg == ConstraintAlgorithm::Shake
1618 && ir->shake_tol <= 1.1 * shakeToleranceThreshold);
1619 if (ir->cutoff_scheme == CutoffScheme::Verlet && ir->verletbuf_tol <= 1.1 * bufferToleranceThreshold
1620 && (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1625 bool haveDecoupledMode = false;
1626 for (const gmx_moltype_t& molt : mtop->moltype)
1628 if (haveDecoupledModeInMol(molt, mtop->ffparams.iparams, massFactorThreshold))
1630 haveDecoupledMode = true;
1634 if (haveDecoupledMode)
1636 std::string message = gmx::formatString(
1637 "There are atoms at both ends of an angle, connected by constraints "
1638 "and with masses that differ by more than a factor of %g. This means "
1639 "that there are likely dynamic modes that are only very weakly coupled.",
1640 massFactorThreshold);
1641 if (ir->cutoff_scheme == CutoffScheme::Verlet)
1643 message += gmx::formatString(
1644 " To ensure good equipartitioning, you need to either not use "
1645 "constraints on all bonds (but, if possible, only on bonds involving "
1646 "hydrogens) or use integrator = %s or decrease one or more tolerances: "
1647 "verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order "
1648 ">= %d or SHAKE tolerance <= %g",
1649 enumValueToString(IntegrationAlgorithm::SD1),
1650 bufferToleranceThreshold,
1651 lincsIterationThreshold,
1652 lincsOrderThreshold,
1653 shakeToleranceThreshold);
1657 message += gmx::formatString(
1658 " To ensure good equipartitioning, we suggest to switch to the %s "
1659 "cutoff-scheme, since that allows for better control over the Verlet "
1660 "buffer size and thus over the energy drift.",
1661 enumValueToString(CutoffScheme::Verlet));
1663 warning(wi, message);
1667 static void set_verlet_buffer(const gmx_mtop_t* mtop,
1672 const gmx::MDLogger& logger)
1674 GMX_LOG(logger.info)
1676 .appendTextFormatted(
1677 "Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K",
1681 /* Calculate the buffer size for simple atom vs atoms list */
1682 VerletbufListSetup listSetup1x1;
1683 listSetup1x1.cluster_size_i = 1;
1684 listSetup1x1.cluster_size_j = 1;
1685 const real rlist_1x1 = calcVerletBufferSize(
1686 *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup1x1);
1688 /* Set the pair-list buffer size in ir */
1689 VerletbufListSetup listSetup4x4 = verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1690 ir->rlist = calcVerletBufferSize(
1691 *mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, buffer_temp, listSetup4x4);
1693 const int n_nonlin_vsite = gmx::countNonlinearVsites(*mtop);
1694 if (n_nonlin_vsite > 0)
1696 std::string warningMessage = gmx::formatString(
1697 "There are %d non-linear virtual site constructions. Their contribution to the "
1698 "energy error is approximated. In most cases this does not affect the error "
1701 warning_note(wi, warningMessage);
1704 GMX_LOG(logger.info)
1706 .appendTextFormatted(
1707 "Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm",
1711 rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
1713 GMX_LOG(logger.info)
1715 .appendTextFormatted(
1716 "Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm",
1717 listSetup4x4.cluster_size_i,
1718 listSetup4x4.cluster_size_j,
1720 ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
1722 GMX_LOG(logger.info)
1724 .appendTextFormatted(
1725 "Note that mdrun will redetermine rlist based on the actual pair-list setup");
1727 if (gmx::square(ir->rlist) >= max_cutoff2(ir->pbcType, box))
1730 "The pair-list cut-off (%g nm) is longer than half the shortest box vector or "
1731 "longer than the smallest box diagonal element (%g nm). Increase the box size or "
1732 "decrease nstlist or increase verlet-buffer-tolerance.",
1734 std::sqrt(max_cutoff2(ir->pbcType, box)));
1738 int gmx_grompp(int argc, char* argv[])
1740 const char* desc[] = {
1741 "[THISMODULE] (the gromacs preprocessor)",
1742 "reads a molecular topology file, checks the validity of the",
1743 "file, expands the topology from a molecular description to an atomic",
1744 "description. The topology file contains information about",
1745 "molecule types and the number of molecules, the preprocessor",
1746 "copies each molecule as needed. ",
1747 "There is no limitation on the number of molecule types. ",
1748 "Bonds and bond-angles can be converted into constraints, separately",
1749 "for hydrogens and heavy atoms.",
1750 "Then a coordinate file is read and velocities can be generated",
1751 "from a Maxwellian distribution if requested.",
1752 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1753 "(eg. number of MD steps, time step, cut-off).",
1754 "Eventually a binary file is produced that can serve as the sole input",
1755 "file for the MD program.[PAR]",
1757 "[THISMODULE] uses the atom names from the topology file. The atom names",
1758 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1759 "warnings when they do not match the atom names in the topology.",
1760 "Note that the atom names are irrelevant for the simulation as",
1761 "only the atom types are used for generating interaction parameters.[PAR]",
1763 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1764 "etc. The preprocessor supports the following keywords::",
1767 " #ifndef VARIABLE",
1770 " #define VARIABLE",
1772 " #include \"filename\"",
1773 " #include <filename>",
1775 "The functioning of these statements in your topology may be modulated by",
1776 "using the following two flags in your [REF].mdp[ref] file::",
1778 " define = -DVARIABLE1 -DVARIABLE2",
1779 " include = -I/home/john/doe",
1781 "For further information a C-programming textbook may help you out.",
1782 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1783 "topology file written out so that you can verify its contents.[PAR]",
1785 "When using position restraints, a file with restraint coordinates",
1786 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1787 "for [TT]-c[tt]). For free energy calculations, separate reference",
1788 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1789 "otherwise they will be equal to those of the A topology.[PAR]",
1791 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1792 "The last frame with coordinates and velocities will be read,",
1793 "unless the [TT]-time[tt] option is used. Only if this information",
1794 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1795 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1796 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1797 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1800 "[THISMODULE] can be used to restart simulations (preserving",
1801 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1802 "However, for simply changing the number of run steps to extend",
1803 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1804 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1805 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1806 "like output frequency, then supplying the checkpoint file to",
1807 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1808 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1809 "the ensemble (if possible) still requires passing the checkpoint",
1810 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1812 "By default, all bonded interactions which have constant energy due to",
1813 "virtual site constructions will be removed. If this constant energy is",
1814 "not zero, this will result in a shift in the total energy. All bonded",
1815 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1816 "all constraints for distances which will be constant anyway because",
1817 "of virtual site constructions will be removed. If any constraints remain",
1818 "which involve virtual sites, a fatal error will result.[PAR]",
1820 "To verify your run input file, please take note of all warnings",
1821 "on the screen, and correct where necessary. Do also look at the contents",
1822 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1823 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1824 "with the [TT]-debug[tt] option which will give you more information",
1825 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1826 "can see the contents of the run input file with the [gmx-dump]",
1827 "program. [gmx-check] can be used to compare the contents of two",
1828 "run input files.[PAR]",
1830 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1831 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1832 "harmless, but usually they are not. The user is advised to carefully",
1833 "interpret the output messages before attempting to bypass them with",
1836 std::vector<MoleculeInformation> mi;
1837 std::unique_ptr<MoleculeInformation> intermolecular_interactions;
1839 CombinationRule comb;
1842 const char* mdparin;
1843 bool bNeedVel, bGenVel;
1844 gmx_output_env_t* oenv;
1845 gmx_bool bVerbose = FALSE;
1848 t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
1849 { efMDP, "-po", "mdout", ffWRITE },
1850 { efSTX, "-c", nullptr, ffREAD },
1851 { efSTX, "-r", "restraint", ffOPTRD },
1852 { efSTX, "-rb", "restraint", ffOPTRD },
1853 { efNDX, nullptr, nullptr, ffOPTRD },
1854 { efTOP, nullptr, nullptr, ffREAD },
1855 { efTOP, "-pp", "processed", ffOPTWR },
1856 { efTPR, "-o", nullptr, ffWRITE },
1857 { efTRN, "-t", nullptr, ffOPTRD },
1858 { efEDR, "-e", nullptr, ffOPTRD },
1859 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1860 { efGRO, "-imd", "imdgroup", ffOPTWR },
1861 { efTRN, "-ref", "rotref", ffOPTRW | ffALLOW_MISSING },
1862 /* This group is needed by the QMMM MDModule: */
1863 { efQMI, "-qmi", nullptr, ffOPTRD } };
1864 #define NFILE asize(fnm)
1866 /* Command line options */
1867 gmx_bool bRenum = TRUE;
1868 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1872 { "-v", FALSE, etBOOL, { &bVerbose }, "Be loud and noisy" },
1873 { "-time", FALSE, etREAL, { &fr_time }, "Take frame at or first after this time." },
1878 "Remove constant bonded interactions with virtual sites" },
1883 "Number of allowed warnings during input processing. Not for normal use and may "
1884 "generate unstable systems" },
1889 "Set parameters for bonded interactions without defaults to zero instead of "
1890 "generating an error" },
1895 "Renumber atomtypes and minimize number of atomtypes" }
1898 /* Parse the command line */
1899 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
1904 /* Initiate some variables */
1905 gmx::MDModules mdModules;
1906 t_inputrec irInstance;
1907 t_inputrec* ir = &irInstance;
1908 t_gromppopts optsInstance;
1909 t_gromppopts* opts = &optsInstance;
1910 snew(opts->include, STRLEN);
1911 snew(opts->define, STRLEN);
1913 gmx::LoggerBuilder builder;
1914 builder.addTargetStream(gmx::MDLogger::LogLevel::Info, &gmx::TextOutputFile::standardOutput());
1915 builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
1916 gmx::LoggerOwner logOwner(builder.build());
1917 const gmx::MDLogger logger(logOwner.logger());
1920 wi = init_warning(TRUE, maxwarn);
1922 /* PARAMETER file processing */
1923 mdparin = opt2fn("-f", NFILE, fnm);
1924 set_warning_line(wi, mdparin, -1);
1927 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1929 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1931 // Now that the MDModules have their options assigned from get_ir, subscribe
1932 // to eventual notifications during pre-processing their data
1933 mdModules.subscribeToPreProcessingNotifications();
1935 // Notify MDModules of existing logger
1936 mdModules.notifiers().preProcessingNotifier_.notify(logger);
1938 // Notify MDModules of existing warninp
1939 mdModules.notifiers().preProcessingNotifier_.notify(wi);
1941 // Notify QMMM MDModule of external QM input file command-line option
1943 gmx::QMInputFileName qmInputFileName = { ftp2bSet(efQMI, NFILE, fnm), ftp2fn(efQMI, NFILE, fnm) };
1944 mdModules.notifiers().preProcessingNotifier_.notify(qmInputFileName);
1949 GMX_LOG(logger.info)
1951 .appendTextFormatted("checking input for internal consistency...");
1953 check_ir(mdparin, mdModules.notifiers(), ir, opts, wi);
1955 if (ir->ld_seed == -1)
1957 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1958 GMX_LOG(logger.info)
1960 .appendTextFormatted("Setting the LD random seed to %" PRId64 "", ir->ld_seed);
1963 if (ir->expandedvals->lmc_seed == -1)
1965 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1966 GMX_LOG(logger.info)
1968 .appendTextFormatted("Setting the lambda MC random seed to %d", ir->expandedvals->lmc_seed);
1971 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1972 bGenVel = (bNeedVel && opts->bGenVel);
1973 if (bGenVel && ir->bContinuation)
1975 std::string warningMessage = gmx::formatString(
1976 "Generating velocities is inconsistent with attempting "
1977 "to continue a previous run. Choose only one of "
1978 "gen-vel = yes and continuation = yes.");
1979 warning_error(wi, warningMessage);
1982 std::array<InteractionsOfType, F_NRE> interactions;
1984 PreprocessingAtomTypes atypes;
1987 pr_symtab(debug, 0, "Just opened", &sys.symtab);
1990 const char* fn = ftp2fn(efTOP, NFILE, fnm);
1991 if (!gmx_fexist(fn))
1993 gmx_fatal(FARGS, "%s does not exist", fn);
1998 opt2fn_null("-pp", NFILE, fnm),
1999 opt2fn("-c", NFILE, fnm),
2009 &intermolecular_interactions,
2020 pr_symtab(debug, 0, "After new_status", &sys.symtab);
2024 /* set parameters for virtual site construction (not for vsiten) */
2025 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2027 nvsite += set_vsites(bVerbose, &sys.moltype[mt].atoms, &atypes, mi[mt].interactions, logger);
2029 /* now throw away all obsolete bonds, angles and dihedrals: */
2030 /* note: constraints are ALWAYS removed */
2033 for (size_t mt = 0; mt < sys.moltype.size(); mt++)
2035 clean_vsite_bondeds(mi[mt].interactions, sys.moltype[mt].atoms.nr, bRmVSBds, logger);
2039 if ((count_constraints(&sys, mi, wi) != 0) && (ir->eConstrAlg == ConstraintAlgorithm::Shake))
2041 if (ir->eI == IntegrationAlgorithm::CG || ir->eI == IntegrationAlgorithm::LBFGS)
2043 std::string warningMessage =
2044 gmx::formatString("Can not do %s with %s, use %s",
2045 enumValueToString(ir->eI),
2046 enumValueToString(ConstraintAlgorithm::Shake),
2047 enumValueToString(ConstraintAlgorithm::Lincs));
2048 warning_error(wi, warningMessage);
2050 if (ir->bPeriodicMols)
2052 std::string warningMessage =
2053 gmx::formatString("Can not do periodic molecules with %s, use %s",
2054 enumValueToString(ConstraintAlgorithm::Shake),
2055 enumValueToString(ConstraintAlgorithm::Lincs));
2056 warning_error(wi, warningMessage);
2060 if (EI_SD(ir->eI) && ir->etc != TemperatureCoupling::No)
2062 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
2065 /* Check for errors in the input now, since they might cause problems
2066 * during processing further down.
2068 check_warning_error(wi, FARGS);
2070 if (nint_ftype(&sys, mi, F_POSRES) > 0 || nint_ftype(&sys, mi, F_FBPOSRES) > 0)
2072 if (ir->epc == PressureCoupling::ParrinelloRahman || ir->epc == PressureCoupling::Mttk)
2074 std::string warningMessage = gmx::formatString(
2075 "You are combining position restraints with %s pressure coupling, which can "
2076 "lead to instabilities. If you really want to combine position restraints with "
2077 "pressure coupling, we suggest to use %s pressure coupling instead.",
2078 enumValueToString(ir->epc),
2079 enumValueToString(PressureCoupling::Berendsen));
2080 warning_note(wi, warningMessage);
2083 const char* fn = opt2fn("-r", NFILE, fnm);
2086 if (!gmx_fexist(fn))
2089 "Cannot find position restraint file %s (option -r).\n"
2090 "From GROMACS-2018, you need to specify the position restraint "
2091 "coordinate files explicitly to avoid mistakes, although you can "
2092 "still use the same file as you specify for the -c option.",
2096 if (opt2bSet("-rb", NFILE, fnm))
2098 fnB = opt2fn("-rb", NFILE, fnm);
2099 if (!gmx_fexist(fnB))
2102 "Cannot find B-state position restraint file %s (option -rb).\n"
2103 "From GROMACS-2018, you need to specify the position restraint "
2104 "coordinate files explicitly to avoid mistakes, although you can "
2105 "still use the same file as you specify for the -c option.",
2116 std::string message = gmx::formatString("Reading position restraint coords from %s", fn);
2117 if (strcmp(fn, fnB) != 0)
2119 message += gmx::formatString(" and %s", fnB);
2121 GMX_LOG(logger.info).asParagraph().appendText(message);
2123 gen_posres(&sys, mi, fn, fnB, ir->refcoord_scaling, ir->pbcType, ir->posres_com, ir->posres_comB, wi, logger);
2126 /* If we are using CMAP, setup the pre-interpolation grid */
2127 if (interactions[F_CMAP].ncmap() > 0)
2129 init_cmap_grid(&sys.ffparams.cmap_grid,
2130 interactions[F_CMAP].cmapAngles,
2131 interactions[F_CMAP].cmakeGridSpacing);
2132 setup_cmap(interactions[F_CMAP].cmakeGridSpacing,
2133 interactions[F_CMAP].cmapAngles,
2134 interactions[F_CMAP].cmap,
2135 &sys.ffparams.cmap_grid);
2138 set_wall_atomtype(&atypes, opts, ir, wi, logger);
2141 atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
2144 if (ir->implicit_solvent)
2146 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
2149 /* PELA: Copy the atomtype data to the topology atomtype list */
2150 atypes.copyTot_atomtypes(&(sys.atomtypes));
2154 pr_symtab(debug, 0, "After atype.renumberTypes", &sys.symtab);
2159 GMX_LOG(logger.info).asParagraph().appendTextFormatted("converting bonded parameters...");
2162 const int ntype = atypes.size();
2163 convertInteractionsOfType(
2164 ntype, interactions, mi, intermolecular_interactions.get(), comb, reppow, fudgeQQ, &sys);
2168 pr_symtab(debug, 0, "After converInteractionsOfType", &sys.symtab);
2171 /* set ptype to VSite for virtual sites */
2172 for (gmx_moltype_t& moltype : sys.moltype)
2174 set_vsites_ptype(FALSE, &moltype, logger);
2178 pr_symtab(debug, 0, "After virtual sites", &sys.symtab);
2180 /* Check velocity for virtual sites and shells */
2183 check_vel(&sys, state.v.rvec_array());
2186 /* check for shells and inpurecs */
2187 check_shells_inputrec(&sys, ir, wi);
2190 check_mol(&sys, wi);
2192 if (haveFepPerturbedMassesInSettles(sys))
2195 "SETTLE is not implemented for atoms whose mass is perturbed. "
2196 "You might instead use normal constraints.");
2199 checkForUnboundAtoms(&sys, bVerbose, wi, logger);
2201 if (EI_DYNAMICS(ir->eI) && ir->eI != IntegrationAlgorithm::BD)
2203 check_bonds_timestep(&sys, ir->delta_t, wi);
2206 checkDecoupledModeAccuracy(&sys, ir, wi);
2208 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2212 "Zero-step energy minimization will alter the coordinates before calculating the "
2213 "energy. If you just want the energy of a single point, try zero-step MD (with "
2214 "unconstrained_start = yes). To do multiple single-point energy evaluations of "
2215 "different configurations of the same topology, use mdrun -rerun.");
2218 check_warning_error(wi, FARGS);
2222 GMX_LOG(logger.info).asParagraph().appendTextFormatted("initialising group options...");
2224 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifiers(), ir, wi);
2226 // Notify topology to MdModules for pre-processing after all indexes were built
2227 mdModules.notifiers().preProcessingNotifier_.notify(&sys);
2229 if (ir->cutoff_scheme == CutoffScheme::Verlet && ir->verletbuf_tol > 0)
2231 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2235 if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No)
2239 buffer_temp = opts->tempi;
2243 buffer_temp = calc_temp(&sys, ir, state.v.rvec_array());
2245 if (buffer_temp > 0)
2247 std::string warningMessage = gmx::formatString(
2248 "NVE simulation: will use the initial temperature of %.3f K for "
2249 "determining the Verlet buffer size",
2251 warning_note(wi, warningMessage);
2255 std::string warningMessage = gmx::formatString(
2256 "NVE simulation with an initial temperature of zero: will use a Verlet "
2257 "buffer of %d%%. Check your energy drift!",
2258 gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
2259 warning_note(wi, warningMessage);
2264 buffer_temp = get_max_reference_temp(ir, wi);
2267 if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No && buffer_temp == 0)
2269 /* NVE with initial T=0: we add a fixed ratio to rlist.
2270 * Since we don't actually use verletbuf_tol, we set it to -1
2271 * so it can't be misused later.
2273 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2274 ir->verletbuf_tol = -1;
2278 /* We warn for NVE simulations with a drift tolerance that
2279 * might result in a 1(.1)% drift over the total run-time.
2280 * Note that we can't warn when nsteps=0, since we don't
2281 * know how many steps the user intends to run.
2283 if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No && ir->nstlist > 1 && ir->nsteps > 0)
2285 const real driftTolerance = 0.01;
2286 /* We use 2 DOF per atom = 2kT pot+kin energy,
2287 * to be on the safe side with constraints.
2289 const real totalEnergyDriftPerAtomPerPicosecond =
2290 2 * gmx::c_boltz * buffer_temp / (ir->nsteps * ir->delta_t);
2292 if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
2294 std::string warningMessage = gmx::formatString(
2295 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
2296 "NVE simulation of length %g ps, which can give a final drift of "
2297 "%d%%. For conserving energy to %d%% when using constraints, you "
2298 "might need to set verlet-buffer-tolerance to %.1e.",
2300 ir->nsteps * ir->delta_t,
2301 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
2302 gmx::roundToInt(100 * driftTolerance),
2303 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
2304 warning_note(wi, warningMessage);
2308 set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi, logger);
2313 /* Init the temperature coupling state */
2314 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2318 pr_symtab(debug, 0, "After index", &sys.symtab);
2321 triple_check(mdparin, ir, &sys, wi);
2322 close_symtab(&sys.symtab);
2325 pr_symtab(debug, 0, "After close", &sys.symtab);
2328 if (ir->eI == IntegrationAlgorithm::Mimic)
2330 generate_qmexcl(&sys, ir, logger);
2333 if (ftp2bSet(efTRN, NFILE, fnm))
2337 GMX_LOG(logger.info)
2339 .appendTextFormatted("getting data from old trajectory ...");
2341 cont_status(ftp2fn(efTRN, NFILE, fnm),
2342 ftp2fn_null(efEDR, NFILE, fnm),
2353 if (ir->pbcType == PbcType::XY && ir->nwall != 2)
2355 clear_rvec(state.box[ZZ]);
2358 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2360 /* Calculate the optimal grid dimensions */
2362 EwaldBoxZScaler boxScaler(inputrecPbcXY2Walls(ir), ir->wall_ewald_zfac);
2363 boxScaler.scaleBox(state.box, scaledBox);
2365 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2367 /* Mark fourier_spacing as not used */
2368 ir->fourier_spacing = 0;
2370 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2372 set_warning_line(wi, mdparin, -1);
2374 wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2376 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2377 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize, &(ir->nkx), &(ir->nky), &(ir->nkz));
2378 if (ir->nkx < minGridSize || ir->nky < minGridSize || ir->nkz < minGridSize)
2381 "The PME grid size should be >= 2*(pme-order - 1); either manually "
2382 "increase the grid size or decrease pme-order");
2386 /* MRS: eventually figure out better logic for initializing the fep
2387 values that makes declaring the lambda and declaring the state not
2388 potentially conflict if not handled correctly. */
2389 if (ir->efep != FreeEnergyPerturbationType::No)
2391 state.fep_state = ir->fepvals->init_fep_state;
2392 for (i = 0; i < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); i++)
2394 /* init_lambda trumps state definitions*/
2395 if (ir->fepvals->init_lambda >= 0)
2397 state.lambda[i] = ir->fepvals->init_lambda;
2401 if (ir->fepvals->all_lambda[i].empty())
2403 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2407 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2413 pull_t* pull = nullptr;
2417 pull = set_pull_init(
2418 ir, sys, state.x, state.box, state.lambda[FreeEnergyPerturbationCouplingType::Mass], wi);
2421 /* Modules that supply external potential for pull coordinates
2422 * should register those potentials here. finish_pull() will check
2423 * that providers have been registerd for all external potentials.
2428 tensor compressibility = { { 0 } };
2429 if (ir->epc != PressureCoupling::No)
2431 copy_mat(ir->compress, compressibility);
2433 setStateDependentAwhParams(
2434 ir->awhParams.get(),
2441 ir->efep != FreeEnergyPerturbationType::No ? ir->fepvals->all_lambda[static_cast<int>(
2442 FreeEnergyPerturbationCouplingType::Fep)][ir->fepvals->init_fep_state]
2455 set_reference_positions(ir->rot,
2456 state.x.rvec_array(),
2458 opt2fn("-ref", NFILE, fnm),
2459 opt2bSet("-ref", NFILE, fnm),
2463 /* reset_multinr(sys); */
2465 if (EEL_PME(ir->coulombtype))
2467 float ratio = pme_load_estimate(sys, *ir, state.box);
2468 GMX_LOG(logger.info)
2470 .appendTextFormatted(
2471 "Estimate for the relative computational load of the PME mesh part: %.2f", ratio);
2472 /* With free energy we might need to do PME both for the A and B state
2473 * charges. This will double the cost, but the optimal performance will
2474 * then probably be at a slightly larger cut-off and grid spacing.
2476 if ((ir->efep == FreeEnergyPerturbationType::No && ratio > 1.0 / 2.0)
2477 || (ir->efep != FreeEnergyPerturbationType::No && ratio > 2.0 / 3.0))
2481 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2482 "and for highly parallel simulations between 0.25 and 0.33,\n"
2483 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2484 if (ir->efep != FreeEnergyPerturbationType::No)
2487 "For free energy simulations, the optimal load limit increases from "
2494 double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
2495 std::string warningMessage =
2496 gmx::formatString("This run will generate roughly %.0f Mb of data", cio);
2499 set_warning_line(wi, mdparin, -1);
2500 warning_note(wi, warningMessage);
2504 GMX_LOG(logger.info).asParagraph().appendText(warningMessage);
2508 // Hand over box and coordiantes to MdModules before they evaluate their final parameters
2510 gmx::CoordinatesAndBoxPreprocessed coordinatesAndBoxPreprocessed;
2511 coordinatesAndBoxPreprocessed.coordinates_ = state.x.arrayRefWithPadding();
2512 copy_mat(state.box, coordinatesAndBoxPreprocessed.box_);
2513 coordinatesAndBoxPreprocessed.pbc_ = ir->pbcType;
2514 mdModules.notifiers().preProcessingNotifier_.notify(coordinatesAndBoxPreprocessed);
2517 // Add the md modules internal parameters that are not mdp options
2518 // e.g., atom indices
2521 gmx::KeyValueTreeBuilder internalParameterBuilder;
2522 mdModules.notifiers().preProcessingNotifier_.notify(internalParameterBuilder.rootObject());
2523 ir->internalParameters =
2524 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
2527 if (ir->comm_mode != ComRemovalAlgorithm::No)
2529 const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
2530 if (ir->nstcomm % nstglobalcomm != 0)
2535 "COM removal frequency is set to (%d).\n"
2536 "Other settings require a global communication frequency of %d.\n"
2537 "Note that this will require additional global communication steps,\n"
2538 "which will reduce performance when using multiple ranks.\n"
2539 "Consider setting nstcomm to a multiple of %d.",
2548 GMX_LOG(logger.info).asParagraph().appendTextFormatted("writing run input file...");
2551 done_warning(wi, FARGS);
2552 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, sys);
2554 /* Output IMD group, if bIMD is TRUE */
2555 gmx::write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2557 sfree(opts->define);
2558 sfree(opts->wall_atomtype[0]);
2559 sfree(opts->wall_atomtype[1]);
2560 sfree(opts->include);
2561 sfree(opts->couple_moltype);
2563 for (auto& mol : mi)
2565 // Some of the contents of molinfo have been stolen, so
2566 // fullCleanUp can't be called.
2567 mol.partialCleanUp();
2569 done_inputrec_strings();
2570 output_env_done(oenv);