2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
49 #include <sys/types.h>
51 #include "gromacs/awh/read-params.h"
52 #include "gromacs/commandline/pargs.h"
53 #include "gromacs/ewald/ewald-utils.h"
54 #include "gromacs/ewald/pme.h"
55 #include "gromacs/fft/calcgrid.h"
56 #include "gromacs/fileio/confio.h"
57 #include "gromacs/fileio/enxio.h"
58 #include "gromacs/fileio/tpxio.h"
59 #include "gromacs/fileio/trxio.h"
60 #include "gromacs/fileio/warninp.h"
61 #include "gromacs/gmxpreprocess/add_par.h"
62 #include "gromacs/gmxpreprocess/convparm.h"
63 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
64 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
65 #include "gromacs/gmxpreprocess/grompp-impl.h"
66 #include "gromacs/gmxpreprocess/notset.h"
67 #include "gromacs/gmxpreprocess/readir.h"
68 #include "gromacs/gmxpreprocess/tomorse.h"
69 #include "gromacs/gmxpreprocess/topio.h"
70 #include "gromacs/gmxpreprocess/toputil.h"
71 #include "gromacs/gmxpreprocess/vsite_parm.h"
72 #include "gromacs/imd/imd.h"
73 #include "gromacs/math/functions.h"
74 #include "gromacs/math/invertmatrix.h"
75 #include "gromacs/math/units.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/mdlib/calc_verletbuf.h"
78 #include "gromacs/mdlib/compute_io.h"
79 #include "gromacs/mdlib/constr.h"
80 #include "gromacs/mdlib/perf_est.h"
81 #include "gromacs/mdlib/sim_util.h"
82 #include "gromacs/mdrunutility/mdmodules.h"
83 #include "gromacs/mdtypes/inputrec.h"
84 #include "gromacs/mdtypes/md_enums.h"
85 #include "gromacs/mdtypes/nblist.h"
86 #include "gromacs/mdtypes/state.h"
87 #include "gromacs/pbcutil/boxutilities.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/random/seed.h"
91 #include "gromacs/topology/ifunc.h"
92 #include "gromacs/topology/mtop_util.h"
93 #include "gromacs/topology/symtab.h"
94 #include "gromacs/topology/topology.h"
95 #include "gromacs/trajectory/trajectoryframe.h"
96 #include "gromacs/utility/arraysize.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/futil.h"
101 #include "gromacs/utility/gmxassert.h"
102 #include "gromacs/utility/smalloc.h"
103 #include "gromacs/utility/snprintf.h"
105 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
110 /* For all the molecule types */
111 for (i = 0; i < nrmols; i++)
113 n += mols[i].plist[ifunc].nr;
114 mols[i].plist[ifunc].nr = 0;
119 static int check_atom_names(const char *fn1, const char *fn2,
120 gmx_mtop_t *mtop, const t_atoms *at)
122 int mb, m, i, j, nmismatch;
124 #define MAXMISMATCH 20
126 if (mtop->natoms != at->nr)
128 gmx_incons("comparing atom names");
133 for (mb = 0; mb < mtop->nmolblock; mb++)
135 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
136 for (m = 0; m < mtop->molblock[mb].nmol; m++)
138 for (j = 0; j < tat->nr; j++)
140 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
142 if (nmismatch < MAXMISMATCH)
145 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
146 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
148 else if (nmismatch == MAXMISMATCH)
150 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
162 static void check_eg_vs_cg(gmx_mtop_t *mtop)
164 int astart, mb, m, cg, j, firstj;
165 unsigned char firsteg, eg;
168 /* Go through all the charge groups and make sure all their
169 * atoms are in the same energy group.
173 for (mb = 0; mb < mtop->nmolblock; mb++)
175 molt = &mtop->moltype[mtop->molblock[mb].type];
176 for (m = 0; m < mtop->molblock[mb].nmol; m++)
178 for (cg = 0; cg < molt->cgs.nr; cg++)
180 /* Get the energy group of the first atom in this charge group */
181 firstj = astart + molt->cgs.index[cg];
182 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
183 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
185 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
188 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
189 firstj+1, astart+j+1, cg+1, *molt->name);
193 astart += molt->atoms.nr;
198 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
201 char warn_buf[STRLEN];
204 for (cg = 0; cg < cgs->nr; cg++)
206 maxsize = std::max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
209 if (maxsize > MAX_CHARGEGROUP_SIZE)
211 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
213 else if (maxsize > 10)
215 set_warning_line(wi, topfn, -1);
217 "The largest charge group contains %d atoms.\n"
218 "Since atoms only see each other when the centers of geometry of the charge groups they belong to are within the cut-off distance, too large charge groups can lead to serious cut-off artifacts.\n"
219 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
220 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
222 warning_note(wi, warn_buf);
226 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
228 /* This check is not intended to ensure accurate integration,
229 * rather it is to signal mistakes in the mdp settings.
230 * A common mistake is to forget to turn on constraints
231 * for MD after energy minimization with flexible bonds.
232 * This check can also detect too large time steps for flexible water
233 * models, but such errors will often be masked by the constraints
234 * mdp options, which turns flexible water into water with bond constraints,
235 * but without an angle constraint. Unfortunately such incorrect use
236 * of water models can not easily be detected without checking
237 * for specific model names.
239 * The stability limit of leap-frog or velocity verlet is 4.44 steps
240 * per oscillational period.
241 * But accurate bonds distributions are lost far before that limit.
242 * To allow relatively common schemes (although not common with Gromacs)
243 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
244 * we set the note limit to 10.
246 int min_steps_warn = 5;
247 int min_steps_note = 10;
250 gmx_moltype_t *moltype, *w_moltype;
252 t_ilist *ilist, *ilb, *ilc, *ils;
254 int i, a1, a2, w_a1, w_a2, j;
255 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
256 gmx_bool bFound, bWater, bWarn;
257 char warn_buf[STRLEN];
259 ip = mtop->ffparams.iparams;
261 twopi2 = gmx::square(2*M_PI);
263 limit2 = gmx::square(min_steps_note*dt);
269 for (molt = 0; molt < mtop->nmoltype; molt++)
271 moltype = &mtop->moltype[molt];
272 atom = moltype->atoms.atom;
273 ilist = moltype->ilist;
274 ilc = &ilist[F_CONSTR];
275 ils = &ilist[F_SETTLE];
276 for (ftype = 0; ftype < F_NRE; ftype++)
278 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
284 for (i = 0; i < ilb->nr; i += 3)
286 fc = ip[ilb->iatoms[i]].harmonic.krA;
287 re = ip[ilb->iatoms[i]].harmonic.rA;
288 if (ftype == F_G96BONDS)
290 /* Convert squared sqaure fc to harmonic fc */
293 a1 = ilb->iatoms[i+1];
294 a2 = ilb->iatoms[i+2];
297 if (fc > 0 && m1 > 0 && m2 > 0)
299 period2 = twopi2*m1*m2/((m1 + m2)*fc);
303 period2 = GMX_FLOAT_MAX;
307 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
308 fc, m1, m2, std::sqrt(period2));
310 if (period2 < limit2)
313 for (j = 0; j < ilc->nr; j += 3)
315 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
316 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
321 for (j = 0; j < ils->nr; j += 4)
323 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
324 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
330 (w_moltype == nullptr || period2 < w_period2))
342 if (w_moltype != nullptr)
344 bWarn = (w_period2 < gmx::square(min_steps_warn*dt));
345 /* A check that would recognize most water models */
346 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
347 w_moltype->atoms.nr <= 5);
348 sprintf(warn_buf, "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated oscillational period of %.1e ps, which is less than %d times the time step of %.1e ps.\n"
351 w_a1+1, *w_moltype->atoms.atomname[w_a1],
352 w_a2+1, *w_moltype->atoms.atomname[w_a2],
353 std::sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
355 "Maybe you asked for fexible water." :
356 "Maybe you forgot to change the constraints mdp option.");
359 warning(wi, warn_buf);
363 warning_note(wi, warn_buf);
368 static void check_vel(gmx_mtop_t *mtop, rvec v[])
370 gmx_mtop_atomloop_all_t aloop;
374 aloop = gmx_mtop_atomloop_all_init(mtop);
375 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
377 if (atom->ptype == eptShell ||
378 atom->ptype == eptBond ||
379 atom->ptype == eptVSite)
386 static void check_shells_inputrec(gmx_mtop_t *mtop,
390 gmx_mtop_atomloop_all_t aloop;
393 char warn_buf[STRLEN];
395 aloop = gmx_mtop_atomloop_all_init(mtop);
396 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
398 if (atom->ptype == eptShell ||
399 atom->ptype == eptBond)
404 if ((nshells > 0) && (ir->nstcalcenergy != 1))
406 set_warning_line(wi, "unknown", -1);
407 snprintf(warn_buf, STRLEN,
408 "There are %d shells, changing nstcalcenergy from %d to 1",
409 nshells, ir->nstcalcenergy);
410 ir->nstcalcenergy = 1;
411 warning(wi, warn_buf);
415 /* TODO Decide whether this function can be consolidated with
416 * gmx_mtop_ftype_count */
417 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
422 for (mb = 0; mb < mtop->nmolblock; mb++)
424 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
430 /* This routine reorders the molecule type array
431 * in the order of use in the molblocks,
432 * unused molecule types are deleted.
434 static void renumber_moltypes(gmx_mtop_t *sys,
435 int *nmolinfo, t_molinfo **molinfo)
437 int *order, norder, i;
441 snew(order, *nmolinfo);
443 for (mb = 0; mb < sys->nmolblock; mb++)
445 for (i = 0; i < norder; i++)
447 if (order[i] == sys->molblock[mb].type)
454 /* This type did not occur yet, add it */
455 order[norder] = sys->molblock[mb].type;
456 /* Renumber the moltype in the topology */
459 sys->molblock[mb].type = i;
462 /* We still need to reorder the molinfo structs */
464 for (mi = 0; mi < *nmolinfo; mi++)
466 for (i = 0; i < norder; i++)
475 done_mi(&(*molinfo)[mi]);
479 minew[i] = (*molinfo)[mi];
488 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
493 mtop->nmoltype = nmi;
494 snew(mtop->moltype, nmi);
495 for (m = 0; m < nmi; m++)
497 molt = &mtop->moltype[m];
498 molt->name = mi[m].name;
499 molt->atoms = mi[m].atoms;
500 /* ilists are copied later */
501 molt->cgs = mi[m].cgs;
502 molt->excls = mi[m].excls;
507 new_status(const char *topfile, const char *topppfile, const char *confin,
508 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
509 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
510 gpp_atomtype_t atype, gmx_mtop_t *sys,
511 int *nmi, t_molinfo **mi, t_molinfo **intermolecular_interactions,
513 int *comb, double *reppow, real *fudgeQQ,
517 t_molinfo *molinfo = nullptr;
519 gmx_molblock_t *molblock, *molbs;
520 int mb, i, nrmols, nmismatch;
522 char warn_buf[STRLEN];
526 /* TOPOLOGY processing */
527 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
528 plist, comb, reppow, fudgeQQ,
529 atype, &nrmols, &molinfo, intermolecular_interactions,
531 &nmolblock, &molblock,
535 snew(sys->molblock, nmolblock);
538 for (mb = 0; mb < nmolblock; mb++)
540 if (sys->nmolblock > 0 &&
541 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
543 /* Merge consecutive blocks with the same molecule type */
544 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
545 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
547 else if (molblock[mb].nmol > 0)
549 /* Add a new molblock to the topology */
550 molbs = &sys->molblock[sys->nmolblock];
551 *molbs = molblock[mb];
552 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
553 molbs->nposres_xA = 0;
554 molbs->nposres_xB = 0;
555 sys->natoms += molbs->nmol*molbs->natoms_mol;
559 if (sys->nmolblock == 0)
561 gmx_fatal(FARGS, "No molecules were defined in the system");
564 renumber_moltypes(sys, &nrmols, &molinfo);
568 convert_harmonics(nrmols, molinfo, atype);
571 if (ir->eDisre == edrNone)
573 i = rm_interactions(F_DISRES, nrmols, molinfo);
576 set_warning_line(wi, "unknown", -1);
577 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
578 warning_note(wi, warn_buf);
581 if (opts->bOrire == FALSE)
583 i = rm_interactions(F_ORIRES, nrmols, molinfo);
586 set_warning_line(wi, "unknown", -1);
587 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
588 warning_note(wi, warn_buf);
592 /* Copy structures from msys to sys */
593 molinfo2mtop(nrmols, molinfo, sys);
595 gmx_mtop_finalize(sys);
597 /* COORDINATE file processing */
600 fprintf(stderr, "processing coordinates...\n");
607 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
608 state->natoms = conftop->atoms.nr;
609 if (state->natoms != sys->natoms)
611 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
612 " does not match topology (%s, %d)",
613 confin, state->natoms, topfile, sys->natoms);
615 /* It would be nice to get rid of the copies below, but we don't know
616 * a priori if the number of atoms in confin matches what we expect.
618 state->flags |= (1 << estX);
621 state->flags |= (1 << estV);
623 state_change_natoms(state, state->natoms);
624 for (int i = 0; i < state->natoms; i++)
626 copy_rvec(x[i], state->x[i]);
631 for (int i = 0; i < state->natoms; i++)
633 copy_rvec(v[i], state->v[i]);
637 /* This call fixes the box shape for runs with pressure scaling */
638 set_box_rel(ir, state);
640 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
646 sprintf(buf, "%d non-matching atom name%s\n"
647 "atom names from %s will be used\n"
648 "atom names from %s will be ignored\n",
649 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
653 /* If using the group scheme, make sure charge groups are made whole to avoid errors
654 * in calculating charge group size later on
656 if (ir->cutoff_scheme == ecutsGROUP && ir->ePBC != epbcNONE)
658 // Need temporary rvec for coordinates
659 do_pbc_first_mtop(nullptr, ir->ePBC, state->box, sys, as_rvec_array(state->x.data()));
662 /* Do more checks, mostly related to constraints */
665 fprintf(stderr, "double-checking input for internal consistency...\n");
668 int bHasNormalConstraints = 0 < (nint_ftype(sys, molinfo, F_CONSTR) +
669 nint_ftype(sys, molinfo, F_CONSTRNC));
670 int bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, molinfo, F_SETTLE);
671 double_check(ir, state->box,
672 bHasNormalConstraints,
680 gmx_mtop_atomloop_all_t aloop;
683 snew(mass, state->natoms);
684 aloop = gmx_mtop_atomloop_all_init(sys);
685 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
690 if (opts->seed == -1)
692 opts->seed = static_cast<int>(gmx::makeRandomSeed());
693 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
695 state->flags |= (1 << estV);
696 maxwell_speed(opts->tempi, opts->seed, sys, as_rvec_array(state->v.data()));
698 stop_cm(stdout, state->natoms, mass, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()));
706 static void copy_state(const char *slog, t_trxframe *fr,
707 gmx_bool bReadVel, t_state *state,
712 if (fr->not_ok & FRAME_NOT_OK)
714 gmx_fatal(FARGS, "Can not start from an incomplete frame");
718 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
722 for (i = 0; i < state->natoms; i++)
724 copy_rvec(fr->x[i], state->x[i]);
730 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
732 for (i = 0; i < state->natoms; i++)
734 copy_rvec(fr->v[i], state->v[i]);
739 copy_mat(fr->box, state->box);
742 *use_time = fr->time;
745 static void cont_status(const char *slog, const char *ener,
746 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
747 t_inputrec *ir, t_state *state,
749 const gmx_output_env_t *oenv)
750 /* If fr_time == -1 read the last frame available which is complete */
758 bReadVel = (bNeedVel && !bGenVel);
761 "Reading Coordinates%s and Box size from old trajectory\n",
762 bReadVel ? ", Velocities" : "");
765 fprintf(stderr, "Will read whole trajectory\n");
769 fprintf(stderr, "Will read till time %g\n", fr_time);
775 fprintf(stderr, "Velocities generated: "
776 "ignoring velocities in input trajectory\n");
778 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
782 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
788 "WARNING: Did not find a frame with velocities in file %s,\n"
789 " all velocities will be set to zero!\n\n", slog);
790 for (i = 0; i < sys->natoms; i++)
792 clear_rvec(state->v[i]);
795 /* Search for a frame without velocities */
797 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
801 state->natoms = fr.natoms;
803 if (sys->natoms != state->natoms)
805 gmx_fatal(FARGS, "Number of atoms in Topology "
806 "is not the same as in Trajectory");
808 copy_state(slog, &fr, bReadVel, state, &use_time);
810 /* Find the appropriate frame */
811 while ((fr_time == -1 || fr.time < fr_time) &&
812 read_next_frame(oenv, fp, &fr))
814 copy_state(slog, &fr, bReadVel, state, &use_time);
819 /* Set the relative box lengths for preserving the box shape.
820 * Note that this call can lead to differences in the last bit
821 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
823 set_box_rel(ir, state);
825 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
826 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
828 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
830 get_enx_state(ener, use_time, &sys->groups, ir, state);
831 preserve_box_shape(ir, state->box_rel, state->boxv);
835 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
837 int rc_scaling, int ePBC,
847 int natoms, npbcdim = 0;
848 char warn_buf[STRLEN];
849 int a, i, ai, j, k, mb, nat_molb;
850 gmx_molblock_t *molb;
855 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
856 natoms = top->atoms.nr;
859 if (natoms != mtop->natoms)
861 sprintf(warn_buf, "The number of atoms in %s (%d) does not match the number of atoms in the topology (%d). Will assume that the first %d atoms in the topology and %s match.", fn, natoms, mtop->natoms, std::min(mtop->natoms, natoms), fn);
862 warning(wi, warn_buf);
865 npbcdim = ePBC2npbcdim(ePBC);
867 if (rc_scaling != erscNO)
869 copy_mat(box, invbox);
870 for (j = npbcdim; j < DIM; j++)
872 clear_rvec(invbox[j]);
875 gmx::invertBoxMatrix(invbox, invbox);
878 /* Copy the reference coordinates to mtop */
882 snew(hadAtom, natoms);
883 for (mb = 0; mb < mtop->nmolblock; mb++)
885 molb = &mtop->molblock[mb];
886 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
887 pr = &(molinfo[molb->type].plist[F_POSRES]);
888 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
889 if (pr->nr > 0 || prfb->nr > 0)
891 atom = mtop->moltype[molb->type].atoms.atom;
892 for (i = 0; (i < pr->nr); i++)
894 ai = pr->param[i].ai();
897 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
898 ai+1, *molinfo[molb->type].name, fn, natoms);
901 if (rc_scaling == erscCOM)
903 /* Determine the center of mass of the posres reference coordinates */
904 for (j = 0; j < npbcdim; j++)
906 sum[j] += atom[ai].m*x[a+ai][j];
908 totmass += atom[ai].m;
911 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
912 for (i = 0; (i < prfb->nr); i++)
914 ai = prfb->param[i].ai();
917 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
918 ai+1, *molinfo[molb->type].name, fn, natoms);
920 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
922 /* Determine the center of mass of the posres reference coordinates */
923 for (j = 0; j < npbcdim; j++)
925 sum[j] += atom[ai].m*x[a+ai][j];
927 totmass += atom[ai].m;
932 molb->nposres_xA = nat_molb;
933 snew(molb->posres_xA, molb->nposres_xA);
934 for (i = 0; i < nat_molb; i++)
936 copy_rvec(x[a+i], molb->posres_xA[i]);
941 molb->nposres_xB = nat_molb;
942 snew(molb->posres_xB, molb->nposres_xB);
943 for (i = 0; i < nat_molb; i++)
945 copy_rvec(x[a+i], molb->posres_xB[i]);
951 if (rc_scaling == erscCOM)
955 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
957 for (j = 0; j < npbcdim; j++)
959 com[j] = sum[j]/totmass;
961 fprintf(stderr, "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n", com[XX], com[YY], com[ZZ]);
964 if (rc_scaling != erscNO)
966 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
968 for (mb = 0; mb < mtop->nmolblock; mb++)
970 molb = &mtop->molblock[mb];
971 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
972 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
974 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
975 for (i = 0; i < nat_molb; i++)
977 for (j = 0; j < npbcdim; j++)
979 if (rc_scaling == erscALL)
981 /* Convert from Cartesian to crystal coordinates */
982 xp[i][j] *= invbox[j][j];
983 for (k = j+1; k < npbcdim; k++)
985 xp[i][j] += invbox[k][j]*xp[i][k];
988 else if (rc_scaling == erscCOM)
990 /* Subtract the center of mass */
998 if (rc_scaling == erscCOM)
1000 /* Convert the COM from Cartesian to crystal coordinates */
1001 for (j = 0; j < npbcdim; j++)
1003 com[j] *= invbox[j][j];
1004 for (k = j+1; k < npbcdim; k++)
1006 com[j] += invbox[k][j]*com[k];
1017 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
1018 const char *fnA, const char *fnB,
1019 int rc_scaling, int ePBC,
1020 rvec com, rvec comB,
1023 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1024 /* It is safer to simply read the b-state posres rather than trying
1025 * to be smart and copy the positions.
1027 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1030 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
1031 t_inputrec *ir, warninp_t wi)
1034 char warn_buf[STRLEN];
1038 fprintf(stderr, "Searching the wall atom type(s)\n");
1040 for (i = 0; i < ir->nwall; i++)
1042 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1043 if (ir->wall_atomtype[i] == NOTSET)
1045 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1046 warning_error(wi, warn_buf);
1051 static int nrdf_internal(t_atoms *atoms)
1056 for (i = 0; i < atoms->nr; i++)
1058 /* Vsite ptype might not be set here yet, so also check the mass */
1059 if ((atoms->atom[i].ptype == eptAtom ||
1060 atoms->atom[i].ptype == eptNucleus)
1061 && atoms->atom[i].m > 0)
1068 case 0: nrdf = 0; break;
1069 case 1: nrdf = 0; break;
1070 case 2: nrdf = 1; break;
1071 default: nrdf = nmass*3 - 6; break;
1078 spline1d( double dx,
1090 for (i = 1; i < n-1; i++)
1092 p = 0.5*y2[i-1]+2.0;
1094 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1095 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1100 for (i = n-2; i >= 0; i--)
1102 y2[i] = y2[i]*y2[i+1]+u[i];
1108 interpolate1d( double xmin,
1119 ix = static_cast<int>((x-xmin)/dx);
1121 a = (xmin+(ix+1)*dx-x)/dx;
1122 b = (x-xmin-ix*dx)/dx;
1124 *y = a*ya[ix]+b*ya[ix+1]+((a*a*a-a)*y2a[ix]+(b*b*b-b)*y2a[ix+1])*(dx*dx)/6.0;
1125 *y1 = (ya[ix+1]-ya[ix])/dx-(3.0*a*a-1.0)/6.0*dx*y2a[ix]+(3.0*b*b-1.0)/6.0*dx*y2a[ix+1];
1130 setup_cmap (int grid_spacing,
1133 gmx_cmap_t * cmap_grid)
1135 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1137 int i, j, k, ii, jj, kk, idx;
1139 double dx, xmin, v, v1, v2, v12;
1142 snew(tmp_u, 2*grid_spacing);
1143 snew(tmp_u2, 2*grid_spacing);
1144 snew(tmp_yy, 2*grid_spacing);
1145 snew(tmp_y1, 2*grid_spacing);
1146 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1147 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1149 dx = 360.0/grid_spacing;
1150 xmin = -180.0-dx*grid_spacing/2;
1152 for (kk = 0; kk < nc; kk++)
1154 /* Compute an offset depending on which cmap we are using
1155 * Offset will be the map number multiplied with the
1156 * grid_spacing * grid_spacing * 2
1158 offset = kk * grid_spacing * grid_spacing * 2;
1160 for (i = 0; i < 2*grid_spacing; i++)
1162 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1164 for (j = 0; j < 2*grid_spacing; j++)
1166 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1167 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1171 for (i = 0; i < 2*grid_spacing; i++)
1173 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1176 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1178 ii = i-grid_spacing/2;
1181 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1183 jj = j-grid_spacing/2;
1186 for (k = 0; k < 2*grid_spacing; k++)
1188 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1189 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1192 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1193 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1194 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1195 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1197 idx = ii*grid_spacing+jj;
1198 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1199 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1200 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1201 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1207 static void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1211 cmap_grid->ngrid = ngrid;
1212 cmap_grid->grid_spacing = grid_spacing;
1213 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1215 snew(cmap_grid->cmapdata, ngrid);
1217 for (i = 0; i < cmap_grid->ngrid; i++)
1219 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1224 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1226 int count, count_mol, i, mb;
1227 gmx_molblock_t *molb;
1232 for (mb = 0; mb < mtop->nmolblock; mb++)
1235 molb = &mtop->molblock[mb];
1236 plist = mi[molb->type].plist;
1238 for (i = 0; i < F_NRE; i++)
1242 count_mol += 3*plist[i].nr;
1244 else if (interaction_function[i].flags & IF_CONSTRAINT)
1246 count_mol += plist[i].nr;
1250 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1253 "Molecule type '%s' has %d constraints.\n"
1254 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1255 *mi[molb->type].name, count_mol,
1256 nrdf_internal(&mi[molb->type].atoms));
1259 count += molb->nmol*count_mol;
1265 static real calc_temp(const gmx_mtop_t *mtop,
1266 const t_inputrec *ir,
1269 gmx_mtop_atomloop_all_t aloop;
1274 aloop = gmx_mtop_atomloop_all_init(mtop);
1275 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1277 sum_mv2 += atom->m*norm2(v[a]);
1281 for (int g = 0; g < ir->opts.ngtc; g++)
1283 nrdf += ir->opts.nrdf[g];
1286 return sum_mv2/(nrdf*BOLTZ);
1289 static real get_max_reference_temp(const t_inputrec *ir,
1298 for (i = 0; i < ir->opts.ngtc; i++)
1300 if (ir->opts.tau_t[i] < 0)
1306 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1314 sprintf(buf, "Some temperature coupling groups do not use temperature coupling. We will assume their temperature is not more than %.3f K. If their temperature is higher, the energy error and the Verlet buffer might be underestimated.",
1322 /* Checks if there are unbound atoms in moleculetype molt.
1323 * Prints a note for each unbound atoms and a warning if any is present.
1325 static void checkForUnboundAtoms(const gmx_moltype_t *molt,
1329 const t_atoms *atoms = &molt->atoms;
1333 /* Only one atom, there can't be unbound atoms */
1337 std::vector<int> count(atoms->nr, 0);
1339 for (int ftype = 0; ftype < F_NRE; ftype++)
1341 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS) ||
1342 (interaction_function[ftype].flags & IF_CONSTRAINT) ||
1345 const t_ilist *il = &molt->ilist[ftype];
1346 int nral = NRAL(ftype);
1348 for (int i = 0; i < il->nr; i += 1 + nral)
1350 for (int j = 0; j < nral; j++)
1352 count[il->iatoms[i + 1 + j]]++;
1358 int numDanglingAtoms = 0;
1359 for (int a = 0; a < atoms->nr; a++)
1361 if (atoms->atom[a].ptype != eptVSite &&
1366 fprintf(stderr, "\nAtom %d '%s' in moleculetype '%s' is not bound by a potential or constraint to any other atom in the same moleculetype.\n",
1367 a + 1, *atoms->atomname[a], *molt->name);
1373 if (numDanglingAtoms > 0)
1376 sprintf(buf, "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any other atom in the same moleculetype. Although technically this might not cause issues in a simulation, this often means that the user forgot to add a bond/potential/constraint or put multiple molecules in the same moleculetype definition by mistake. Run with -v to get information for each atom.",
1377 *molt->name, numDanglingAtoms);
1378 warning_note(wi, buf);
1382 /* Checks all moleculetypes for unbound atoms */
1383 static void checkForUnboundAtoms(const gmx_mtop_t *mtop,
1387 for (int mt = 0; mt < mtop->nmoltype; mt++)
1389 checkForUnboundAtoms(&mtop->moltype[mt], bVerbose, wi);
1393 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1395 * The specific decoupled modes this routine check for are angle modes
1396 * where the two bonds are constrained and the atoms a both ends are only
1397 * involved in a single constraint; the mass of the two atoms needs to
1398 * differ by more than \p massFactorThreshold.
1400 static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
1401 const t_iparams *iparams,
1402 real massFactorThreshold)
1404 if (molt->ilist[F_CONSTR].nr == 0 && molt->ilist[F_CONSTRNC].nr == 0)
1409 const t_atom * atom = molt->atoms.atom;
1411 int numFlexibleConstraints;
1412 t_blocka atomToConstraints = make_at2con(0, molt->atoms.nr,
1413 molt->ilist, iparams,
1415 &numFlexibleConstraints);
1417 bool haveDecoupledMode = false;
1418 for (int ftype = 0; ftype < F_NRE; ftype++)
1420 if (interaction_function[ftype].flags & IF_ATYPE)
1422 const int nral = NRAL(ftype);
1423 const t_ilist *il = &molt->ilist[ftype];
1424 for (int i = 0; i < il->nr; i += 1 + nral)
1426 /* Here we check for the mass difference between the atoms
1427 * at both ends of the angle, that the atoms at the ends have
1428 * 1 contraint and the atom in the middle at least 3; we check
1429 * that the 3 atoms are linked by constraints below.
1430 * We check for at least three constraints for the middle atom,
1431 * since with only the two bonds in the angle, we have 3-atom
1432 * molecule, which has much more thermal exhange in this single
1433 * angle mode than molecules with more atoms.
1434 * Note that this check also catches molecules where more atoms
1435 * are connected to one or more atoms in the angle, but by
1436 * bond potentials instead of angles. But such cases will not
1437 * occur in "normal" molecules and it doesn't hurt running
1438 * those with higher accuracy settings as well.
1440 int a0 = il->iatoms[1 + i];
1441 int a1 = il->iatoms[1 + i + 1];
1442 int a2 = il->iatoms[1 + i + 2];
1443 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
1444 atom[a2].m > atom[a0].m*massFactorThreshold) &&
1445 atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
1446 atomToConstraints.index[a2 + 1] - atomToConstraints.index[a2] == 1 &&
1447 atomToConstraints.index[a1 + 1] - atomToConstraints.index[a1] >= 3)
1449 int constraint0 = atomToConstraints.a[atomToConstraints.index[a0]];
1450 int constraint2 = atomToConstraints.a[atomToConstraints.index[a2]];
1452 bool foundAtom0 = false;
1453 bool foundAtom2 = false;
1454 for (int conIndex = atomToConstraints.index[a1]; conIndex < atomToConstraints.index[a1 + 1]; conIndex++)
1456 if (atomToConstraints.a[conIndex] == constraint0)
1460 if (atomToConstraints.a[conIndex] == constraint2)
1465 if (foundAtom0 && foundAtom2)
1467 haveDecoupledMode = true;
1474 done_blocka(&atomToConstraints);
1476 return haveDecoupledMode;
1479 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1481 * When decoupled modes are present and the accuray in insufficient,
1482 * this routine issues a warning if the accuracy is insufficient.
1484 static void checkDecoupledModeAccuracy(const gmx_mtop_t *mtop,
1485 const t_inputrec *ir,
1488 /* We only have issues with decoupled modes with normal MD.
1489 * With stochastic dynamics equipartitioning is enforced strongly.
1496 /* When atoms of very different mass are involved in an angle potential
1497 * and both bonds in the angle are constrained, the dynamic modes in such
1498 * angles have very different periods and significant energy exchange
1499 * takes several nanoseconds. Thus even a small amount of error in
1500 * different algorithms can lead to violation of equipartitioning.
1501 * The parameters below are mainly based on an all-atom chloroform model
1502 * with all bonds constrained. Equipartitioning is off by more than 1%
1503 * (need to run 5-10 ns) when the difference in mass between H and Cl
1504 * is more than a factor 13 and the accuracy is less than the thresholds
1505 * given below. This has been verified on some other molecules.
1507 * Note that the buffer and shake parameters have unit length and
1508 * energy/time, respectively, so they will "only" work correctly
1509 * for atomistic force fields using MD units.
1511 const real massFactorThreshold = 13.0;
1512 const real bufferToleranceThreshold = 1e-4;
1513 const int lincsIterationThreshold = 2;
1514 const int lincsOrderThreshold = 4;
1515 const real shakeToleranceThreshold = 0.005*ir->delta_t;
1517 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1518 bool shakeWithSufficientTolerance = (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1*shakeToleranceThreshold);
1519 if (ir->cutoff_scheme == ecutsVERLET &&
1520 ir->verletbuf_tol <= 1.1*bufferToleranceThreshold &&
1521 (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1526 bool haveDecoupledMode = false;
1527 for (int mt = 0; mt < mtop->nmoltype; mt++)
1529 if (haveDecoupledModeInMol(&mtop->moltype[mt], mtop->ffparams.iparams,
1530 massFactorThreshold))
1532 haveDecoupledMode = true;
1536 if (haveDecoupledMode)
1538 char modeMessage[STRLEN];
1539 sprintf(modeMessage, "There are atoms at both ends of an angle, connected by constraints and with masses that differ by more than a factor of %g. This means that there are likely dynamic modes that are only very weakly coupled.",
1540 massFactorThreshold);
1542 if (ir->cutoff_scheme == ecutsVERLET)
1544 sprintf(buf, "%s To ensure good equipartitioning, you need to either not use constraints on all bonds (but, if possible, only on bonds involving hydrogens) or use integrator = %s or decrease one or more tolerances: verlet-buffer-tolerance <= %g, LINCS iterations >= %d, LINCS order >= %d or SHAKE tolerance <= %g",
1547 bufferToleranceThreshold,
1548 lincsIterationThreshold, lincsOrderThreshold,
1549 shakeToleranceThreshold);
1553 sprintf(buf, "%s To ensure good equipartitioning, we suggest to switch to the %s cutoff-scheme, since that allows for better control over the Verlet buffer size and thus over the energy drift.",
1555 ecutscheme_names[ecutsVERLET]);
1561 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1569 char warn_buf[STRLEN];
1571 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1573 /* Calculate the buffer size for simple atom vs atoms list */
1574 VerletbufListSetup listSetup1x1;
1575 listSetup1x1.cluster_size_i = 1;
1576 listSetup1x1.cluster_size_j = 1;
1577 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1,
1578 buffer_temp, &listSetup1x1,
1579 &n_nonlin_vsite, &rlist_1x1);
1581 /* Set the pair-list buffer size in ir */
1582 VerletbufListSetup listSetup4x4 =
1583 verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1584 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1,
1585 buffer_temp, &listSetup4x4,
1586 &n_nonlin_vsite, &ir->rlist);
1588 if (n_nonlin_vsite > 0)
1590 sprintf(warn_buf, "There are %d non-linear virtual site constructions. Their contribution to the energy error is approximated. In most cases this does not affect the error significantly.", n_nonlin_vsite);
1591 warning_note(wi, warn_buf);
1594 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1595 1, 1, rlist_1x1, rlist_1x1-std::max(ir->rvdw, ir->rcoulomb));
1597 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1598 listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j,
1599 ir->rlist, ir->rlist-std::max(ir->rvdw, ir->rcoulomb));
1601 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1603 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1605 gmx_fatal(FARGS, "The pair-list cut-off (%g nm) is longer than half the shortest box vector or longer than the smallest box diagonal element (%g nm). Increase the box size or decrease nstlist or increase verlet-buffer-tolerance.", ir->rlist, std::sqrt(max_cutoff2(ir->ePBC, box)));
1609 int gmx_grompp(int argc, char *argv[])
1611 const char *desc[] = {
1612 "[THISMODULE] (the gromacs preprocessor)",
1613 "reads a molecular topology file, checks the validity of the",
1614 "file, expands the topology from a molecular description to an atomic",
1615 "description. The topology file contains information about",
1616 "molecule types and the number of molecules, the preprocessor",
1617 "copies each molecule as needed. ",
1618 "There is no limitation on the number of molecule types. ",
1619 "Bonds and bond-angles can be converted into constraints, separately",
1620 "for hydrogens and heavy atoms.",
1621 "Then a coordinate file is read and velocities can be generated",
1622 "from a Maxwellian distribution if requested.",
1623 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1624 "(eg. number of MD steps, time step, cut-off), and others such as",
1625 "NEMD parameters, which are corrected so that the net acceleration",
1627 "Eventually a binary file is produced that can serve as the sole input",
1628 "file for the MD program.[PAR]",
1630 "[THISMODULE] uses the atom names from the topology file. The atom names",
1631 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1632 "warnings when they do not match the atom names in the topology.",
1633 "Note that the atom names are irrelevant for the simulation as",
1634 "only the atom types are used for generating interaction parameters.[PAR]",
1636 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1637 "etc. The preprocessor supports the following keywords::",
1640 " #ifndef VARIABLE",
1643 " #define VARIABLE",
1645 " #include \"filename\"",
1646 " #include <filename>",
1648 "The functioning of these statements in your topology may be modulated by",
1649 "using the following two flags in your [REF].mdp[ref] file::",
1651 " define = -DVARIABLE1 -DVARIABLE2",
1652 " include = -I/home/john/doe",
1654 "For further information a C-programming textbook may help you out.",
1655 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1656 "topology file written out so that you can verify its contents.[PAR]",
1658 "When using position restraints, a file with restraint coordinates",
1659 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1660 "for [TT]-c[tt]). For free energy calculations, separate reference",
1661 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1662 "otherwise they will be equal to those of the A topology.[PAR]",
1664 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1665 "The last frame with coordinates and velocities will be read,",
1666 "unless the [TT]-time[tt] option is used. Only if this information",
1667 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1668 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1669 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1670 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1673 "[THISMODULE] can be used to restart simulations (preserving",
1674 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1675 "However, for simply changing the number of run steps to extend",
1676 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1677 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1678 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1679 "like output frequency, then supplying the checkpoint file to",
1680 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1681 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1682 "the ensemble (if possible) still requires passing the checkpoint",
1683 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1685 "By default, all bonded interactions which have constant energy due to",
1686 "virtual site constructions will be removed. If this constant energy is",
1687 "not zero, this will result in a shift in the total energy. All bonded",
1688 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1689 "all constraints for distances which will be constant anyway because",
1690 "of virtual site constructions will be removed. If any constraints remain",
1691 "which involve virtual sites, a fatal error will result.[PAR]"
1693 "To verify your run input file, please take note of all warnings",
1694 "on the screen, and correct where necessary. Do also look at the contents",
1695 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1696 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1697 "with the [TT]-debug[tt] option which will give you more information",
1698 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1699 "can see the contents of the run input file with the [gmx-dump]",
1700 "program. [gmx-check] can be used to compare the contents of two",
1701 "run input files.[PAR]"
1703 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1704 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1705 "harmless, but usually they are not. The user is advised to carefully",
1706 "interpret the output messages before attempting to bypass them with",
1712 t_molinfo *mi, *intermolecular_interactions;
1713 gpp_atomtype_t atype;
1714 int nvsite, comb, mt;
1718 const char *mdparin;
1720 gmx_bool bNeedVel, bGenVel;
1721 gmx_bool have_atomnumber;
1722 gmx_output_env_t *oenv;
1723 gmx_bool bVerbose = FALSE;
1725 char warn_buf[STRLEN];
1728 { efMDP, nullptr, nullptr, ffREAD },
1729 { efMDP, "-po", "mdout", ffWRITE },
1730 { efSTX, "-c", nullptr, ffREAD },
1731 { efSTX, "-r", "restraint", ffOPTRD },
1732 { efSTX, "-rb", "restraint", ffOPTRD },
1733 { efNDX, nullptr, nullptr, ffOPTRD },
1734 { efTOP, nullptr, nullptr, ffREAD },
1735 { efTOP, "-pp", "processed", ffOPTWR },
1736 { efTPR, "-o", nullptr, ffWRITE },
1737 { efTRN, "-t", nullptr, ffOPTRD },
1738 { efEDR, "-e", nullptr, ffOPTRD },
1739 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1740 { efGRO, "-imd", "imdgroup", ffOPTWR },
1741 { efTRN, "-ref", "rotref", ffOPTRW }
1743 #define NFILE asize(fnm)
1745 /* Command line options */
1746 gmx_bool bRenum = TRUE;
1747 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1751 { "-v", FALSE, etBOOL, {&bVerbose},
1752 "Be loud and noisy" },
1753 { "-time", FALSE, etREAL, {&fr_time},
1754 "Take frame at or first after this time." },
1755 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1756 "Remove constant bonded interactions with virtual sites" },
1757 { "-maxwarn", FALSE, etINT, {&maxwarn},
1758 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1759 { "-zero", FALSE, etBOOL, {&bZero},
1760 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1761 { "-renum", FALSE, etBOOL, {&bRenum},
1762 "Renumber atomtypes and minimize number of atomtypes" }
1765 /* Parse the command line */
1766 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1767 asize(desc), desc, 0, nullptr, &oenv))
1772 /* Initiate some variables */
1773 gmx::MDModules mdModules;
1774 t_inputrec irInstance;
1775 t_inputrec *ir = &irInstance;
1777 snew(opts->include, STRLEN);
1778 snew(opts->define, STRLEN);
1780 wi = init_warning(TRUE, maxwarn);
1782 /* PARAMETER file processing */
1783 mdparin = opt2fn("-f", NFILE, fnm);
1784 set_warning_line(wi, mdparin, -1);
1787 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1789 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1793 fprintf(stderr, "checking input for internal consistency...\n");
1795 check_ir(mdparin, ir, opts, wi);
1797 if (ir->ld_seed == -1)
1799 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1800 fprintf(stderr, "Setting the LD random seed to %" GMX_PRId64 "\n", ir->ld_seed);
1803 if (ir->expandedvals->lmc_seed == -1)
1805 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1806 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1809 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1810 bGenVel = (bNeedVel && opts->bGenVel);
1811 if (bGenVel && ir->bContinuation)
1814 "Generating velocities is inconsistent with attempting "
1815 "to continue a previous run. Choose only one of "
1816 "gen-vel = yes and continuation = yes.");
1817 warning_error(wi, warn_buf);
1823 atype = init_atomtype();
1826 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1829 const char *fn = ftp2fn(efTOP, NFILE, fnm);
1830 if (!gmx_fexist(fn))
1832 gmx_fatal(FARGS, "%s does not exist", fn);
1836 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1837 opts, ir, bZero, bGenVel, bVerbose, &state,
1838 atype, sys, &nmi, &mi, &intermolecular_interactions,
1839 plist, &comb, &reppow, &fudgeQQ,
1845 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1849 /* set parameters for virtual site construction (not for vsiten) */
1850 for (mt = 0; mt < sys->nmoltype; mt++)
1853 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1855 /* now throw away all obsolete bonds, angles and dihedrals: */
1856 /* note: constraints are ALWAYS removed */
1859 for (mt = 0; mt < sys->nmoltype; mt++)
1861 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1865 if (ir->cutoff_scheme == ecutsVERLET)
1867 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1868 ecutscheme_names[ir->cutoff_scheme]);
1870 /* Remove all charge groups */
1871 gmx_mtop_remove_chargegroups(sys);
1874 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1876 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1878 sprintf(warn_buf, "Can not do %s with %s, use %s",
1879 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1880 warning_error(wi, warn_buf);
1882 if (ir->bPeriodicMols)
1884 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1885 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1886 warning_error(wi, warn_buf);
1890 if (EI_SD (ir->eI) && ir->etc != etcNO)
1892 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1895 /* If we are doing QM/MM, check that we got the atom numbers */
1896 have_atomnumber = TRUE;
1897 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1899 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1901 if (!have_atomnumber && ir->bQMMM)
1905 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1906 "field you are using does not contain atom numbers fields. This is an\n"
1907 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1908 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1909 "an integer just before the mass column in ffXXXnb.itp.\n"
1910 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1913 /* Check for errors in the input now, since they might cause problems
1914 * during processing further down.
1916 check_warning_error(wi, FARGS);
1918 if (nint_ftype(sys, mi, F_POSRES) > 0 ||
1919 nint_ftype(sys, mi, F_FBPOSRES) > 0)
1921 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1923 sprintf(warn_buf, "You are combining position restraints with %s pressure coupling, which can lead to instabilities. If you really want to combine position restraints with pressure coupling, we suggest to use %s pressure coupling instead.",
1924 EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
1925 warning_note(wi, warn_buf);
1928 const char *fn = opt2fn("-r", NFILE, fnm);
1931 if (!gmx_fexist(fn))
1934 "Cannot find position restraint file %s (option -r).\n"
1935 "From GROMACS-2018, you need to specify the position restraint "
1936 "coordinate files explicitly to avoid mistakes, although you can "
1937 "still use the same file as you specify for the -c option.", fn);
1940 if (opt2bSet("-rb", NFILE, fnm))
1942 fnB = opt2fn("-rb", NFILE, fnm);
1943 if (!gmx_fexist(fnB))
1946 "Cannot find B-state position restraint file %s (option -rb).\n"
1947 "From GROMACS-2018, you need to specify the position restraint "
1948 "coordinate files explicitly to avoid mistakes, although you can "
1949 "still use the same file as you specify for the -c option.", fn);
1959 fprintf(stderr, "Reading position restraint coords from %s", fn);
1960 if (strcmp(fn, fnB) == 0)
1962 fprintf(stderr, "\n");
1966 fprintf(stderr, " and %s\n", fnB);
1969 gen_posres(sys, mi, fn, fnB,
1970 ir->refcoord_scaling, ir->ePBC,
1971 ir->posres_com, ir->posres_comB,
1975 /* If we are using CMAP, setup the pre-interpolation grid */
1976 if (plist[F_CMAP].ncmap > 0)
1978 init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
1979 setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
1982 set_wall_atomtype(atype, opts, ir, wi);
1985 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1986 get_atomtype_ntypes(atype);
1989 if (ir->implicit_solvent)
1991 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
1994 /* PELA: Copy the atomtype data to the topology atomtype list */
1995 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1999 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
2004 fprintf(stderr, "converting bonded parameters...\n");
2007 ntype = get_atomtype_ntypes(atype);
2008 convert_params(ntype, plist, mi, intermolecular_interactions,
2009 comb, reppow, fudgeQQ, sys);
2013 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
2016 /* set ptype to VSite for virtual sites */
2017 for (mt = 0; mt < sys->nmoltype; mt++)
2019 set_vsites_ptype(FALSE, &sys->moltype[mt]);
2023 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
2025 /* Check velocity for virtual sites and shells */
2028 check_vel(sys, as_rvec_array(state.v.data()));
2031 /* check for shells and inpurecs */
2032 check_shells_inputrec(sys, ir, wi);
2037 checkForUnboundAtoms(sys, bVerbose, wi);
2039 for (i = 0; i < sys->nmoltype; i++)
2041 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
2044 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2046 check_bonds_timestep(sys, ir->delta_t, wi);
2049 checkDecoupledModeAccuracy(sys, ir, wi);
2051 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2053 warning_note(wi, "Zero-step energy minimization will alter the coordinates before calculating the energy. If you just want the energy of a single point, try zero-step MD (with unconstrained_start = yes). To do multiple single-point energy evaluations of different configurations of the same topology, use mdrun -rerun.");
2056 check_warning_error(wi, FARGS);
2060 fprintf(stderr, "initialising group options...\n");
2062 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
2066 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2068 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2072 if (EI_MD(ir->eI) && ir->etc == etcNO)
2076 buffer_temp = opts->tempi;
2080 buffer_temp = calc_temp(sys, ir, as_rvec_array(state.v.data()));
2082 if (buffer_temp > 0)
2084 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
2085 warning_note(wi, warn_buf);
2089 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
2090 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
2091 warning_note(wi, warn_buf);
2096 buffer_temp = get_max_reference_temp(ir, wi);
2099 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2101 /* NVE with initial T=0: we add a fixed ratio to rlist.
2102 * Since we don't actually use verletbuf_tol, we set it to -1
2103 * so it can't be misused later.
2105 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2106 ir->verletbuf_tol = -1;
2110 /* We warn for NVE simulations with a drift tolerance that
2111 * might result in a 1(.1)% drift over the total run-time.
2112 * Note that we can't warn when nsteps=0, since we don't
2113 * know how many steps the user intends to run.
2115 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
2118 const real driftTolerance = 0.01;
2119 /* We use 2 DOF per atom = 2kT pot+kin energy,
2120 * to be on the safe side with constraints.
2122 const real totalEnergyDriftPerAtomPerPicosecond = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
2124 if (ir->verletbuf_tol > 1.1*driftTolerance*totalEnergyDriftPerAtomPerPicosecond)
2126 sprintf(warn_buf, "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an NVE simulation of length %g ps, which can give a final drift of %d%%. For conserving energy to %d%% when using constraints, you might need to set verlet-buffer-tolerance to %.1e.",
2127 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
2128 (int)(ir->verletbuf_tol/totalEnergyDriftPerAtomPerPicosecond*100 + 0.5),
2129 (int)(100*driftTolerance + 0.5),
2130 driftTolerance*totalEnergyDriftPerAtomPerPicosecond);
2131 warning_note(wi, warn_buf);
2135 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
2140 /* Init the temperature coupling state */
2141 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2145 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
2147 check_eg_vs_cg(sys);
2151 pr_symtab(debug, 0, "After index", &sys->symtab);
2154 triple_check(mdparin, ir, sys, wi);
2155 close_symtab(&sys->symtab);
2158 pr_symtab(debug, 0, "After close", &sys->symtab);
2161 /* make exclusions between QM atoms */
2164 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
2166 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
2170 generate_qmexcl(sys, ir, wi);
2174 if (ftp2bSet(efTRN, NFILE, fnm))
2178 fprintf(stderr, "getting data from old trajectory ...\n");
2180 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
2181 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
2184 if (ir->ePBC == epbcXY && ir->nwall != 2)
2186 clear_rvec(state.box[ZZ]);
2189 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
2191 set_warning_line(wi, mdparin, -1);
2192 check_chargegroup_radii(sys, ir, as_rvec_array(state.x.data()), wi);
2195 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2197 /* Calculate the optimal grid dimensions */
2199 EwaldBoxZScaler boxScaler(*ir);
2200 boxScaler.scaleBox(state.box, scaledBox);
2202 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2204 /* Mark fourier_spacing as not used */
2205 ir->fourier_spacing = 0;
2207 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2209 set_warning_line(wi, mdparin, -1);
2210 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2212 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2213 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize,
2214 &(ir->nkx), &(ir->nky), &(ir->nkz));
2215 if (ir->nkx < minGridSize ||
2216 ir->nky < minGridSize ||
2217 ir->nkz < minGridSize)
2219 warning_error(wi, "The PME grid size should be >= 2*(pme-order - 1); either manually increase the grid size or decrease pme-order");
2223 /* MRS: eventually figure out better logic for initializing the fep
2224 values that makes declaring the lambda and declaring the state not
2225 potentially conflict if not handled correctly. */
2226 if (ir->efep != efepNO)
2228 state.fep_state = ir->fepvals->init_fep_state;
2229 for (i = 0; i < efptNR; i++)
2231 /* init_lambda trumps state definitions*/
2232 if (ir->fepvals->init_lambda >= 0)
2234 state.lambda[i] = ir->fepvals->init_lambda;
2238 if (ir->fepvals->all_lambda[i] == nullptr)
2240 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2244 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2250 struct pull_t *pull = nullptr;
2254 pull = set_pull_init(ir, sys, as_rvec_array(state.x.data()), state.box, state.lambda[efptMASS], oenv);
2257 /* Modules that supply external potential for pull coordinates
2258 * should register those potentials here. finish_pull() will check
2259 * that providers have been registerd for all external potentials.
2264 setStateDependentAwhParams(ir->awhParams, ir->pull, pull,
2265 state.box, ir->ePBC, &ir->opts, wi);
2275 set_reference_positions(ir->rot, as_rvec_array(state.x.data()), state.box,
2276 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2280 /* reset_multinr(sys); */
2282 if (EEL_PME(ir->coulombtype))
2284 float ratio = pme_load_estimate(sys, ir, state.box);
2285 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2286 /* With free energy we might need to do PME both for the A and B state
2287 * charges. This will double the cost, but the optimal performance will
2288 * then probably be at a slightly larger cut-off and grid spacing.
2290 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2291 (ir->efep != efepNO && ratio > 2.0/3.0))
2294 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2295 "and for highly parallel simulations between 0.25 and 0.33,\n"
2296 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2297 if (ir->efep != efepNO)
2300 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2306 char warn_buf[STRLEN];
2307 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2308 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2311 set_warning_line(wi, mdparin, -1);
2312 warning_note(wi, warn_buf);
2316 printf("%s\n", warn_buf);
2322 fprintf(stderr, "writing run input file...\n");
2325 done_warning(wi, FARGS);
2326 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, sys);
2328 /* Output IMD group, if bIMD is TRUE */
2329 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2331 done_atomtype(atype);
2333 done_inputrec_strings();