Free more memory in grompp and mdrun
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.cpp
blob394bf80086ab539049d14c4da27261f5c3c48e28
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,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.
37 #include "gmxpre.h"
39 #include "grompp.h"
41 #include <errno.h>
42 #include <limits.h>
43 #include <string.h>
45 #include <cmath>
47 #include <algorithm>
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[])
107 int i, n;
109 n = 0;
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;
116 return n;
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;
123 t_atoms *tat;
124 #define MAXMISMATCH 20
126 if (mtop->natoms != at->nr)
128 gmx_incons("comparing atom names");
131 nmismatch = 0;
132 i = 0;
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)
144 fprintf(stderr,
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);
152 nmismatch++;
154 i++;
159 return nmismatch;
162 static void check_eg_vs_cg(gmx_mtop_t *mtop)
164 int astart, mb, m, cg, j, firstj;
165 unsigned char firsteg, eg;
166 gmx_moltype_t *molt;
168 /* Go through all the charge groups and make sure all their
169 * atoms are in the same energy group.
172 astart = 0;
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);
186 if (eg != firsteg)
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)
200 int maxsize, cg;
201 char warn_buf[STRLEN];
203 maxsize = 0;
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);
216 sprintf(warn_buf,
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.",
221 maxsize);
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;
248 t_iparams *ip;
249 int molt;
250 gmx_moltype_t *moltype, *w_moltype;
251 t_atom *atom;
252 t_ilist *ilist, *ilb, *ilc, *ils;
253 int ftype;
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);
265 w_a1 = w_a2 = -1;
266 w_period2 = -1.0;
268 w_moltype = nullptr;
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))
280 continue;
283 ilb = &ilist[ftype];
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 */
291 fc = 2*fc*re;
293 a1 = ilb->iatoms[i+1];
294 a2 = ilb->iatoms[i+2];
295 m1 = atom[a1].m;
296 m2 = atom[a2].m;
297 if (fc > 0 && m1 > 0 && m2 > 0)
299 period2 = twopi2*m1*m2/((m1 + m2)*fc);
301 else
303 period2 = GMX_FLOAT_MAX;
305 if (debug)
307 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
308 fc, m1, m2, std::sqrt(period2));
310 if (period2 < limit2)
312 bFound = FALSE;
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))
318 bFound = TRUE;
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]))
326 bFound = TRUE;
329 if (!bFound &&
330 (w_moltype == nullptr || period2 < w_period2))
332 w_moltype = moltype;
333 w_a1 = a1;
334 w_a2 = a2;
335 w_period2 = 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"
349 "%s",
350 *w_moltype->name,
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,
354 bWater ?
355 "Maybe you asked for fexible water." :
356 "Maybe you forgot to change the constraints mdp option.");
357 if (bWarn)
359 warning(wi, warn_buf);
361 else
363 warning_note(wi, warn_buf);
368 static void check_vel(gmx_mtop_t *mtop, rvec v[])
370 gmx_mtop_atomloop_all_t aloop;
371 const t_atom *atom;
372 int a;
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)
381 clear_rvec(v[a]);
386 static void check_shells_inputrec(gmx_mtop_t *mtop,
387 t_inputrec *ir,
388 warninp_t wi)
390 gmx_mtop_atomloop_all_t aloop;
391 const t_atom *atom;
392 int a, nshells = 0;
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)
401 nshells++;
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)
419 int nint, mb;
421 nint = 0;
422 for (mb = 0; mb < mtop->nmolblock; mb++)
424 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
427 return nint;
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;
438 int mb, mi;
439 t_molinfo *minew;
441 snew(order, *nmolinfo);
442 norder = 0;
443 for (mb = 0; mb < sys->nmolblock; mb++)
445 for (i = 0; i < norder; i++)
447 if (order[i] == sys->molblock[mb].type)
449 break;
452 if (i == norder)
454 /* This type did not occur yet, add it */
455 order[norder] = sys->molblock[mb].type;
456 /* Renumber the moltype in the topology */
457 norder++;
459 sys->molblock[mb].type = i;
462 /* We still need to reorder the molinfo structs */
463 snew(minew, norder);
464 for (mi = 0; mi < *nmolinfo; mi++)
466 for (i = 0; i < norder; i++)
468 if (order[i] == mi)
470 break;
473 if (i == norder)
475 done_mi(&(*molinfo)[mi]);
477 else
479 minew[i] = (*molinfo)[mi];
482 sfree(order);
483 sfree(*molinfo);
485 *nmolinfo = norder;
486 *molinfo = minew;
489 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
491 int m;
492 gmx_moltype_t *molt;
494 mtop->nmoltype = nmi;
495 snew(mtop->moltype, nmi);
496 for (m = 0; m < nmi; m++)
498 molt = &mtop->moltype[m];
499 molt->name = mi[m].name;
500 molt->atoms = mi[m].atoms;
501 /* ilists are copied later */
502 molt->cgs = mi[m].cgs;
503 molt->excls = mi[m].excls;
507 static void
508 new_status(const char *topfile, const char *topppfile, const char *confin,
509 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
510 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
511 gpp_atomtype_t atype, gmx_mtop_t *sys,
512 int *nmi, t_molinfo **mi, t_molinfo **intermolecular_interactions,
513 t_params plist[],
514 int *comb, double *reppow, real *fudgeQQ,
515 gmx_bool bMorse,
516 warninp_t wi)
518 t_molinfo *molinfo = nullptr;
519 int nmolblock;
520 gmx_molblock_t *molblock, *molbs;
521 int mb, i, nrmols, nmismatch;
522 char buf[STRLEN];
523 char warn_buf[STRLEN];
525 init_mtop(sys);
527 /* TOPOLOGY processing */
528 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
529 plist, comb, reppow, fudgeQQ,
530 atype, &nrmols, &molinfo, intermolecular_interactions,
532 &nmolblock, &molblock,
533 wi);
535 sys->nmolblock = 0;
536 snew(sys->molblock, nmolblock);
538 sys->natoms = 0;
539 for (mb = 0; mb < nmolblock; mb++)
541 if (sys->nmolblock > 0 &&
542 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
544 /* Merge consecutive blocks with the same molecule type */
545 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
546 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
548 else if (molblock[mb].nmol > 0)
550 /* Add a new molblock to the topology */
551 molbs = &sys->molblock[sys->nmolblock];
552 *molbs = molblock[mb];
553 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
554 molbs->nposres_xA = 0;
555 molbs->nposres_xB = 0;
556 sys->natoms += molbs->nmol*molbs->natoms_mol;
557 sys->nmolblock++;
560 if (sys->nmolblock == 0)
562 gmx_fatal(FARGS, "No molecules were defined in the system");
565 renumber_moltypes(sys, &nrmols, &molinfo);
567 if (bMorse)
569 convert_harmonics(nrmols, molinfo, atype);
572 if (ir->eDisre == edrNone)
574 i = rm_interactions(F_DISRES, nrmols, molinfo);
575 if (i > 0)
577 set_warning_line(wi, "unknown", -1);
578 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
579 warning_note(wi, warn_buf);
582 if (opts->bOrire == FALSE)
584 i = rm_interactions(F_ORIRES, nrmols, molinfo);
585 if (i > 0)
587 set_warning_line(wi, "unknown", -1);
588 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
589 warning_note(wi, warn_buf);
593 /* Copy structures from msys to sys */
594 molinfo2mtop(nrmols, molinfo, sys);
596 gmx_mtop_finalize(sys);
598 /* COORDINATE file processing */
599 if (bVerbose)
601 fprintf(stderr, "processing coordinates...\n");
604 t_topology *conftop;
605 rvec *x = nullptr;
606 rvec *v = nullptr;
607 snew(conftop, 1);
608 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
609 state->natoms = conftop->atoms.nr;
610 if (state->natoms != sys->natoms)
612 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
613 " does not match topology (%s, %d)",
614 confin, state->natoms, topfile, sys->natoms);
616 /* It would be nice to get rid of the copies below, but we don't know
617 * a priori if the number of atoms in confin matches what we expect.
619 state->flags |= (1 << estX);
620 if (v != NULL)
622 state->flags |= (1 << estV);
624 state_change_natoms(state, state->natoms);
625 for (int i = 0; i < state->natoms; i++)
627 copy_rvec(x[i], state->x[i]);
629 sfree(x);
630 if (v != nullptr)
632 for (int i = 0; i < state->natoms; i++)
634 copy_rvec(v[i], state->v[i]);
636 sfree(v);
638 /* This call fixes the box shape for runs with pressure scaling */
639 set_box_rel(ir, state);
641 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
642 done_top(conftop);
643 sfree(conftop);
645 if (nmismatch)
647 sprintf(buf, "%d non-matching atom name%s\n"
648 "atom names from %s will be used\n"
649 "atom names from %s will be ignored\n",
650 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
651 warning(wi, buf);
654 /* If using the group scheme, make sure charge groups are made whole to avoid errors
655 * in calculating charge group size later on
657 if (ir->cutoff_scheme == ecutsGROUP && ir->ePBC != epbcNONE)
659 // Need temporary rvec for coordinates
660 do_pbc_first_mtop(nullptr, ir->ePBC, state->box, sys, as_rvec_array(state->x.data()));
663 /* Do more checks, mostly related to constraints */
664 if (bVerbose)
666 fprintf(stderr, "double-checking input for internal consistency...\n");
669 int bHasNormalConstraints = 0 < (nint_ftype(sys, molinfo, F_CONSTR) +
670 nint_ftype(sys, molinfo, F_CONSTRNC));
671 int bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, molinfo, F_SETTLE);
672 double_check(ir, state->box,
673 bHasNormalConstraints,
674 bHasAnyConstraints,
675 wi);
678 if (bGenVel)
680 real *mass;
681 gmx_mtop_atomloop_all_t aloop;
682 const t_atom *atom;
684 snew(mass, state->natoms);
685 aloop = gmx_mtop_atomloop_all_init(sys);
686 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
688 mass[i] = atom->m;
691 if (opts->seed == -1)
693 opts->seed = static_cast<int>(gmx::makeRandomSeed());
694 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
696 state->flags |= (1 << estV);
697 maxwell_speed(opts->tempi, opts->seed, sys, as_rvec_array(state->v.data()));
699 stop_cm(stdout, state->natoms, mass, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()));
700 sfree(mass);
703 *nmi = nrmols;
704 *mi = molinfo;
707 static void copy_state(const char *slog, t_trxframe *fr,
708 gmx_bool bReadVel, t_state *state,
709 double *use_time)
711 int i;
713 if (fr->not_ok & FRAME_NOT_OK)
715 gmx_fatal(FARGS, "Can not start from an incomplete frame");
717 if (!fr->bX)
719 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
720 slog);
723 for (i = 0; i < state->natoms; i++)
725 copy_rvec(fr->x[i], state->x[i]);
727 if (bReadVel)
729 if (!fr->bV)
731 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
733 for (i = 0; i < state->natoms; i++)
735 copy_rvec(fr->v[i], state->v[i]);
738 if (fr->bBox)
740 copy_mat(fr->box, state->box);
743 *use_time = fr->time;
746 static void cont_status(const char *slog, const char *ener,
747 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
748 t_inputrec *ir, t_state *state,
749 gmx_mtop_t *sys,
750 const gmx_output_env_t *oenv)
751 /* If fr_time == -1 read the last frame available which is complete */
753 gmx_bool bReadVel;
754 t_trxframe fr;
755 t_trxstatus *fp;
756 int i;
757 double use_time;
759 bReadVel = (bNeedVel && !bGenVel);
761 fprintf(stderr,
762 "Reading Coordinates%s and Box size from old trajectory\n",
763 bReadVel ? ", Velocities" : "");
764 if (fr_time == -1)
766 fprintf(stderr, "Will read whole trajectory\n");
768 else
770 fprintf(stderr, "Will read till time %g\n", fr_time);
772 if (!bReadVel)
774 if (bGenVel)
776 fprintf(stderr, "Velocities generated: "
777 "ignoring velocities in input trajectory\n");
779 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
781 else
783 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
785 if (!fr.bV)
787 fprintf(stderr,
788 "\n"
789 "WARNING: Did not find a frame with velocities in file %s,\n"
790 " all velocities will be set to zero!\n\n", slog);
791 for (i = 0; i < sys->natoms; i++)
793 clear_rvec(state->v[i]);
795 close_trx(fp);
796 /* Search for a frame without velocities */
797 bReadVel = FALSE;
798 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
802 state->natoms = fr.natoms;
804 if (sys->natoms != state->natoms)
806 gmx_fatal(FARGS, "Number of atoms in Topology "
807 "is not the same as in Trajectory");
809 copy_state(slog, &fr, bReadVel, state, &use_time);
811 /* Find the appropriate frame */
812 while ((fr_time == -1 || fr.time < fr_time) &&
813 read_next_frame(oenv, fp, &fr))
815 copy_state(slog, &fr, bReadVel, state, &use_time);
818 close_trx(fp);
820 /* Set the relative box lengths for preserving the box shape.
821 * Note that this call can lead to differences in the last bit
822 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
824 set_box_rel(ir, state);
826 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
827 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
829 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
831 get_enx_state(ener, use_time, &sys->groups, ir, state);
832 preserve_box_shape(ir, state->box_rel, state->boxv);
836 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
837 const char *fn,
838 int rc_scaling, int ePBC,
839 rvec com,
840 warninp_t wi)
842 gmx_bool *hadAtom;
843 rvec *x, *v, *xp;
844 dvec sum;
845 double totmass;
846 t_topology *top;
847 matrix box, invbox;
848 int natoms, npbcdim = 0;
849 char warn_buf[STRLEN];
850 int a, i, ai, j, k, mb, nat_molb;
851 gmx_molblock_t *molb;
852 t_params *pr, *prfb;
853 t_atom *atom;
855 snew(top, 1);
856 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
857 natoms = top->atoms.nr;
858 done_top(top);
859 sfree(top);
860 if (natoms != mtop->natoms)
862 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);
863 warning(wi, warn_buf);
866 npbcdim = ePBC2npbcdim(ePBC);
867 clear_rvec(com);
868 if (rc_scaling != erscNO)
870 copy_mat(box, invbox);
871 for (j = npbcdim; j < DIM; j++)
873 clear_rvec(invbox[j]);
874 invbox[j][j] = 1;
876 gmx::invertBoxMatrix(invbox, invbox);
879 /* Copy the reference coordinates to mtop */
880 clear_dvec(sum);
881 totmass = 0;
882 a = 0;
883 snew(hadAtom, natoms);
884 for (mb = 0; mb < mtop->nmolblock; mb++)
886 molb = &mtop->molblock[mb];
887 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
888 pr = &(molinfo[molb->type].plist[F_POSRES]);
889 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
890 if (pr->nr > 0 || prfb->nr > 0)
892 atom = mtop->moltype[molb->type].atoms.atom;
893 for (i = 0; (i < pr->nr); i++)
895 ai = pr->param[i].ai();
896 if (ai >= natoms)
898 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
899 ai+1, *molinfo[molb->type].name, fn, natoms);
901 hadAtom[ai] = TRUE;
902 if (rc_scaling == erscCOM)
904 /* Determine the center of mass of the posres reference coordinates */
905 for (j = 0; j < npbcdim; j++)
907 sum[j] += atom[ai].m*x[a+ai][j];
909 totmass += atom[ai].m;
912 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
913 for (i = 0; (i < prfb->nr); i++)
915 ai = prfb->param[i].ai();
916 if (ai >= natoms)
918 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
919 ai+1, *molinfo[molb->type].name, fn, natoms);
921 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
923 /* Determine the center of mass of the posres reference coordinates */
924 for (j = 0; j < npbcdim; j++)
926 sum[j] += atom[ai].m*x[a+ai][j];
928 totmass += atom[ai].m;
931 if (!bTopB)
933 molb->nposres_xA = nat_molb;
934 snew(molb->posres_xA, molb->nposres_xA);
935 for (i = 0; i < nat_molb; i++)
937 copy_rvec(x[a+i], molb->posres_xA[i]);
940 else
942 molb->nposres_xB = nat_molb;
943 snew(molb->posres_xB, molb->nposres_xB);
944 for (i = 0; i < nat_molb; i++)
946 copy_rvec(x[a+i], molb->posres_xB[i]);
950 a += nat_molb;
952 if (rc_scaling == erscCOM)
954 if (totmass == 0)
956 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
958 for (j = 0; j < npbcdim; j++)
960 com[j] = sum[j]/totmass;
962 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]);
965 if (rc_scaling != erscNO)
967 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
969 for (mb = 0; mb < mtop->nmolblock; mb++)
971 molb = &mtop->molblock[mb];
972 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
973 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
975 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
976 for (i = 0; i < nat_molb; i++)
978 for (j = 0; j < npbcdim; j++)
980 if (rc_scaling == erscALL)
982 /* Convert from Cartesian to crystal coordinates */
983 xp[i][j] *= invbox[j][j];
984 for (k = j+1; k < npbcdim; k++)
986 xp[i][j] += invbox[k][j]*xp[i][k];
989 else if (rc_scaling == erscCOM)
991 /* Subtract the center of mass */
992 xp[i][j] -= com[j];
999 if (rc_scaling == erscCOM)
1001 /* Convert the COM from Cartesian to crystal coordinates */
1002 for (j = 0; j < npbcdim; j++)
1004 com[j] *= invbox[j][j];
1005 for (k = j+1; k < npbcdim; k++)
1007 com[j] += invbox[k][j]*com[k];
1013 sfree(x);
1014 sfree(v);
1015 sfree(hadAtom);
1018 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
1019 const char *fnA, const char *fnB,
1020 int rc_scaling, int ePBC,
1021 rvec com, rvec comB,
1022 warninp_t wi)
1024 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1025 /* It is safer to simply read the b-state posres rather than trying
1026 * to be smart and copy the positions.
1028 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1031 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
1032 t_inputrec *ir, warninp_t wi)
1034 int i;
1035 char warn_buf[STRLEN];
1037 if (ir->nwall > 0)
1039 fprintf(stderr, "Searching the wall atom type(s)\n");
1041 for (i = 0; i < ir->nwall; i++)
1043 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1044 if (ir->wall_atomtype[i] == NOTSET)
1046 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1047 warning_error(wi, warn_buf);
1052 static int nrdf_internal(t_atoms *atoms)
1054 int i, nmass, nrdf;
1056 nmass = 0;
1057 for (i = 0; i < atoms->nr; i++)
1059 /* Vsite ptype might not be set here yet, so also check the mass */
1060 if ((atoms->atom[i].ptype == eptAtom ||
1061 atoms->atom[i].ptype == eptNucleus)
1062 && atoms->atom[i].m > 0)
1064 nmass++;
1067 switch (nmass)
1069 case 0: nrdf = 0; break;
1070 case 1: nrdf = 0; break;
1071 case 2: nrdf = 1; break;
1072 default: nrdf = nmass*3 - 6; break;
1075 return nrdf;
1078 static void
1079 spline1d( double dx,
1080 double * y,
1081 int n,
1082 double * u,
1083 double * y2 )
1085 int i;
1086 double p, q;
1088 y2[0] = 0.0;
1089 u[0] = 0.0;
1091 for (i = 1; i < n-1; i++)
1093 p = 0.5*y2[i-1]+2.0;
1094 y2[i] = -0.5/p;
1095 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1096 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1099 y2[n-1] = 0.0;
1101 for (i = n-2; i >= 0; i--)
1103 y2[i] = y2[i]*y2[i+1]+u[i];
1108 static void
1109 interpolate1d( double xmin,
1110 double dx,
1111 double * ya,
1112 double * y2a,
1113 double x,
1114 double * y,
1115 double * y1)
1117 int ix;
1118 double a, b;
1120 ix = static_cast<int>((x-xmin)/dx);
1122 a = (xmin+(ix+1)*dx-x)/dx;
1123 b = (x-xmin-ix*dx)/dx;
1125 *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;
1126 *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 static void
1131 setup_cmap (int grid_spacing,
1132 int nc,
1133 real * grid,
1134 gmx_cmap_t * cmap_grid)
1136 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1138 int i, j, k, ii, jj, kk, idx;
1139 int offset;
1140 double dx, xmin, v, v1, v2, v12;
1141 double phi, psi;
1143 snew(tmp_u, 2*grid_spacing);
1144 snew(tmp_u2, 2*grid_spacing);
1145 snew(tmp_yy, 2*grid_spacing);
1146 snew(tmp_y1, 2*grid_spacing);
1147 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1148 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1150 dx = 360.0/grid_spacing;
1151 xmin = -180.0-dx*grid_spacing/2;
1153 for (kk = 0; kk < nc; kk++)
1155 /* Compute an offset depending on which cmap we are using
1156 * Offset will be the map number multiplied with the
1157 * grid_spacing * grid_spacing * 2
1159 offset = kk * grid_spacing * grid_spacing * 2;
1161 for (i = 0; i < 2*grid_spacing; i++)
1163 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1165 for (j = 0; j < 2*grid_spacing; j++)
1167 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1168 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1172 for (i = 0; i < 2*grid_spacing; i++)
1174 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1177 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1179 ii = i-grid_spacing/2;
1180 phi = ii*dx-180.0;
1182 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1184 jj = j-grid_spacing/2;
1185 psi = jj*dx-180.0;
1187 for (k = 0; k < 2*grid_spacing; k++)
1189 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1190 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1193 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1194 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1195 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1196 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1198 idx = ii*grid_spacing+jj;
1199 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1200 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1201 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1202 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1208 static void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1210 int i, nelem;
1212 cmap_grid->ngrid = ngrid;
1213 cmap_grid->grid_spacing = grid_spacing;
1214 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1216 snew(cmap_grid->cmapdata, ngrid);
1218 for (i = 0; i < cmap_grid->ngrid; i++)
1220 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1225 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1227 int count, count_mol, i, mb;
1228 gmx_molblock_t *molb;
1229 t_params *plist;
1230 char buf[STRLEN];
1232 count = 0;
1233 for (mb = 0; mb < mtop->nmolblock; mb++)
1235 count_mol = 0;
1236 molb = &mtop->molblock[mb];
1237 plist = mi[molb->type].plist;
1239 for (i = 0; i < F_NRE; i++)
1241 if (i == F_SETTLE)
1243 count_mol += 3*plist[i].nr;
1245 else if (interaction_function[i].flags & IF_CONSTRAINT)
1247 count_mol += plist[i].nr;
1251 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1253 sprintf(buf,
1254 "Molecule type '%s' has %d constraints.\n"
1255 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1256 *mi[molb->type].name, count_mol,
1257 nrdf_internal(&mi[molb->type].atoms));
1258 warning(wi, buf);
1260 count += molb->nmol*count_mol;
1263 return count;
1266 static real calc_temp(const gmx_mtop_t *mtop,
1267 const t_inputrec *ir,
1268 rvec *v)
1270 gmx_mtop_atomloop_all_t aloop;
1271 const t_atom *atom;
1272 int a;
1274 double sum_mv2 = 0;
1275 aloop = gmx_mtop_atomloop_all_init(mtop);
1276 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1278 sum_mv2 += atom->m*norm2(v[a]);
1281 double nrdf = 0;
1282 for (int g = 0; g < ir->opts.ngtc; g++)
1284 nrdf += ir->opts.nrdf[g];
1287 return sum_mv2/(nrdf*BOLTZ);
1290 static real get_max_reference_temp(const t_inputrec *ir,
1291 warninp_t wi)
1293 real ref_t;
1294 int i;
1295 gmx_bool bNoCoupl;
1297 ref_t = 0;
1298 bNoCoupl = FALSE;
1299 for (i = 0; i < ir->opts.ngtc; i++)
1301 if (ir->opts.tau_t[i] < 0)
1303 bNoCoupl = TRUE;
1305 else
1307 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1311 if (bNoCoupl)
1313 char buf[STRLEN];
1315 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.",
1316 ref_t);
1317 warning(wi, buf);
1320 return ref_t;
1323 /* Checks if there are unbound atoms in moleculetype molt.
1324 * Prints a note for each unbound atoms and a warning if any is present.
1326 static void checkForUnboundAtoms(const gmx_moltype_t *molt,
1327 gmx_bool bVerbose,
1328 warninp_t wi)
1330 const t_atoms *atoms = &molt->atoms;
1332 if (atoms->nr == 1)
1334 /* Only one atom, there can't be unbound atoms */
1335 return;
1338 std::vector<int> count(atoms->nr, 0);
1340 for (int ftype = 0; ftype < F_NRE; ftype++)
1342 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS) ||
1343 (interaction_function[ftype].flags & IF_CONSTRAINT) ||
1344 ftype == F_SETTLE)
1346 const t_ilist *il = &molt->ilist[ftype];
1347 int nral = NRAL(ftype);
1349 for (int i = 0; i < il->nr; i += 1 + nral)
1351 for (int j = 0; j < nral; j++)
1353 count[il->iatoms[i + 1 + j]]++;
1359 int numDanglingAtoms = 0;
1360 for (int a = 0; a < atoms->nr; a++)
1362 if (atoms->atom[a].ptype != eptVSite &&
1363 count[a] == 0)
1365 if (bVerbose)
1367 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",
1368 a + 1, *atoms->atomname[a], *molt->name);
1370 numDanglingAtoms++;
1374 if (numDanglingAtoms > 0)
1376 char buf[STRLEN];
1377 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.",
1378 *molt->name, numDanglingAtoms);
1379 warning_note(wi, buf);
1383 /* Checks all moleculetypes for unbound atoms */
1384 static void checkForUnboundAtoms(const gmx_mtop_t *mtop,
1385 gmx_bool bVerbose,
1386 warninp_t wi)
1388 for (int mt = 0; mt < mtop->nmoltype; mt++)
1390 checkForUnboundAtoms(&mtop->moltype[mt], bVerbose, wi);
1394 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1396 * The specific decoupled modes this routine check for are angle modes
1397 * where the two bonds are constrained and the atoms a both ends are only
1398 * involved in a single constraint; the mass of the two atoms needs to
1399 * differ by more than \p massFactorThreshold.
1401 static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
1402 const t_iparams *iparams,
1403 real massFactorThreshold)
1405 if (molt->ilist[F_CONSTR].nr == 0 && molt->ilist[F_CONSTRNC].nr == 0)
1407 return false;
1410 const t_atom * atom = molt->atoms.atom;
1412 int numFlexibleConstraints;
1413 t_blocka atomToConstraints = make_at2con(0, molt->atoms.nr,
1414 molt->ilist, iparams,
1415 FALSE,
1416 &numFlexibleConstraints);
1418 bool haveDecoupledMode = false;
1419 for (int ftype = 0; ftype < F_NRE; ftype++)
1421 if (interaction_function[ftype].flags & IF_ATYPE)
1423 const int nral = NRAL(ftype);
1424 const t_ilist *il = &molt->ilist[ftype];
1425 for (int i = 0; i < il->nr; i += 1 + nral)
1427 /* Here we check for the mass difference between the atoms
1428 * at both ends of the angle, that the atoms at the ends have
1429 * 1 contraint and the atom in the middle at least 3; we check
1430 * that the 3 atoms are linked by constraints below.
1431 * We check for at least three constraints for the middle atom,
1432 * since with only the two bonds in the angle, we have 3-atom
1433 * molecule, which has much more thermal exhange in this single
1434 * angle mode than molecules with more atoms.
1435 * Note that this check also catches molecules where more atoms
1436 * are connected to one or more atoms in the angle, but by
1437 * bond potentials instead of angles. But such cases will not
1438 * occur in "normal" molecules and it doesn't hurt running
1439 * those with higher accuracy settings as well.
1441 int a0 = il->iatoms[1 + i];
1442 int a1 = il->iatoms[1 + i + 1];
1443 int a2 = il->iatoms[1 + i + 2];
1444 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
1445 atom[a2].m > atom[a0].m*massFactorThreshold) &&
1446 atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
1447 atomToConstraints.index[a2 + 1] - atomToConstraints.index[a2] == 1 &&
1448 atomToConstraints.index[a1 + 1] - atomToConstraints.index[a1] >= 3)
1450 int constraint0 = atomToConstraints.a[atomToConstraints.index[a0]];
1451 int constraint2 = atomToConstraints.a[atomToConstraints.index[a2]];
1453 bool foundAtom0 = false;
1454 bool foundAtom2 = false;
1455 for (int conIndex = atomToConstraints.index[a1]; conIndex < atomToConstraints.index[a1 + 1]; conIndex++)
1457 if (atomToConstraints.a[conIndex] == constraint0)
1459 foundAtom0 = true;
1461 if (atomToConstraints.a[conIndex] == constraint2)
1463 foundAtom2 = true;
1466 if (foundAtom0 && foundAtom2)
1468 haveDecoupledMode = true;
1475 done_blocka(&atomToConstraints);
1477 return haveDecoupledMode;
1480 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1482 * When decoupled modes are present and the accuray in insufficient,
1483 * this routine issues a warning if the accuracy is insufficient.
1485 static void checkDecoupledModeAccuracy(const gmx_mtop_t *mtop,
1486 const t_inputrec *ir,
1487 warninp_t wi)
1489 /* We only have issues with decoupled modes with normal MD.
1490 * With stochastic dynamics equipartitioning is enforced strongly.
1492 if (!EI_MD(ir->eI))
1494 return;
1497 /* When atoms of very different mass are involved in an angle potential
1498 * and both bonds in the angle are constrained, the dynamic modes in such
1499 * angles have very different periods and significant energy exchange
1500 * takes several nanoseconds. Thus even a small amount of error in
1501 * different algorithms can lead to violation of equipartitioning.
1502 * The parameters below are mainly based on an all-atom chloroform model
1503 * with all bonds constrained. Equipartitioning is off by more than 1%
1504 * (need to run 5-10 ns) when the difference in mass between H and Cl
1505 * is more than a factor 13 and the accuracy is less than the thresholds
1506 * given below. This has been verified on some other molecules.
1508 * Note that the buffer and shake parameters have unit length and
1509 * energy/time, respectively, so they will "only" work correctly
1510 * for atomistic force fields using MD units.
1512 const real massFactorThreshold = 13.0;
1513 const real bufferToleranceThreshold = 1e-4;
1514 const int lincsIterationThreshold = 2;
1515 const int lincsOrderThreshold = 4;
1516 const real shakeToleranceThreshold = 0.005*ir->delta_t;
1518 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1519 bool shakeWithSufficientTolerance = (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1*shakeToleranceThreshold);
1520 if (ir->cutoff_scheme == ecutsVERLET &&
1521 ir->verletbuf_tol <= 1.1*bufferToleranceThreshold &&
1522 (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1524 return;
1527 bool haveDecoupledMode = false;
1528 for (int mt = 0; mt < mtop->nmoltype; mt++)
1530 if (haveDecoupledModeInMol(&mtop->moltype[mt], mtop->ffparams.iparams,
1531 massFactorThreshold))
1533 haveDecoupledMode = true;
1537 if (haveDecoupledMode)
1539 char modeMessage[STRLEN];
1540 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.",
1541 massFactorThreshold);
1542 char buf[STRLEN];
1543 if (ir->cutoff_scheme == ecutsVERLET)
1545 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",
1546 modeMessage,
1547 ei_names[eiSD1],
1548 bufferToleranceThreshold,
1549 lincsIterationThreshold, lincsOrderThreshold,
1550 shakeToleranceThreshold);
1552 else
1554 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 modeMessage,
1556 ecutscheme_names[ecutsVERLET]);
1558 warning(wi, buf);
1562 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1563 t_inputrec *ir,
1564 real buffer_temp,
1565 matrix box,
1566 warninp_t wi)
1568 real rlist_1x1;
1569 int n_nonlin_vsite;
1570 char warn_buf[STRLEN];
1572 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1574 /* Calculate the buffer size for simple atom vs atoms list */
1575 VerletbufListSetup listSetup1x1;
1576 listSetup1x1.cluster_size_i = 1;
1577 listSetup1x1.cluster_size_j = 1;
1578 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1,
1579 buffer_temp, &listSetup1x1,
1580 &n_nonlin_vsite, &rlist_1x1);
1582 /* Set the pair-list buffer size in ir */
1583 VerletbufListSetup listSetup4x4 =
1584 verletbufGetSafeListSetup(ListSetupType::CpuNoSimd);
1585 calc_verlet_buffer_size(mtop, det(box), ir, ir->nstlist, ir->nstlist - 1,
1586 buffer_temp, &listSetup4x4,
1587 &n_nonlin_vsite, &ir->rlist);
1589 if (n_nonlin_vsite > 0)
1591 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);
1592 warning_note(wi, warn_buf);
1595 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1596 1, 1, rlist_1x1, rlist_1x1-std::max(ir->rvdw, ir->rcoulomb));
1598 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1599 listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j,
1600 ir->rlist, ir->rlist-std::max(ir->rvdw, ir->rcoulomb));
1602 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1604 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1606 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)));
1610 int gmx_grompp(int argc, char *argv[])
1612 const char *desc[] = {
1613 "[THISMODULE] (the gromacs preprocessor)",
1614 "reads a molecular topology file, checks the validity of the",
1615 "file, expands the topology from a molecular description to an atomic",
1616 "description. The topology file contains information about",
1617 "molecule types and the number of molecules, the preprocessor",
1618 "copies each molecule as needed. ",
1619 "There is no limitation on the number of molecule types. ",
1620 "Bonds and bond-angles can be converted into constraints, separately",
1621 "for hydrogens and heavy atoms.",
1622 "Then a coordinate file is read and velocities can be generated",
1623 "from a Maxwellian distribution if requested.",
1624 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1625 "(eg. number of MD steps, time step, cut-off), and others such as",
1626 "NEMD parameters, which are corrected so that the net acceleration",
1627 "is zero.",
1628 "Eventually a binary file is produced that can serve as the sole input",
1629 "file for the MD program.[PAR]",
1631 "[THISMODULE] uses the atom names from the topology file. The atom names",
1632 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1633 "warnings when they do not match the atom names in the topology.",
1634 "Note that the atom names are irrelevant for the simulation as",
1635 "only the atom types are used for generating interaction parameters.[PAR]",
1637 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1638 "etc. The preprocessor supports the following keywords::",
1640 " #ifdef VARIABLE",
1641 " #ifndef VARIABLE",
1642 " #else",
1643 " #endif",
1644 " #define VARIABLE",
1645 " #undef VARIABLE",
1646 " #include \"filename\"",
1647 " #include <filename>",
1649 "The functioning of these statements in your topology may be modulated by",
1650 "using the following two flags in your [REF].mdp[ref] file::",
1652 " define = -DVARIABLE1 -DVARIABLE2",
1653 " include = -I/home/john/doe",
1655 "For further information a C-programming textbook may help you out.",
1656 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1657 "topology file written out so that you can verify its contents.[PAR]",
1659 "When using position restraints, a file with restraint coordinates",
1660 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1661 "for [TT]-c[tt]). For free energy calculations, separate reference",
1662 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1663 "otherwise they will be equal to those of the A topology.[PAR]",
1665 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1666 "The last frame with coordinates and velocities will be read,",
1667 "unless the [TT]-time[tt] option is used. Only if this information",
1668 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1669 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1670 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1671 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1672 "variables.[PAR]",
1674 "[THISMODULE] can be used to restart simulations (preserving",
1675 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1676 "However, for simply changing the number of run steps to extend",
1677 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1678 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1679 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1680 "like output frequency, then supplying the checkpoint file to",
1681 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1682 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1683 "the ensemble (if possible) still requires passing the checkpoint",
1684 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1686 "By default, all bonded interactions which have constant energy due to",
1687 "virtual site constructions will be removed. If this constant energy is",
1688 "not zero, this will result in a shift in the total energy. All bonded",
1689 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1690 "all constraints for distances which will be constant anyway because",
1691 "of virtual site constructions will be removed. If any constraints remain",
1692 "which involve virtual sites, a fatal error will result.[PAR]"
1694 "To verify your run input file, please take note of all warnings",
1695 "on the screen, and correct where necessary. Do also look at the contents",
1696 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1697 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1698 "with the [TT]-debug[tt] option which will give you more information",
1699 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1700 "can see the contents of the run input file with the [gmx-dump]",
1701 "program. [gmx-check] can be used to compare the contents of two",
1702 "run input files.[PAR]"
1704 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1705 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1706 "harmless, but usually they are not. The user is advised to carefully",
1707 "interpret the output messages before attempting to bypass them with",
1708 "this option."
1710 t_gromppopts *opts;
1711 gmx_mtop_t *sys;
1712 int nmi;
1713 t_molinfo *mi, *intermolecular_interactions;
1714 gpp_atomtype_t atype;
1715 int nvsite, comb, mt;
1716 t_params *plist;
1717 real fudgeQQ;
1718 double reppow;
1719 const char *mdparin;
1720 int ntype;
1721 gmx_bool bNeedVel, bGenVel;
1722 gmx_bool have_atomnumber;
1723 gmx_output_env_t *oenv;
1724 gmx_bool bVerbose = FALSE;
1725 warninp_t wi;
1726 char warn_buf[STRLEN];
1728 t_filenm fnm[] = {
1729 { efMDP, nullptr, nullptr, ffREAD },
1730 { efMDP, "-po", "mdout", ffWRITE },
1731 { efSTX, "-c", nullptr, ffREAD },
1732 { efSTX, "-r", "restraint", ffOPTRD },
1733 { efSTX, "-rb", "restraint", ffOPTRD },
1734 { efNDX, nullptr, nullptr, ffOPTRD },
1735 { efTOP, nullptr, nullptr, ffREAD },
1736 { efTOP, "-pp", "processed", ffOPTWR },
1737 { efTPR, "-o", nullptr, ffWRITE },
1738 { efTRN, "-t", nullptr, ffOPTRD },
1739 { efEDR, "-e", nullptr, ffOPTRD },
1740 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1741 { efGRO, "-imd", "imdgroup", ffOPTWR },
1742 { efTRN, "-ref", "rotref", ffOPTRW }
1744 #define NFILE asize(fnm)
1746 /* Command line options */
1747 gmx_bool bRenum = TRUE;
1748 gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1749 int i, maxwarn = 0;
1750 real fr_time = -1;
1751 t_pargs pa[] = {
1752 { "-v", FALSE, etBOOL, {&bVerbose},
1753 "Be loud and noisy" },
1754 { "-time", FALSE, etREAL, {&fr_time},
1755 "Take frame at or first after this time." },
1756 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1757 "Remove constant bonded interactions with virtual sites" },
1758 { "-maxwarn", FALSE, etINT, {&maxwarn},
1759 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1760 { "-zero", FALSE, etBOOL, {&bZero},
1761 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1762 { "-renum", FALSE, etBOOL, {&bRenum},
1763 "Renumber atomtypes and minimize number of atomtypes" }
1766 /* Parse the command line */
1767 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1768 asize(desc), desc, 0, nullptr, &oenv))
1770 return 0;
1773 /* Initiate some variables */
1774 gmx::MDModules mdModules;
1775 t_inputrec irInstance;
1776 t_inputrec *ir = &irInstance;
1777 snew(opts, 1);
1778 snew(opts->include, STRLEN);
1779 snew(opts->define, STRLEN);
1781 wi = init_warning(TRUE, maxwarn);
1783 /* PARAMETER file processing */
1784 mdparin = opt2fn("-f", NFILE, fnm);
1785 set_warning_line(wi, mdparin, -1);
1788 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1790 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1792 if (bVerbose)
1794 fprintf(stderr, "checking input for internal consistency...\n");
1796 check_ir(mdparin, ir, opts, wi);
1798 if (ir->ld_seed == -1)
1800 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1801 fprintf(stderr, "Setting the LD random seed to %" GMX_PRId64 "\n", ir->ld_seed);
1804 if (ir->expandedvals->lmc_seed == -1)
1806 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1807 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1810 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1811 bGenVel = (bNeedVel && opts->bGenVel);
1812 if (bGenVel && ir->bContinuation)
1814 sprintf(warn_buf,
1815 "Generating velocities is inconsistent with attempting "
1816 "to continue a previous run. Choose only one of "
1817 "gen-vel = yes and continuation = yes.");
1818 warning_error(wi, warn_buf);
1821 snew(plist, F_NRE);
1822 init_plist(plist);
1823 snew(sys, 1);
1824 atype = init_atomtype();
1825 if (debug)
1827 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1830 const char *fn = ftp2fn(efTOP, NFILE, fnm);
1831 if (!gmx_fexist(fn))
1833 gmx_fatal(FARGS, "%s does not exist", fn);
1836 t_state state;
1837 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1838 opts, ir, bZero, bGenVel, bVerbose, &state,
1839 atype, sys, &nmi, &mi, &intermolecular_interactions,
1840 plist, &comb, &reppow, &fudgeQQ,
1841 opts->bMorse,
1842 wi);
1844 if (debug)
1846 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1849 nvsite = 0;
1850 /* set parameters for virtual site construction (not for vsiten) */
1851 for (mt = 0; mt < sys->nmoltype; mt++)
1853 nvsite +=
1854 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1856 /* now throw away all obsolete bonds, angles and dihedrals: */
1857 /* note: constraints are ALWAYS removed */
1858 if (nvsite)
1860 for (mt = 0; mt < sys->nmoltype; mt++)
1862 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1866 if (ir->cutoff_scheme == ecutsVERLET)
1868 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1869 ecutscheme_names[ir->cutoff_scheme]);
1871 /* Remove all charge groups */
1872 gmx_mtop_remove_chargegroups(sys);
1875 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1877 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1879 sprintf(warn_buf, "Can not do %s with %s, use %s",
1880 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1881 warning_error(wi, warn_buf);
1883 if (ir->bPeriodicMols)
1885 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1886 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1887 warning_error(wi, warn_buf);
1891 if (EI_SD (ir->eI) && ir->etc != etcNO)
1893 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1896 /* If we are doing QM/MM, check that we got the atom numbers */
1897 have_atomnumber = TRUE;
1898 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1900 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1902 if (!have_atomnumber && ir->bQMMM)
1904 warning_error(wi,
1905 "\n"
1906 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1907 "field you are using does not contain atom numbers fields. This is an\n"
1908 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1909 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1910 "an integer just before the mass column in ffXXXnb.itp.\n"
1911 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1914 /* Check for errors in the input now, since they might cause problems
1915 * during processing further down.
1917 check_warning_error(wi, FARGS);
1919 if (nint_ftype(sys, mi, F_POSRES) > 0 ||
1920 nint_ftype(sys, mi, F_FBPOSRES) > 0)
1922 if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
1924 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.",
1925 EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
1926 warning_note(wi, warn_buf);
1929 const char *fn = opt2fn("-r", NFILE, fnm);
1930 const char *fnB;
1932 if (!gmx_fexist(fn))
1934 gmx_fatal(FARGS,
1935 "Cannot find position restraint file %s (option -r).\n"
1936 "From GROMACS-2018, you need to specify the position restraint "
1937 "coordinate files explicitly to avoid mistakes, although you can "
1938 "still use the same file as you specify for the -c option.", fn);
1941 if (opt2bSet("-rb", NFILE, fnm))
1943 fnB = opt2fn("-rb", NFILE, fnm);
1944 if (!gmx_fexist(fnB))
1946 gmx_fatal(FARGS,
1947 "Cannot find B-state position restraint file %s (option -rb).\n"
1948 "From GROMACS-2018, you need to specify the position restraint "
1949 "coordinate files explicitly to avoid mistakes, although you can "
1950 "still use the same file as you specify for the -c option.", fn);
1953 else
1955 fnB = fn;
1958 if (bVerbose)
1960 fprintf(stderr, "Reading position restraint coords from %s", fn);
1961 if (strcmp(fn, fnB) == 0)
1963 fprintf(stderr, "\n");
1965 else
1967 fprintf(stderr, " and %s\n", fnB);
1970 gen_posres(sys, mi, fn, fnB,
1971 ir->refcoord_scaling, ir->ePBC,
1972 ir->posres_com, ir->posres_comB,
1973 wi);
1976 /* If we are using CMAP, setup the pre-interpolation grid */
1977 if (plist[F_CMAP].ncmap > 0)
1979 init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
1980 setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
1983 set_wall_atomtype(atype, opts, ir, wi);
1984 if (bRenum)
1986 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1987 get_atomtype_ntypes(atype);
1990 if (ir->implicit_solvent)
1992 gmx_fatal(FARGS, "Implicit solvation is no longer supported");
1995 /* PELA: Copy the atomtype data to the topology atomtype list */
1996 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1998 if (debug)
2000 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
2003 if (bVerbose)
2005 fprintf(stderr, "converting bonded parameters...\n");
2008 ntype = get_atomtype_ntypes(atype);
2009 convert_params(ntype, plist, mi, intermolecular_interactions,
2010 comb, reppow, fudgeQQ, sys);
2012 if (debug)
2014 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
2017 /* set ptype to VSite for virtual sites */
2018 for (mt = 0; mt < sys->nmoltype; mt++)
2020 set_vsites_ptype(FALSE, &sys->moltype[mt]);
2022 if (debug)
2024 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
2026 /* Check velocity for virtual sites and shells */
2027 if (bGenVel)
2029 check_vel(sys, as_rvec_array(state.v.data()));
2032 /* check for shells and inpurecs */
2033 check_shells_inputrec(sys, ir, wi);
2035 /* check masses */
2036 check_mol(sys, wi);
2038 checkForUnboundAtoms(sys, bVerbose, wi);
2040 for (i = 0; i < sys->nmoltype; i++)
2042 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
2045 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2047 check_bonds_timestep(sys, ir->delta_t, wi);
2050 checkDecoupledModeAccuracy(sys, ir, wi);
2052 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2054 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.");
2057 check_warning_error(wi, FARGS);
2059 if (bVerbose)
2061 fprintf(stderr, "initialising group options...\n");
2063 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
2064 sys, bVerbose, ir,
2065 wi);
2067 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2069 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2071 real buffer_temp;
2073 if (EI_MD(ir->eI) && ir->etc == etcNO)
2075 if (bGenVel)
2077 buffer_temp = opts->tempi;
2079 else
2081 buffer_temp = calc_temp(sys, ir, as_rvec_array(state.v.data()));
2083 if (buffer_temp > 0)
2085 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
2086 warning_note(wi, warn_buf);
2088 else
2090 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
2091 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
2092 warning_note(wi, warn_buf);
2095 else
2097 buffer_temp = get_max_reference_temp(ir, wi);
2100 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2102 /* NVE with initial T=0: we add a fixed ratio to rlist.
2103 * Since we don't actually use verletbuf_tol, we set it to -1
2104 * so it can't be misused later.
2106 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2107 ir->verletbuf_tol = -1;
2109 else
2111 /* We warn for NVE simulations with a drift tolerance that
2112 * might result in a 1(.1)% drift over the total run-time.
2113 * Note that we can't warn when nsteps=0, since we don't
2114 * know how many steps the user intends to run.
2116 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
2117 ir->nsteps > 0)
2119 const real driftTolerance = 0.01;
2120 /* We use 2 DOF per atom = 2kT pot+kin energy,
2121 * to be on the safe side with constraints.
2123 const real totalEnergyDriftPerAtomPerPicosecond = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
2125 if (ir->verletbuf_tol > 1.1*driftTolerance*totalEnergyDriftPerAtomPerPicosecond)
2127 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.",
2128 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
2129 (int)(ir->verletbuf_tol/totalEnergyDriftPerAtomPerPicosecond*100 + 0.5),
2130 (int)(100*driftTolerance + 0.5),
2131 driftTolerance*totalEnergyDriftPerAtomPerPicosecond);
2132 warning_note(wi, warn_buf);
2136 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
2141 /* Init the temperature coupling state */
2142 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2144 if (bVerbose)
2146 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
2148 check_eg_vs_cg(sys);
2150 if (debug)
2152 pr_symtab(debug, 0, "After index", &sys->symtab);
2155 triple_check(mdparin, ir, sys, wi);
2156 close_symtab(&sys->symtab);
2157 if (debug)
2159 pr_symtab(debug, 0, "After close", &sys->symtab);
2162 /* make exclusions between QM atoms */
2163 if (ir->bQMMM)
2165 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
2167 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
2169 else
2171 generate_qmexcl(sys, ir, wi);
2175 if (ftp2bSet(efTRN, NFILE, fnm))
2177 if (bVerbose)
2179 fprintf(stderr, "getting data from old trajectory ...\n");
2181 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
2182 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
2185 if (ir->ePBC == epbcXY && ir->nwall != 2)
2187 clear_rvec(state.box[ZZ]);
2190 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
2192 set_warning_line(wi, mdparin, -1);
2193 check_chargegroup_radii(sys, ir, as_rvec_array(state.x.data()), wi);
2196 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2198 /* Calculate the optimal grid dimensions */
2199 matrix scaledBox;
2200 EwaldBoxZScaler boxScaler(*ir);
2201 boxScaler.scaleBox(state.box, scaledBox);
2203 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2205 /* Mark fourier_spacing as not used */
2206 ir->fourier_spacing = 0;
2208 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2210 set_warning_line(wi, mdparin, -1);
2211 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2213 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2214 calcFftGrid(stdout, scaledBox, ir->fourier_spacing, minGridSize,
2215 &(ir->nkx), &(ir->nky), &(ir->nkz));
2216 if (ir->nkx < minGridSize ||
2217 ir->nky < minGridSize ||
2218 ir->nkz < minGridSize)
2220 warning_error(wi, "The PME grid size should be >= 2*(pme-order - 1); either manually increase the grid size or decrease pme-order");
2224 /* MRS: eventually figure out better logic for initializing the fep
2225 values that makes declaring the lambda and declaring the state not
2226 potentially conflict if not handled correctly. */
2227 if (ir->efep != efepNO)
2229 state.fep_state = ir->fepvals->init_fep_state;
2230 for (i = 0; i < efptNR; i++)
2232 /* init_lambda trumps state definitions*/
2233 if (ir->fepvals->init_lambda >= 0)
2235 state.lambda[i] = ir->fepvals->init_lambda;
2237 else
2239 if (ir->fepvals->all_lambda[i] == nullptr)
2241 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2243 else
2245 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2251 struct pull_t *pull = nullptr;
2253 if (ir->bPull)
2255 pull = set_pull_init(ir, sys, as_rvec_array(state.x.data()), state.box, state.lambda[efptMASS], oenv);
2258 /* Modules that supply external potential for pull coordinates
2259 * should register those potentials here. finish_pull() will check
2260 * that providers have been registerd for all external potentials.
2263 if (ir->bDoAwh)
2265 setStateDependentAwhParams(ir->awhParams, ir->pull, pull,
2266 state.box, ir->ePBC, &ir->opts, wi);
2269 if (ir->bPull)
2271 finish_pull(pull);
2274 if (ir->bRot)
2276 set_reference_positions(ir->rot, as_rvec_array(state.x.data()), state.box,
2277 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2278 wi);
2281 /* reset_multinr(sys); */
2283 if (EEL_PME(ir->coulombtype))
2285 float ratio = pme_load_estimate(sys, ir, state.box);
2286 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2287 /* With free energy we might need to do PME both for the A and B state
2288 * charges. This will double the cost, but the optimal performance will
2289 * then probably be at a slightly larger cut-off and grid spacing.
2291 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2292 (ir->efep != efepNO && ratio > 2.0/3.0))
2294 warning_note(wi,
2295 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2296 "and for highly parallel simulations between 0.25 and 0.33,\n"
2297 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2298 if (ir->efep != efepNO)
2300 warning_note(wi,
2301 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2307 char warn_buf[STRLEN];
2308 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2309 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2310 if (cio > 2000)
2312 set_warning_line(wi, mdparin, -1);
2313 warning_note(wi, warn_buf);
2315 else
2317 printf("%s\n", warn_buf);
2321 if (bVerbose)
2323 fprintf(stderr, "writing run input file...\n");
2326 done_warning(wi, FARGS);
2327 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, sys);
2329 /* Output IMD group, if bIMD is TRUE */
2330 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2332 sfree(opts->define);
2333 sfree(opts->include);
2334 sfree(opts);
2335 done_plist(plist);
2336 sfree(plist);
2337 sfree(mi);
2338 done_atomtype(atype);
2339 done_mtop(sys);
2340 sfree(sys);
2341 done_inputrec_strings();
2343 return 0;