Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.cpp
blobb9916e5f86cacac5c0852876884cfac62b67a91a
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, 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/commandline/pargs.h"
52 #include "gromacs/ewald/pme.h"
53 #include "gromacs/fft/calcgrid.h"
54 #include "gromacs/fileio/confio.h"
55 #include "gromacs/fileio/enxio.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/fileio/warninp.h"
59 #include "gromacs/gmxpreprocess/add_par.h"
60 #include "gromacs/gmxpreprocess/convparm.h"
61 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
62 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
63 #include "gromacs/gmxpreprocess/grompp-impl.h"
64 #include "gromacs/gmxpreprocess/notset.h"
65 #include "gromacs/gmxpreprocess/readir.h"
66 #include "gromacs/gmxpreprocess/tomorse.h"
67 #include "gromacs/gmxpreprocess/topio.h"
68 #include "gromacs/gmxpreprocess/toputil.h"
69 #include "gromacs/gmxpreprocess/vsite_parm.h"
70 #include "gromacs/imd/imd.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/invertmatrix.h"
73 #include "gromacs/math/units.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/calc_verletbuf.h"
76 #include "gromacs/mdlib/compute_io.h"
77 #include "gromacs/mdlib/constr.h"
78 #include "gromacs/mdlib/genborn.h"
79 #include "gromacs/mdlib/perf_est.h"
80 #include "gromacs/mdrunutility/mdmodules.h"
81 #include "gromacs/mdtypes/inputrec.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/mdtypes/nblist.h"
84 #include "gromacs/mdtypes/state.h"
85 #include "gromacs/pbcutil/boxutilities.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/pulling/pull.h"
88 #include "gromacs/random/seed.h"
89 #include "gromacs/topology/ifunc.h"
90 #include "gromacs/topology/mtop_util.h"
91 #include "gromacs/topology/symtab.h"
92 #include "gromacs/topology/topology.h"
93 #include "gromacs/trajectory/trajectoryframe.h"
94 #include "gromacs/utility/arraysize.h"
95 #include "gromacs/utility/cstringutil.h"
96 #include "gromacs/utility/exceptions.h"
97 #include "gromacs/utility/fatalerror.h"
98 #include "gromacs/utility/futil.h"
99 #include "gromacs/utility/gmxassert.h"
100 #include "gromacs/utility/smalloc.h"
101 #include "gromacs/utility/snprintf.h"
103 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
105 int i, n;
107 n = 0;
108 /* For all the molecule types */
109 for (i = 0; i < nrmols; i++)
111 n += mols[i].plist[ifunc].nr;
112 mols[i].plist[ifunc].nr = 0;
114 return n;
117 static int check_atom_names(const char *fn1, const char *fn2,
118 gmx_mtop_t *mtop, const t_atoms *at)
120 int mb, m, i, j, nmismatch;
121 t_atoms *tat;
122 #define MAXMISMATCH 20
124 if (mtop->natoms != at->nr)
126 gmx_incons("comparing atom names");
129 nmismatch = 0;
130 i = 0;
131 for (mb = 0; mb < mtop->nmolblock; mb++)
133 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
134 for (m = 0; m < mtop->molblock[mb].nmol; m++)
136 for (j = 0; j < tat->nr; j++)
138 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
140 if (nmismatch < MAXMISMATCH)
142 fprintf(stderr,
143 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
144 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
146 else if (nmismatch == MAXMISMATCH)
148 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
150 nmismatch++;
152 i++;
157 return nmismatch;
160 static void check_eg_vs_cg(gmx_mtop_t *mtop)
162 int astart, mb, m, cg, j, firstj;
163 unsigned char firsteg, eg;
164 gmx_moltype_t *molt;
166 /* Go through all the charge groups and make sure all their
167 * atoms are in the same energy group.
170 astart = 0;
171 for (mb = 0; mb < mtop->nmolblock; mb++)
173 molt = &mtop->moltype[mtop->molblock[mb].type];
174 for (m = 0; m < mtop->molblock[mb].nmol; m++)
176 for (cg = 0; cg < molt->cgs.nr; cg++)
178 /* Get the energy group of the first atom in this charge group */
179 firstj = astart + molt->cgs.index[cg];
180 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
181 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
183 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
184 if (eg != firsteg)
186 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
187 firstj+1, astart+j+1, cg+1, *molt->name);
191 astart += molt->atoms.nr;
196 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
198 int maxsize, cg;
199 char warn_buf[STRLEN];
201 maxsize = 0;
202 for (cg = 0; cg < cgs->nr; cg++)
204 maxsize = std::max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
207 if (maxsize > MAX_CHARGEGROUP_SIZE)
209 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
211 else if (maxsize > 10)
213 set_warning_line(wi, topfn, -1);
214 sprintf(warn_buf,
215 "The largest charge group contains %d atoms.\n"
216 "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"
217 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
218 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
219 maxsize);
220 warning_note(wi, warn_buf);
224 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
226 /* This check is not intended to ensure accurate integration,
227 * rather it is to signal mistakes in the mdp settings.
228 * A common mistake is to forget to turn on constraints
229 * for MD after energy minimization with flexible bonds.
230 * This check can also detect too large time steps for flexible water
231 * models, but such errors will often be masked by the constraints
232 * mdp options, which turns flexible water into water with bond constraints,
233 * but without an angle constraint. Unfortunately such incorrect use
234 * of water models can not easily be detected without checking
235 * for specific model names.
237 * The stability limit of leap-frog or velocity verlet is 4.44 steps
238 * per oscillational period.
239 * But accurate bonds distributions are lost far before that limit.
240 * To allow relatively common schemes (although not common with Gromacs)
241 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
242 * we set the note limit to 10.
244 int min_steps_warn = 5;
245 int min_steps_note = 10;
246 t_iparams *ip;
247 int molt;
248 gmx_moltype_t *moltype, *w_moltype;
249 t_atom *atom;
250 t_ilist *ilist, *ilb, *ilc, *ils;
251 int ftype;
252 int i, a1, a2, w_a1, w_a2, j;
253 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
254 gmx_bool bFound, bWater, bWarn;
255 char warn_buf[STRLEN];
257 ip = mtop->ffparams.iparams;
259 twopi2 = gmx::square(2*M_PI);
261 limit2 = gmx::square(min_steps_note*dt);
263 w_a1 = w_a2 = -1;
264 w_period2 = -1.0;
266 w_moltype = nullptr;
267 for (molt = 0; molt < mtop->nmoltype; molt++)
269 moltype = &mtop->moltype[molt];
270 atom = moltype->atoms.atom;
271 ilist = moltype->ilist;
272 ilc = &ilist[F_CONSTR];
273 ils = &ilist[F_SETTLE];
274 for (ftype = 0; ftype < F_NRE; ftype++)
276 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
278 continue;
281 ilb = &ilist[ftype];
282 for (i = 0; i < ilb->nr; i += 3)
284 fc = ip[ilb->iatoms[i]].harmonic.krA;
285 re = ip[ilb->iatoms[i]].harmonic.rA;
286 if (ftype == F_G96BONDS)
288 /* Convert squared sqaure fc to harmonic fc */
289 fc = 2*fc*re;
291 a1 = ilb->iatoms[i+1];
292 a2 = ilb->iatoms[i+2];
293 m1 = atom[a1].m;
294 m2 = atom[a2].m;
295 if (fc > 0 && m1 > 0 && m2 > 0)
297 period2 = twopi2*m1*m2/((m1 + m2)*fc);
299 else
301 period2 = GMX_FLOAT_MAX;
303 if (debug)
305 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
306 fc, m1, m2, std::sqrt(period2));
308 if (period2 < limit2)
310 bFound = FALSE;
311 for (j = 0; j < ilc->nr; j += 3)
313 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
314 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
316 bFound = TRUE;
319 for (j = 0; j < ils->nr; j += 4)
321 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
322 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
324 bFound = TRUE;
327 if (!bFound &&
328 (w_moltype == nullptr || period2 < w_period2))
330 w_moltype = moltype;
331 w_a1 = a1;
332 w_a2 = a2;
333 w_period2 = period2;
340 if (w_moltype != nullptr)
342 bWarn = (w_period2 < gmx::square(min_steps_warn*dt));
343 /* A check that would recognize most water models */
344 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
345 w_moltype->atoms.nr <= 5);
346 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"
347 "%s",
348 *w_moltype->name,
349 w_a1+1, *w_moltype->atoms.atomname[w_a1],
350 w_a2+1, *w_moltype->atoms.atomname[w_a2],
351 std::sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
352 bWater ?
353 "Maybe you asked for fexible water." :
354 "Maybe you forgot to change the constraints mdp option.");
355 if (bWarn)
357 warning(wi, warn_buf);
359 else
361 warning_note(wi, warn_buf);
366 static void check_vel(gmx_mtop_t *mtop, rvec v[])
368 gmx_mtop_atomloop_all_t aloop;
369 const t_atom *atom;
370 int a;
372 aloop = gmx_mtop_atomloop_all_init(mtop);
373 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
375 if (atom->ptype == eptShell ||
376 atom->ptype == eptBond ||
377 atom->ptype == eptVSite)
379 clear_rvec(v[a]);
384 static void check_shells_inputrec(gmx_mtop_t *mtop,
385 t_inputrec *ir,
386 warninp_t wi)
388 gmx_mtop_atomloop_all_t aloop;
389 const t_atom *atom;
390 int a, nshells = 0;
391 char warn_buf[STRLEN];
393 aloop = gmx_mtop_atomloop_all_init(mtop);
394 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
396 if (atom->ptype == eptShell ||
397 atom->ptype == eptBond)
399 nshells++;
402 if ((nshells > 0) && (ir->nstcalcenergy != 1))
404 set_warning_line(wi, "unknown", -1);
405 snprintf(warn_buf, STRLEN,
406 "There are %d shells, changing nstcalcenergy from %d to 1",
407 nshells, ir->nstcalcenergy);
408 ir->nstcalcenergy = 1;
409 warning(wi, warn_buf);
413 /* TODO Decide whether this function can be consolidated with
414 * gmx_mtop_ftype_count */
415 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
417 int nint, mb;
419 nint = 0;
420 for (mb = 0; mb < mtop->nmolblock; mb++)
422 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
425 return nint;
428 /* This routine reorders the molecule type array
429 * in the order of use in the molblocks,
430 * unused molecule types are deleted.
432 static void renumber_moltypes(gmx_mtop_t *sys,
433 int *nmolinfo, t_molinfo **molinfo)
435 int *order, norder, i;
436 int mb, mi;
437 t_molinfo *minew;
439 snew(order, *nmolinfo);
440 norder = 0;
441 for (mb = 0; mb < sys->nmolblock; mb++)
443 for (i = 0; i < norder; i++)
445 if (order[i] == sys->molblock[mb].type)
447 break;
450 if (i == norder)
452 /* This type did not occur yet, add it */
453 order[norder] = sys->molblock[mb].type;
454 /* Renumber the moltype in the topology */
455 norder++;
457 sys->molblock[mb].type = i;
460 /* We still need to reorder the molinfo structs */
461 snew(minew, norder);
462 for (mi = 0; mi < *nmolinfo; mi++)
464 for (i = 0; i < norder; i++)
466 if (order[i] == mi)
468 break;
471 if (i == norder)
473 done_mi(&(*molinfo)[mi]);
475 else
477 minew[i] = (*molinfo)[mi];
480 sfree(*molinfo);
482 *nmolinfo = norder;
483 *molinfo = minew;
486 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
488 int m;
489 gmx_moltype_t *molt;
491 mtop->nmoltype = nmi;
492 snew(mtop->moltype, nmi);
493 for (m = 0; m < nmi; m++)
495 molt = &mtop->moltype[m];
496 molt->name = mi[m].name;
497 molt->atoms = mi[m].atoms;
498 /* ilists are copied later */
499 molt->cgs = mi[m].cgs;
500 molt->excls = mi[m].excls;
504 static void
505 new_status(const char *topfile, const char *topppfile, const char *confin,
506 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
507 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
508 gpp_atomtype_t atype, gmx_mtop_t *sys,
509 int *nmi, t_molinfo **mi, t_molinfo **intermolecular_interactions,
510 t_params plist[],
511 int *comb, double *reppow, real *fudgeQQ,
512 gmx_bool bMorse,
513 warninp_t wi)
515 t_molinfo *molinfo = nullptr;
516 int nmolblock;
517 gmx_molblock_t *molblock, *molbs;
518 int mb, i, nrmols, nmismatch;
519 char buf[STRLEN];
520 gmx_bool bGB = FALSE;
521 char warn_buf[STRLEN];
523 init_mtop(sys);
525 /* Set gmx_boolean for GB */
526 if (ir->implicit_solvent)
528 bGB = TRUE;
531 /* TOPOLOGY processing */
532 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
533 plist, comb, reppow, fudgeQQ,
534 atype, &nrmols, &molinfo, intermolecular_interactions,
536 &nmolblock, &molblock, bGB,
537 wi);
539 sys->nmolblock = 0;
540 snew(sys->molblock, nmolblock);
542 sys->natoms = 0;
543 for (mb = 0; mb < nmolblock; mb++)
545 if (sys->nmolblock > 0 &&
546 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
548 /* Merge consecutive blocks with the same molecule type */
549 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
550 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
552 else if (molblock[mb].nmol > 0)
554 /* Add a new molblock to the topology */
555 molbs = &sys->molblock[sys->nmolblock];
556 *molbs = molblock[mb];
557 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
558 molbs->nposres_xA = 0;
559 molbs->nposres_xB = 0;
560 sys->natoms += molbs->nmol*molbs->natoms_mol;
561 sys->nmolblock++;
564 if (sys->nmolblock == 0)
566 gmx_fatal(FARGS, "No molecules were defined in the system");
569 renumber_moltypes(sys, &nrmols, &molinfo);
571 if (bMorse)
573 convert_harmonics(nrmols, molinfo, atype);
576 if (ir->eDisre == edrNone)
578 i = rm_interactions(F_DISRES, nrmols, molinfo);
579 if (i > 0)
581 set_warning_line(wi, "unknown", -1);
582 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
583 warning_note(wi, warn_buf);
586 if (opts->bOrire == FALSE)
588 i = rm_interactions(F_ORIRES, nrmols, molinfo);
589 if (i > 0)
591 set_warning_line(wi, "unknown", -1);
592 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
593 warning_note(wi, warn_buf);
597 /* Copy structures from msys to sys */
598 molinfo2mtop(nrmols, molinfo, sys);
600 gmx_mtop_finalize(sys);
602 /* COORDINATE file processing */
603 if (bVerbose)
605 fprintf(stderr, "processing coordinates...\n");
608 t_topology *conftop;
609 rvec *x = nullptr;
610 rvec *v = nullptr;
611 snew(conftop, 1);
612 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
613 state->natoms = conftop->atoms.nr;
614 if (state->natoms != sys->natoms)
616 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
617 " does not match topology (%s, %d)",
618 confin, state->natoms, topfile, sys->natoms);
620 /* It would be nice to get rid of the copies below, but we don't know
621 * a priori if the number of atoms in confin matches what we expect.
623 state->flags |= (1 << estX);
624 if (v != NULL)
626 state->flags |= (1 << estV);
628 state_change_natoms(state, state->natoms);
629 for (int i = 0; i < state->natoms; i++)
631 copy_rvec(x[i], state->x[i]);
633 sfree(x);
634 if (v != nullptr)
636 for (int i = 0; i < state->natoms; i++)
638 copy_rvec(v[i], state->v[i]);
640 sfree(v);
642 /* This call fixes the box shape for runs with pressure scaling */
643 set_box_rel(ir, state->box_rel, state->box);
645 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
646 done_top(conftop);
647 sfree(conftop);
649 if (nmismatch)
651 sprintf(buf, "%d non-matching atom name%s\n"
652 "atom names from %s will be used\n"
653 "atom names from %s will be ignored\n",
654 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
655 warning(wi, buf);
658 /* Do more checks, mostly related to constraints */
659 if (bVerbose)
661 fprintf(stderr, "double-checking input for internal consistency...\n");
664 int bHasNormalConstraints = 0 < (nint_ftype(sys, molinfo, F_CONSTR) +
665 nint_ftype(sys, molinfo, F_CONSTRNC));
666 int bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, molinfo, F_SETTLE);
667 double_check(ir, state->box,
668 bHasNormalConstraints,
669 bHasAnyConstraints,
670 wi);
673 if (bGenVel)
675 real *mass;
676 gmx_mtop_atomloop_all_t aloop;
677 const t_atom *atom;
679 snew(mass, state->natoms);
680 aloop = gmx_mtop_atomloop_all_init(sys);
681 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
683 mass[i] = atom->m;
686 if (opts->seed == -1)
688 opts->seed = static_cast<int>(gmx::makeRandomSeed());
689 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
691 state->flags |= (1 << estV);
692 maxwell_speed(opts->tempi, opts->seed, sys, as_rvec_array(state->v.data()));
694 stop_cm(stdout, state->natoms, mass, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()));
695 sfree(mass);
698 *nmi = nrmols;
699 *mi = molinfo;
702 static void copy_state(const char *slog, t_trxframe *fr,
703 gmx_bool bReadVel, t_state *state,
704 double *use_time)
706 int i;
708 if (fr->not_ok & FRAME_NOT_OK)
710 gmx_fatal(FARGS, "Can not start from an incomplete frame");
712 if (!fr->bX)
714 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
715 slog);
718 for (i = 0; i < state->natoms; i++)
720 copy_rvec(fr->x[i], state->x[i]);
722 if (bReadVel)
724 if (!fr->bV)
726 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
728 for (i = 0; i < state->natoms; i++)
730 copy_rvec(fr->v[i], state->v[i]);
733 if (fr->bBox)
735 copy_mat(fr->box, state->box);
738 *use_time = fr->time;
741 static void cont_status(const char *slog, const char *ener,
742 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
743 t_inputrec *ir, t_state *state,
744 gmx_mtop_t *sys,
745 const gmx_output_env_t *oenv)
746 /* If fr_time == -1 read the last frame available which is complete */
748 gmx_bool bReadVel;
749 t_trxframe fr;
750 t_trxstatus *fp;
751 int i;
752 double use_time;
754 bReadVel = (bNeedVel && !bGenVel);
756 fprintf(stderr,
757 "Reading Coordinates%s and Box size from old trajectory\n",
758 bReadVel ? ", Velocities" : "");
759 if (fr_time == -1)
761 fprintf(stderr, "Will read whole trajectory\n");
763 else
765 fprintf(stderr, "Will read till time %g\n", fr_time);
767 if (!bReadVel)
769 if (bGenVel)
771 fprintf(stderr, "Velocities generated: "
772 "ignoring velocities in input trajectory\n");
774 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
776 else
778 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
780 if (!fr.bV)
782 fprintf(stderr,
783 "\n"
784 "WARNING: Did not find a frame with velocities in file %s,\n"
785 " all velocities will be set to zero!\n\n", slog);
786 for (i = 0; i < sys->natoms; i++)
788 clear_rvec(state->v[i]);
790 close_trx(fp);
791 /* Search for a frame without velocities */
792 bReadVel = FALSE;
793 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
797 state->natoms = fr.natoms;
799 if (sys->natoms != state->natoms)
801 gmx_fatal(FARGS, "Number of atoms in Topology "
802 "is not the same as in Trajectory");
804 copy_state(slog, &fr, bReadVel, state, &use_time);
806 /* Find the appropriate frame */
807 while ((fr_time == -1 || fr.time < fr_time) &&
808 read_next_frame(oenv, fp, &fr))
810 copy_state(slog, &fr, bReadVel, state, &use_time);
813 close_trx(fp);
815 /* Set the relative box lengths for preserving the box shape.
816 * Note that this call can lead to differences in the last bit
817 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
819 set_box_rel(ir, state->box_rel, state->box);
821 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
822 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
824 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
826 get_enx_state(ener, use_time, &sys->groups, ir, state);
827 preserve_box_shape(ir, state->box_rel, state->boxv);
831 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
832 const char *fn,
833 int rc_scaling, int ePBC,
834 rvec com,
835 warninp_t wi)
837 gmx_bool *hadAtom;
838 rvec *x, *v, *xp;
839 dvec sum;
840 double totmass;
841 t_topology *top;
842 matrix box, invbox;
843 int natoms, npbcdim = 0;
844 char warn_buf[STRLEN];
845 int a, i, ai, j, k, mb, nat_molb;
846 gmx_molblock_t *molb;
847 t_params *pr, *prfb;
848 t_atom *atom;
850 snew(top, 1);
851 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
852 natoms = top->atoms.nr;
853 done_top(top);
854 sfree(top);
855 if (natoms != mtop->natoms)
857 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);
858 warning(wi, warn_buf);
861 npbcdim = ePBC2npbcdim(ePBC);
862 clear_rvec(com);
863 if (rc_scaling != erscNO)
865 copy_mat(box, invbox);
866 for (j = npbcdim; j < DIM; j++)
868 clear_rvec(invbox[j]);
869 invbox[j][j] = 1;
871 gmx::invertBoxMatrix(invbox, invbox);
874 /* Copy the reference coordinates to mtop */
875 clear_dvec(sum);
876 totmass = 0;
877 a = 0;
878 snew(hadAtom, natoms);
879 for (mb = 0; mb < mtop->nmolblock; mb++)
881 molb = &mtop->molblock[mb];
882 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
883 pr = &(molinfo[molb->type].plist[F_POSRES]);
884 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
885 if (pr->nr > 0 || prfb->nr > 0)
887 atom = mtop->moltype[molb->type].atoms.atom;
888 for (i = 0; (i < pr->nr); i++)
890 ai = pr->param[i].ai();
891 if (ai >= natoms)
893 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
894 ai+1, *molinfo[molb->type].name, fn, natoms);
896 hadAtom[ai] = TRUE;
897 if (rc_scaling == erscCOM)
899 /* Determine the center of mass of the posres reference coordinates */
900 for (j = 0; j < npbcdim; j++)
902 sum[j] += atom[ai].m*x[a+ai][j];
904 totmass += atom[ai].m;
907 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
908 for (i = 0; (i < prfb->nr); i++)
910 ai = prfb->param[i].ai();
911 if (ai >= natoms)
913 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
914 ai+1, *molinfo[molb->type].name, fn, natoms);
916 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
918 /* Determine the center of mass of the posres reference coordinates */
919 for (j = 0; j < npbcdim; j++)
921 sum[j] += atom[ai].m*x[a+ai][j];
923 totmass += atom[ai].m;
926 if (!bTopB)
928 molb->nposres_xA = nat_molb;
929 snew(molb->posres_xA, molb->nposres_xA);
930 for (i = 0; i < nat_molb; i++)
932 copy_rvec(x[a+i], molb->posres_xA[i]);
935 else
937 molb->nposres_xB = nat_molb;
938 snew(molb->posres_xB, molb->nposres_xB);
939 for (i = 0; i < nat_molb; i++)
941 copy_rvec(x[a+i], molb->posres_xB[i]);
945 a += nat_molb;
947 if (rc_scaling == erscCOM)
949 if (totmass == 0)
951 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
953 for (j = 0; j < npbcdim; j++)
955 com[j] = sum[j]/totmass;
957 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]);
960 if (rc_scaling != erscNO)
962 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
964 for (mb = 0; mb < mtop->nmolblock; mb++)
966 molb = &mtop->molblock[mb];
967 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
968 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
970 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
971 for (i = 0; i < nat_molb; i++)
973 for (j = 0; j < npbcdim; j++)
975 if (rc_scaling == erscALL)
977 /* Convert from Cartesian to crystal coordinates */
978 xp[i][j] *= invbox[j][j];
979 for (k = j+1; k < npbcdim; k++)
981 xp[i][j] += invbox[k][j]*xp[i][k];
984 else if (rc_scaling == erscCOM)
986 /* Subtract the center of mass */
987 xp[i][j] -= com[j];
994 if (rc_scaling == erscCOM)
996 /* Convert the COM from Cartesian to crystal coordinates */
997 for (j = 0; j < npbcdim; j++)
999 com[j] *= invbox[j][j];
1000 for (k = j+1; k < npbcdim; k++)
1002 com[j] += invbox[k][j]*com[k];
1008 sfree(x);
1009 sfree(v);
1010 sfree(hadAtom);
1013 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
1014 const char *fnA, const char *fnB,
1015 int rc_scaling, int ePBC,
1016 rvec com, rvec comB,
1017 warninp_t wi)
1019 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1020 /* It is safer to simply read the b-state posres rather than trying
1021 * to be smart and copy the positions.
1023 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1026 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
1027 t_inputrec *ir, warninp_t wi)
1029 int i;
1030 char warn_buf[STRLEN];
1032 if (ir->nwall > 0)
1034 fprintf(stderr, "Searching the wall atom type(s)\n");
1036 for (i = 0; i < ir->nwall; i++)
1038 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1039 if (ir->wall_atomtype[i] == NOTSET)
1041 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1042 warning_error(wi, warn_buf);
1047 static int nrdf_internal(t_atoms *atoms)
1049 int i, nmass, nrdf;
1051 nmass = 0;
1052 for (i = 0; i < atoms->nr; i++)
1054 /* Vsite ptype might not be set here yet, so also check the mass */
1055 if ((atoms->atom[i].ptype == eptAtom ||
1056 atoms->atom[i].ptype == eptNucleus)
1057 && atoms->atom[i].m > 0)
1059 nmass++;
1062 switch (nmass)
1064 case 0: nrdf = 0; break;
1065 case 1: nrdf = 0; break;
1066 case 2: nrdf = 1; break;
1067 default: nrdf = nmass*3 - 6; break;
1070 return nrdf;
1073 void
1074 spline1d( double dx,
1075 double * y,
1076 int n,
1077 double * u,
1078 double * y2 )
1080 int i;
1081 double p, q;
1083 y2[0] = 0.0;
1084 u[0] = 0.0;
1086 for (i = 1; i < n-1; i++)
1088 p = 0.5*y2[i-1]+2.0;
1089 y2[i] = -0.5/p;
1090 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1091 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1094 y2[n-1] = 0.0;
1096 for (i = n-2; i >= 0; i--)
1098 y2[i] = y2[i]*y2[i+1]+u[i];
1103 void
1104 interpolate1d( double xmin,
1105 double dx,
1106 double * ya,
1107 double * y2a,
1108 double x,
1109 double * y,
1110 double * y1)
1112 int ix;
1113 double a, b;
1115 ix = static_cast<int>((x-xmin)/dx);
1117 a = (xmin+(ix+1)*dx-x)/dx;
1118 b = (x-xmin-ix*dx)/dx;
1120 *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;
1121 *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];
1125 void
1126 setup_cmap (int grid_spacing,
1127 int nc,
1128 real * grid,
1129 gmx_cmap_t * cmap_grid)
1131 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1133 int i, j, k, ii, jj, kk, idx;
1134 int offset;
1135 double dx, xmin, v, v1, v2, v12;
1136 double phi, psi;
1138 snew(tmp_u, 2*grid_spacing);
1139 snew(tmp_u2, 2*grid_spacing);
1140 snew(tmp_yy, 2*grid_spacing);
1141 snew(tmp_y1, 2*grid_spacing);
1142 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1143 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1145 dx = 360.0/grid_spacing;
1146 xmin = -180.0-dx*grid_spacing/2;
1148 for (kk = 0; kk < nc; kk++)
1150 /* Compute an offset depending on which cmap we are using
1151 * Offset will be the map number multiplied with the
1152 * grid_spacing * grid_spacing * 2
1154 offset = kk * grid_spacing * grid_spacing * 2;
1156 for (i = 0; i < 2*grid_spacing; i++)
1158 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1160 for (j = 0; j < 2*grid_spacing; j++)
1162 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1163 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1167 for (i = 0; i < 2*grid_spacing; i++)
1169 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1172 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1174 ii = i-grid_spacing/2;
1175 phi = ii*dx-180.0;
1177 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1179 jj = j-grid_spacing/2;
1180 psi = jj*dx-180.0;
1182 for (k = 0; k < 2*grid_spacing; k++)
1184 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1185 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1188 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1189 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1190 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1191 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1193 idx = ii*grid_spacing+jj;
1194 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1195 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1196 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1197 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1203 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1205 int i, nelem;
1207 cmap_grid->ngrid = ngrid;
1208 cmap_grid->grid_spacing = grid_spacing;
1209 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1211 snew(cmap_grid->cmapdata, ngrid);
1213 for (i = 0; i < cmap_grid->ngrid; i++)
1215 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1220 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1222 int count, count_mol, i, mb;
1223 gmx_molblock_t *molb;
1224 t_params *plist;
1225 char buf[STRLEN];
1227 count = 0;
1228 for (mb = 0; mb < mtop->nmolblock; mb++)
1230 count_mol = 0;
1231 molb = &mtop->molblock[mb];
1232 plist = mi[molb->type].plist;
1234 for (i = 0; i < F_NRE; i++)
1236 if (i == F_SETTLE)
1238 count_mol += 3*plist[i].nr;
1240 else if (interaction_function[i].flags & IF_CONSTRAINT)
1242 count_mol += plist[i].nr;
1246 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1248 sprintf(buf,
1249 "Molecule type '%s' has %d constraints.\n"
1250 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1251 *mi[molb->type].name, count_mol,
1252 nrdf_internal(&mi[molb->type].atoms));
1253 warning(wi, buf);
1255 count += molb->nmol*count_mol;
1258 return count;
1261 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1263 int i, nmiss, natoms, mt;
1264 real q;
1265 const t_atoms *atoms;
1267 nmiss = 0;
1268 for (mt = 0; mt < sys->nmoltype; mt++)
1270 atoms = &sys->moltype[mt].atoms;
1271 natoms = atoms->nr;
1273 for (i = 0; i < natoms; i++)
1275 q = atoms->atom[i].q;
1276 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1277 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1278 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1279 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1280 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1281 q != 0)
1283 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1284 get_atomtype_name(atoms->atom[i].type, atype), q);
1285 nmiss++;
1290 if (nmiss > 0)
1292 gmx_fatal(FARGS, "Can't do GB electrostatics; the implicit_genborn_params section of the forcefield has parameters with value zero for %d atomtypes that occur as charged atoms.", nmiss);
1297 static void check_gbsa_params(gpp_atomtype_t atype)
1299 int nmiss, i;
1301 /* If we are doing GBSA, check that we got the parameters we need
1302 * This checking is to see if there are GBSA paratmeters for all
1303 * atoms in the force field. To go around this for testing purposes
1304 * comment out the nerror++ counter temporarily
1306 nmiss = 0;
1307 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1309 if (get_atomtype_radius(i, atype) < 0 ||
1310 get_atomtype_vol(i, atype) < 0 ||
1311 get_atomtype_surftens(i, atype) < 0 ||
1312 get_atomtype_gb_radius(i, atype) < 0 ||
1313 get_atomtype_S_hct(i, atype) < 0)
1315 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1316 get_atomtype_name(i, atype));
1317 nmiss++;
1321 if (nmiss > 0)
1323 gmx_fatal(FARGS, "Can't do GB electrostatics; the implicit_genborn_params section of the forcefield is missing parameters for %d atomtypes or they might be negative.", nmiss);
1328 static real calc_temp(const gmx_mtop_t *mtop,
1329 const t_inputrec *ir,
1330 rvec *v)
1332 gmx_mtop_atomloop_all_t aloop;
1333 const t_atom *atom;
1334 int a;
1336 double sum_mv2 = 0;
1337 aloop = gmx_mtop_atomloop_all_init(mtop);
1338 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1340 sum_mv2 += atom->m*norm2(v[a]);
1343 double nrdf = 0;
1344 for (int g = 0; g < ir->opts.ngtc; g++)
1346 nrdf += ir->opts.nrdf[g];
1349 return sum_mv2/(nrdf*BOLTZ);
1352 static real get_max_reference_temp(const t_inputrec *ir,
1353 warninp_t wi)
1355 real ref_t;
1356 int i;
1357 gmx_bool bNoCoupl;
1359 ref_t = 0;
1360 bNoCoupl = FALSE;
1361 for (i = 0; i < ir->opts.ngtc; i++)
1363 if (ir->opts.tau_t[i] < 0)
1365 bNoCoupl = TRUE;
1367 else
1369 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1373 if (bNoCoupl)
1375 char buf[STRLEN];
1377 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.",
1378 ref_t);
1379 warning(wi, buf);
1382 return ref_t;
1385 /* Checks if there are unbound atoms in moleculetype molt.
1386 * Prints a note for each unbound atoms and a warning if any is present.
1388 static void checkForUnboundAtoms(const gmx_moltype_t *molt,
1389 gmx_bool bVerbose,
1390 warninp_t wi)
1392 const t_atoms *atoms = &molt->atoms;
1394 if (atoms->nr == 1)
1396 /* Only one atom, there can't be unbound atoms */
1397 return;
1400 std::vector<int> count(atoms->nr, 0);
1402 for (int ftype = 0; ftype < F_NRE; ftype++)
1404 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS) ||
1405 (interaction_function[ftype].flags & IF_CONSTRAINT) ||
1406 ftype == F_SETTLE)
1408 const t_ilist *il = &molt->ilist[ftype];
1409 int nral = NRAL(ftype);
1411 for (int i = 0; i < il->nr; i += 1 + nral)
1413 for (int j = 0; j < nral; j++)
1415 count[il->iatoms[i + 1 + j]]++;
1421 int numDanglingAtoms = 0;
1422 for (int a = 0; a < atoms->nr; a++)
1424 if (atoms->atom[a].ptype != eptVSite &&
1425 count[a] == 0)
1427 if (bVerbose)
1429 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",
1430 a + 1, *atoms->atomname[a], *molt->name);
1432 numDanglingAtoms++;
1436 if (numDanglingAtoms > 0)
1438 char buf[STRLEN];
1439 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.",
1440 *molt->name, numDanglingAtoms);
1441 warning_note(wi, buf);
1445 /* Checks all moleculetypes for unbound atoms */
1446 static void checkForUnboundAtoms(const gmx_mtop_t *mtop,
1447 gmx_bool bVerbose,
1448 warninp_t wi)
1450 for (int mt = 0; mt < mtop->nmoltype; mt++)
1452 checkForUnboundAtoms(&mtop->moltype[mt], bVerbose, wi);
1456 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1458 * The specific decoupled modes this routine check for are angle modes
1459 * where the two bonds are constrained and the atoms a both ends are only
1460 * involved in a single constraint; the mass of the two atoms needs to
1461 * differ by more than \p massFactorThreshold.
1463 static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
1464 const t_iparams *iparams,
1465 real massFactorThreshold)
1467 if (molt->ilist[F_CONSTR].nr == 0 && molt->ilist[F_CONSTRNC].nr == 0)
1469 return false;
1472 const t_atom * atom = molt->atoms.atom;
1474 int numFlexibleConstraints;
1475 t_blocka atomToConstraints = make_at2con(0, molt->atoms.nr,
1476 molt->ilist, iparams,
1477 FALSE,
1478 &numFlexibleConstraints);
1480 bool haveDecoupledMode = false;
1481 for (int ftype = 0; ftype < F_NRE; ftype++)
1483 if (interaction_function[ftype].flags & IF_ATYPE)
1485 const int nral = NRAL(ftype);
1486 const t_ilist *il = &molt->ilist[ftype];
1487 for (int i = 0; i < il->nr; i += 1 + nral)
1489 /* Here we check for the mass difference between the atoms
1490 * at both ends of the angle, that the atoms at the ends have
1491 * 1 contraint and the atom in the middle at least 3; we check
1492 * that the 3 atoms are linked by constraints below.
1493 * We check for at least three constraints for the middle atom,
1494 * since with only the two bonds in the angle, we have 3-atom
1495 * molecule, which has much more thermal exhange in this single
1496 * angle mode than molecules with more atoms.
1497 * Note that this check also catches molecules where more atoms
1498 * are connected to one or more atoms in the angle, but by
1499 * bond potentials instead of angles. But such cases will not
1500 * occur in "normal" molecules and it doesn't hurt running
1501 * those with higher accuracy settings as well.
1503 int a0 = il->iatoms[1 + i];
1504 int a1 = il->iatoms[1 + i + 1];
1505 int a2 = il->iatoms[1 + i + 2];
1506 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
1507 atom[a2].m > atom[a0].m*massFactorThreshold) &&
1508 atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
1509 atomToConstraints.index[a2 + 1] - atomToConstraints.index[a2] == 1 &&
1510 atomToConstraints.index[a1 + 1] - atomToConstraints.index[a1] >= 3)
1512 int constraint0 = atomToConstraints.a[atomToConstraints.index[a0]];
1513 int constraint2 = atomToConstraints.a[atomToConstraints.index[a2]];
1515 bool foundAtom0 = false;
1516 bool foundAtom2 = false;
1517 for (int conIndex = atomToConstraints.index[a1]; conIndex < atomToConstraints.index[a1 + 1]; conIndex++)
1519 if (atomToConstraints.a[conIndex] == constraint0)
1521 foundAtom0 = true;
1523 if (atomToConstraints.a[conIndex] == constraint2)
1525 foundAtom2 = true;
1528 if (foundAtom0 && foundAtom2)
1530 haveDecoupledMode = true;
1537 done_blocka(&atomToConstraints);
1539 return haveDecoupledMode;
1542 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1544 * When decoupled modes are present and the accuray in insufficient,
1545 * this routine issues a warning if the accuracy is insufficient.
1547 static void checkDecoupledModeAccuracy(const gmx_mtop_t *mtop,
1548 const t_inputrec *ir,
1549 warninp_t wi)
1551 /* We only have issues with decoupled modes with normal MD.
1552 * With stochastic dynamics equipartitioning is enforced strongly.
1554 if (!EI_MD(ir->eI))
1556 return;
1559 /* When atoms of very different mass are involved in an angle potential
1560 * and both bonds in the angle are constrained, the dynamic modes in such
1561 * angles have very different periods and significant energy exchange
1562 * takes several nanoseconds. Thus even a small amount of error in
1563 * different algorithms can lead to violation of equipartitioning.
1564 * The parameters below are mainly based on an all-atom chloroform model
1565 * with all bonds constrained. Equipartitioning is off by more than 1%
1566 * (need to run 5-10 ns) when the difference in mass between H and Cl
1567 * is more than a factor 13 and the accuracy is less than the thresholds
1568 * given below. This has been verified on some other molecules.
1570 * Note that the buffer and shake parameters have unit length and
1571 * energy/time, respectively, so they will "only" work correctly
1572 * for atomistic force fields using MD units.
1574 const real massFactorThreshold = 13.0;
1575 const real bufferToleranceThreshold = 1e-4;
1576 const int lincsIterationThreshold = 2;
1577 const int lincsOrderThreshold = 4;
1578 const real shakeToleranceThreshold = 0.005*ir->delta_t;
1580 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1581 bool shakeWithSufficientTolerance = (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1*shakeToleranceThreshold);
1582 if (ir->cutoff_scheme == ecutsVERLET &&
1583 ir->verletbuf_tol <= 1.1*bufferToleranceThreshold &&
1584 (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1586 return;
1589 bool haveDecoupledMode = false;
1590 for (int mt = 0; mt < mtop->nmoltype; mt++)
1592 if (haveDecoupledModeInMol(&mtop->moltype[mt], mtop->ffparams.iparams,
1593 massFactorThreshold))
1595 haveDecoupledMode = true;
1599 if (haveDecoupledMode)
1601 char modeMessage[STRLEN];
1602 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.",
1603 massFactorThreshold);
1604 char buf[STRLEN];
1605 if (ir->cutoff_scheme == ecutsVERLET)
1607 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",
1608 modeMessage,
1609 ei_names[eiSD1],
1610 bufferToleranceThreshold,
1611 lincsIterationThreshold, lincsOrderThreshold,
1612 shakeToleranceThreshold);
1614 else
1616 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.",
1617 modeMessage,
1618 ecutscheme_names[ecutsVERLET]);
1620 warning(wi, buf);
1624 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1625 t_inputrec *ir,
1626 real buffer_temp,
1627 matrix box,
1628 warninp_t wi)
1630 verletbuf_list_setup_t ls;
1631 real rlist_1x1;
1632 int n_nonlin_vsite;
1633 char warn_buf[STRLEN];
1635 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1637 /* Calculate the buffer size for simple atom vs atoms list */
1638 ls.cluster_size_i = 1;
1639 ls.cluster_size_j = 1;
1640 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1641 &ls, &n_nonlin_vsite, &rlist_1x1);
1643 /* Set the pair-list buffer size in ir */
1644 verletbuf_get_list_setup(FALSE, FALSE, &ls);
1645 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1646 &ls, &n_nonlin_vsite, &ir->rlist);
1648 if (n_nonlin_vsite > 0)
1650 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);
1651 warning_note(wi, warn_buf);
1654 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1655 1, 1, rlist_1x1, rlist_1x1-std::max(ir->rvdw, ir->rcoulomb));
1657 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1658 ls.cluster_size_i, ls.cluster_size_j,
1659 ir->rlist, ir->rlist-std::max(ir->rvdw, ir->rcoulomb));
1661 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1663 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1665 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)));
1669 int gmx_grompp(int argc, char *argv[])
1671 static const char *desc[] = {
1672 "[THISMODULE] (the gromacs preprocessor)",
1673 "reads a molecular topology file, checks the validity of the",
1674 "file, expands the topology from a molecular description to an atomic",
1675 "description. The topology file contains information about",
1676 "molecule types and the number of molecules, the preprocessor",
1677 "copies each molecule as needed. ",
1678 "There is no limitation on the number of molecule types. ",
1679 "Bonds and bond-angles can be converted into constraints, separately",
1680 "for hydrogens and heavy atoms.",
1681 "Then a coordinate file is read and velocities can be generated",
1682 "from a Maxwellian distribution if requested.",
1683 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1684 "(eg. number of MD steps, time step, cut-off), and others such as",
1685 "NEMD parameters, which are corrected so that the net acceleration",
1686 "is zero.",
1687 "Eventually a binary file is produced that can serve as the sole input",
1688 "file for the MD program.[PAR]",
1690 "[THISMODULE] uses the atom names from the topology file. The atom names",
1691 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1692 "warnings when they do not match the atom names in the topology.",
1693 "Note that the atom names are irrelevant for the simulation as",
1694 "only the atom types are used for generating interaction parameters.[PAR]",
1696 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1697 "etc. The preprocessor supports the following keywords::",
1699 " #ifdef VARIABLE",
1700 " #ifndef VARIABLE",
1701 " #else",
1702 " #endif",
1703 " #define VARIABLE",
1704 " #undef VARIABLE",
1705 " #include \"filename\"",
1706 " #include <filename>",
1708 "The functioning of these statements in your topology may be modulated by",
1709 "using the following two flags in your [REF].mdp[ref] file::",
1711 " define = -DVARIABLE1 -DVARIABLE2",
1712 " include = -I/home/john/doe",
1714 "For further information a C-programming textbook may help you out.",
1715 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1716 "topology file written out so that you can verify its contents.[PAR]",
1718 "When using position restraints, a file with restraint coordinates",
1719 "must be supplied with [TT]-r[tt] (can be the same file as supplied",
1720 "for [TT]-c[tt]). For free energy calculations, separate reference",
1721 "coordinates for the B topology can be supplied with [TT]-rb[tt],",
1722 "otherwise they will be equal to those of the A topology.[PAR]",
1724 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1725 "The last frame with coordinates and velocities will be read,",
1726 "unless the [TT]-time[tt] option is used. Only if this information",
1727 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1728 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1729 "in your [REF].mdp[ref] file. An energy file can be supplied with",
1730 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1731 "variables.[PAR]",
1733 "[THISMODULE] can be used to restart simulations (preserving",
1734 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1735 "However, for simply changing the number of run steps to extend",
1736 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1737 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1738 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1739 "like output frequency, then supplying the checkpoint file to",
1740 "[THISMODULE] with [TT]-t[tt] along with a new [REF].mdp[ref] file",
1741 "with [TT]-f[tt] is the recommended procedure. Actually preserving",
1742 "the ensemble (if possible) still requires passing the checkpoint",
1743 "file to [gmx-mdrun] [TT]-cpi[tt].[PAR]",
1745 "By default, all bonded interactions which have constant energy due to",
1746 "virtual site constructions will be removed. If this constant energy is",
1747 "not zero, this will result in a shift in the total energy. All bonded",
1748 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1749 "all constraints for distances which will be constant anyway because",
1750 "of virtual site constructions will be removed. If any constraints remain",
1751 "which involve virtual sites, a fatal error will result.[PAR]"
1753 "To verify your run input file, please take note of all warnings",
1754 "on the screen, and correct where necessary. Do also look at the contents",
1755 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1756 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1757 "with the [TT]-debug[tt] option which will give you more information",
1758 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1759 "can see the contents of the run input file with the [gmx-dump]",
1760 "program. [gmx-check] can be used to compare the contents of two",
1761 "run input files.[PAR]"
1763 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1764 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1765 "harmless, but usually they are not. The user is advised to carefully",
1766 "interpret the output messages before attempting to bypass them with",
1767 "this option."
1769 t_gromppopts *opts;
1770 gmx_mtop_t *sys;
1771 int nmi;
1772 t_molinfo *mi, *intermolecular_interactions;
1773 gpp_atomtype_t atype;
1774 int nvsite, comb, mt;
1775 t_params *plist;
1776 matrix box;
1777 real fudgeQQ;
1778 double reppow;
1779 const char *mdparin;
1780 int ntype;
1781 gmx_bool bNeedVel, bGenVel;
1782 gmx_bool have_atomnumber;
1783 gmx_output_env_t *oenv;
1784 gmx_bool bVerbose = FALSE;
1785 warninp_t wi;
1786 char warn_buf[STRLEN];
1788 t_filenm fnm[] = {
1789 { efMDP, nullptr, nullptr, ffREAD },
1790 { efMDP, "-po", "mdout", ffWRITE },
1791 { efSTX, "-c", nullptr, ffREAD },
1792 { efSTX, "-r", "restraint", ffOPTRD },
1793 { efSTX, "-rb", "restraint", ffOPTRD },
1794 { efNDX, nullptr, nullptr, ffOPTRD },
1795 { efTOP, nullptr, nullptr, ffREAD },
1796 { efTOP, "-pp", "processed", ffOPTWR },
1797 { efTPR, "-o", nullptr, ffWRITE },
1798 { efTRN, "-t", nullptr, ffOPTRD },
1799 { efEDR, "-e", nullptr, ffOPTRD },
1800 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1801 { efGRO, "-imd", "imdgroup", ffOPTWR },
1802 { efTRN, "-ref", "rotref", ffOPTRW }
1804 #define NFILE asize(fnm)
1806 /* Command line options */
1807 static gmx_bool bRenum = TRUE;
1808 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1809 static int i, maxwarn = 0;
1810 static real fr_time = -1;
1811 t_pargs pa[] = {
1812 { "-v", FALSE, etBOOL, {&bVerbose},
1813 "Be loud and noisy" },
1814 { "-time", FALSE, etREAL, {&fr_time},
1815 "Take frame at or first after this time." },
1816 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1817 "Remove constant bonded interactions with virtual sites" },
1818 { "-maxwarn", FALSE, etINT, {&maxwarn},
1819 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1820 { "-zero", FALSE, etBOOL, {&bZero},
1821 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1822 { "-renum", FALSE, etBOOL, {&bRenum},
1823 "Renumber atomtypes and minimize number of atomtypes" }
1826 /* Parse the command line */
1827 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1828 asize(desc), desc, 0, nullptr, &oenv))
1830 return 0;
1833 /* Initiate some variables */
1834 gmx::MDModules mdModules;
1835 t_inputrec irInstance;
1836 t_inputrec *ir = &irInstance;
1837 snew(opts, 1);
1838 snew(opts->include, STRLEN);
1839 snew(opts->define, STRLEN);
1841 wi = init_warning(TRUE, maxwarn);
1843 /* PARAMETER file processing */
1844 mdparin = opt2fn("-f", NFILE, fnm);
1845 set_warning_line(wi, mdparin, -1);
1848 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, ir, opts, WriteMdpHeader::yes, wi);
1850 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1852 if (bVerbose)
1854 fprintf(stderr, "checking input for internal consistency...\n");
1856 check_ir(mdparin, ir, opts, wi);
1858 if (ir->ld_seed == -1)
1860 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1861 fprintf(stderr, "Setting the LD random seed to %" GMX_PRId64 "\n", ir->ld_seed);
1864 if (ir->expandedvals->lmc_seed == -1)
1866 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1867 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1870 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1871 bGenVel = (bNeedVel && opts->bGenVel);
1872 if (bGenVel && ir->bContinuation)
1874 sprintf(warn_buf,
1875 "Generating velocities is inconsistent with attempting "
1876 "to continue a previous run. Choose only one of "
1877 "gen-vel = yes and continuation = yes.");
1878 warning_error(wi, warn_buf);
1881 snew(plist, F_NRE);
1882 init_plist(plist);
1883 snew(sys, 1);
1884 atype = init_atomtype();
1885 if (debug)
1887 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1890 const char *fn = ftp2fn(efTOP, NFILE, fnm);
1891 if (!gmx_fexist(fn))
1893 gmx_fatal(FARGS, "%s does not exist", fn);
1896 t_state state;
1897 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1898 opts, ir, bZero, bGenVel, bVerbose, &state,
1899 atype, sys, &nmi, &mi, &intermolecular_interactions,
1900 plist, &comb, &reppow, &fudgeQQ,
1901 opts->bMorse,
1902 wi);
1904 if (debug)
1906 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1909 nvsite = 0;
1910 /* set parameters for virtual site construction (not for vsiten) */
1911 for (mt = 0; mt < sys->nmoltype; mt++)
1913 nvsite +=
1914 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1916 /* now throw away all obsolete bonds, angles and dihedrals: */
1917 /* note: constraints are ALWAYS removed */
1918 if (nvsite)
1920 for (mt = 0; mt < sys->nmoltype; mt++)
1922 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1926 if (ir->cutoff_scheme == ecutsVERLET)
1928 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1929 ecutscheme_names[ir->cutoff_scheme]);
1931 /* Remove all charge groups */
1932 gmx_mtop_remove_chargegroups(sys);
1935 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1937 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1939 sprintf(warn_buf, "Can not do %s with %s, use %s",
1940 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1941 warning_error(wi, warn_buf);
1943 if (ir->bPeriodicMols)
1945 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1946 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1947 warning_error(wi, warn_buf);
1951 if (EI_SD (ir->eI) && ir->etc != etcNO)
1953 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1956 /* If we are doing QM/MM, check that we got the atom numbers */
1957 have_atomnumber = TRUE;
1958 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1960 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1962 if (!have_atomnumber && ir->bQMMM)
1964 warning_error(wi,
1965 "\n"
1966 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1967 "field you are using does not contain atom numbers fields. This is an\n"
1968 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1969 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1970 "an integer just before the mass column in ffXXXnb.itp.\n"
1971 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1974 /* Check for errors in the input now, since they might cause problems
1975 * during processing further down.
1977 check_warning_error(wi, FARGS);
1979 if (nint_ftype(sys, mi, F_POSRES) > 0 ||
1980 nint_ftype(sys, mi, F_FBPOSRES) > 0)
1982 const char *fn = opt2fn("-r", NFILE, fnm);
1983 const char *fnB;
1985 if (opt2bSet("-rb", NFILE, fnm))
1987 fnB = opt2fn("-rb", NFILE, fnm);
1989 else
1991 fnB = fn;
1994 if (bVerbose)
1996 fprintf(stderr, "Reading position restraint coords from %s", fn);
1997 if (strcmp(fn, fnB) == 0)
1999 fprintf(stderr, "\n");
2001 else
2003 fprintf(stderr, " and %s\n", fnB);
2006 gen_posres(sys, mi, fn, fnB,
2007 ir->refcoord_scaling, ir->ePBC,
2008 ir->posres_com, ir->posres_comB,
2009 wi);
2012 /* If we are using CMAP, setup the pre-interpolation grid */
2013 if (plist[F_CMAP].ncmap > 0)
2015 init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
2016 setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
2019 set_wall_atomtype(atype, opts, ir, wi);
2020 if (bRenum)
2022 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
2023 get_atomtype_ntypes(atype);
2026 if (ir->implicit_solvent != eisNO)
2028 /* Now we have renumbered the atom types, we can check the GBSA params */
2029 check_gbsa_params(atype);
2031 /* Check that all atoms that have charge and/or LJ-parameters also have
2032 * sensible GB-parameters
2034 check_gbsa_params_charged(sys, atype);
2037 /* PELA: Copy the atomtype data to the topology atomtype list */
2038 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
2040 if (debug)
2042 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
2045 if (bVerbose)
2047 fprintf(stderr, "converting bonded parameters...\n");
2050 ntype = get_atomtype_ntypes(atype);
2051 convert_params(ntype, plist, mi, intermolecular_interactions,
2052 comb, reppow, fudgeQQ, sys);
2054 if (debug)
2056 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
2059 /* set ptype to VSite for virtual sites */
2060 for (mt = 0; mt < sys->nmoltype; mt++)
2062 set_vsites_ptype(FALSE, &sys->moltype[mt]);
2064 if (debug)
2066 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
2068 /* Check velocity for virtual sites and shells */
2069 if (bGenVel)
2071 check_vel(sys, as_rvec_array(state.v.data()));
2074 /* check for shells and inpurecs */
2075 check_shells_inputrec(sys, ir, wi);
2077 /* check masses */
2078 check_mol(sys, wi);
2080 checkForUnboundAtoms(sys, bVerbose, wi);
2082 for (i = 0; i < sys->nmoltype; i++)
2084 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
2087 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
2089 check_bonds_timestep(sys, ir->delta_t, wi);
2092 checkDecoupledModeAccuracy(sys, ir, wi);
2094 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
2096 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.");
2099 check_warning_error(wi, FARGS);
2101 if (bVerbose)
2103 fprintf(stderr, "initialising group options...\n");
2105 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
2106 sys, bVerbose, ir,
2107 wi);
2109 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0)
2111 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
2113 real buffer_temp;
2115 if (EI_MD(ir->eI) && ir->etc == etcNO)
2117 if (bGenVel)
2119 buffer_temp = opts->tempi;
2121 else
2123 buffer_temp = calc_temp(sys, ir, as_rvec_array(state.v.data()));
2125 if (buffer_temp > 0)
2127 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
2128 warning_note(wi, warn_buf);
2130 else
2132 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
2133 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
2134 warning_note(wi, warn_buf);
2137 else
2139 buffer_temp = get_max_reference_temp(ir, wi);
2142 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
2144 /* NVE with initial T=0: we add a fixed ratio to rlist.
2145 * Since we don't actually use verletbuf_tol, we set it to -1
2146 * so it can't be misused later.
2148 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
2149 ir->verletbuf_tol = -1;
2151 else
2153 /* We warn for NVE simulations with a drift tolerance that
2154 * might result in a 1(.1)% drift over the total run-time.
2155 * Note that we can't warn when nsteps=0, since we don't
2156 * know how many steps the user intends to run.
2158 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
2159 ir->nsteps > 0)
2161 const real driftTolerance = 0.01;
2162 /* We use 2 DOF per atom = 2kT pot+kin energy,
2163 * to be on the safe side with constraints.
2165 const real totalEnergyDriftPerAtomPerPicosecond = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
2167 if (ir->verletbuf_tol > 1.1*driftTolerance*totalEnergyDriftPerAtomPerPicosecond)
2169 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.",
2170 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
2171 (int)(ir->verletbuf_tol/totalEnergyDriftPerAtomPerPicosecond*100 + 0.5),
2172 (int)(100*driftTolerance + 0.5),
2173 driftTolerance*totalEnergyDriftPerAtomPerPicosecond);
2174 warning_note(wi, warn_buf);
2178 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
2183 /* Init the temperature coupling state */
2184 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
2186 if (bVerbose)
2188 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
2190 check_eg_vs_cg(sys);
2192 if (debug)
2194 pr_symtab(debug, 0, "After index", &sys->symtab);
2197 triple_check(mdparin, ir, sys, wi);
2198 close_symtab(&sys->symtab);
2199 if (debug)
2201 pr_symtab(debug, 0, "After close", &sys->symtab);
2204 /* make exclusions between QM atoms */
2205 if (ir->bQMMM)
2207 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
2209 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
2211 else
2213 generate_qmexcl(sys, ir, wi);
2217 if (ftp2bSet(efTRN, NFILE, fnm))
2219 if (bVerbose)
2221 fprintf(stderr, "getting data from old trajectory ...\n");
2223 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
2224 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
2227 if (ir->ePBC == epbcXY && ir->nwall != 2)
2229 clear_rvec(state.box[ZZ]);
2232 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
2234 set_warning_line(wi, mdparin, -1);
2235 check_chargegroup_radii(sys, ir, as_rvec_array(state.x.data()), wi);
2238 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
2240 /* Calculate the optimal grid dimensions */
2241 copy_mat(state.box, box);
2242 if (ir->ePBC == epbcXY && ir->nwall == 2)
2244 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
2246 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
2248 /* Mark fourier_spacing as not used */
2249 ir->fourier_spacing = 0;
2251 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
2253 set_warning_line(wi, mdparin, -1);
2254 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
2256 const int minGridSize = minimalPmeGridSize(ir->pme_order);
2257 calcFftGrid(stdout, box, ir->fourier_spacing, minGridSize,
2258 &(ir->nkx), &(ir->nky), &(ir->nkz));
2259 if (ir->nkx < minGridSize ||
2260 ir->nky < minGridSize ||
2261 ir->nkz < minGridSize)
2263 warning_error(wi, "The PME grid size should be >= 2*(pme-order - 1); either manually increase the grid size or decrease pme-order");
2267 /* MRS: eventually figure out better logic for initializing the fep
2268 values that makes declaring the lambda and declaring the state not
2269 potentially conflict if not handled correctly. */
2270 if (ir->efep != efepNO)
2272 state.fep_state = ir->fepvals->init_fep_state;
2273 for (i = 0; i < efptNR; i++)
2275 /* init_lambda trumps state definitions*/
2276 if (ir->fepvals->init_lambda >= 0)
2278 state.lambda[i] = ir->fepvals->init_lambda;
2280 else
2282 if (ir->fepvals->all_lambda[i] == nullptr)
2284 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2286 else
2288 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2294 struct pull_t *pull = nullptr;
2296 if (ir->bPull)
2298 pull = set_pull_init(ir, sys, as_rvec_array(state.x.data()), state.box, state.lambda[efptMASS], oenv);
2301 /* Modules that supply external potential for pull coordinates
2302 * should register those potentials here. finish_pull() will check
2303 * that providers have been registerd for all external potentials.
2306 if (ir->bPull)
2308 finish_pull(pull);
2311 if (ir->bRot)
2313 set_reference_positions(ir->rot, as_rvec_array(state.x.data()), state.box,
2314 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2315 wi);
2318 /* reset_multinr(sys); */
2320 if (EEL_PME(ir->coulombtype))
2322 float ratio = pme_load_estimate(sys, ir, state.box);
2323 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2324 /* With free energy we might need to do PME both for the A and B state
2325 * charges. This will double the cost, but the optimal performance will
2326 * then probably be at a slightly larger cut-off and grid spacing.
2328 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2329 (ir->efep != efepNO && ratio > 2.0/3.0))
2331 warning_note(wi,
2332 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2333 "and for highly parallel simulations between 0.25 and 0.33,\n"
2334 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2335 if (ir->efep != efepNO)
2337 warning_note(wi,
2338 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2344 char warn_buf[STRLEN];
2345 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2346 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2347 if (cio > 2000)
2349 set_warning_line(wi, mdparin, -1);
2350 warning_note(wi, warn_buf);
2352 else
2354 printf("%s\n", warn_buf);
2358 if (bVerbose)
2360 fprintf(stderr, "writing run input file...\n");
2363 done_warning(wi, FARGS);
2364 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, sys);
2366 /* Output IMD group, if bIMD is TRUE */
2367 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2369 done_atomtype(atype);
2370 done_mtop(sys);
2371 done_inputrec_strings();
2373 return 0;