Improved state_change_natoms and used it more
[gromacs.git] / src / gromacs / gmxpreprocess / grompp.cpp
blob17f3a443095d77a58ff24e7ce471fbbe3b7c5993
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/fatalerror.h"
97 #include "gromacs/utility/futil.h"
98 #include "gromacs/utility/gmxassert.h"
99 #include "gromacs/utility/smalloc.h"
100 #include "gromacs/utility/snprintf.h"
102 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
104 int i, n;
106 n = 0;
107 /* For all the molecule types */
108 for (i = 0; i < nrmols; i++)
110 n += mols[i].plist[ifunc].nr;
111 mols[i].plist[ifunc].nr = 0;
113 return n;
116 static int check_atom_names(const char *fn1, const char *fn2,
117 gmx_mtop_t *mtop, const t_atoms *at)
119 int mb, m, i, j, nmismatch;
120 t_atoms *tat;
121 #define MAXMISMATCH 20
123 if (mtop->natoms != at->nr)
125 gmx_incons("comparing atom names");
128 nmismatch = 0;
129 i = 0;
130 for (mb = 0; mb < mtop->nmolblock; mb++)
132 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
133 for (m = 0; m < mtop->molblock[mb].nmol; m++)
135 for (j = 0; j < tat->nr; j++)
137 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
139 if (nmismatch < MAXMISMATCH)
141 fprintf(stderr,
142 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
143 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
145 else if (nmismatch == MAXMISMATCH)
147 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
149 nmismatch++;
151 i++;
156 return nmismatch;
159 static void check_eg_vs_cg(gmx_mtop_t *mtop)
161 int astart, mb, m, cg, j, firstj;
162 unsigned char firsteg, eg;
163 gmx_moltype_t *molt;
165 /* Go through all the charge groups and make sure all their
166 * atoms are in the same energy group.
169 astart = 0;
170 for (mb = 0; mb < mtop->nmolblock; mb++)
172 molt = &mtop->moltype[mtop->molblock[mb].type];
173 for (m = 0; m < mtop->molblock[mb].nmol; m++)
175 for (cg = 0; cg < molt->cgs.nr; cg++)
177 /* Get the energy group of the first atom in this charge group */
178 firstj = astart + molt->cgs.index[cg];
179 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
180 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
182 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
183 if (eg != firsteg)
185 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
186 firstj+1, astart+j+1, cg+1, *molt->name);
190 astart += molt->atoms.nr;
195 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
197 int maxsize, cg;
198 char warn_buf[STRLEN];
200 maxsize = 0;
201 for (cg = 0; cg < cgs->nr; cg++)
203 maxsize = std::max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
206 if (maxsize > MAX_CHARGEGROUP_SIZE)
208 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
210 else if (maxsize > 10)
212 set_warning_line(wi, topfn, -1);
213 sprintf(warn_buf,
214 "The largest charge group contains %d atoms.\n"
215 "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"
216 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
217 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
218 maxsize);
219 warning_note(wi, warn_buf);
223 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
225 /* This check is not intended to ensure accurate integration,
226 * rather it is to signal mistakes in the mdp settings.
227 * A common mistake is to forget to turn on constraints
228 * for MD after energy minimization with flexible bonds.
229 * This check can also detect too large time steps for flexible water
230 * models, but such errors will often be masked by the constraints
231 * mdp options, which turns flexible water into water with bond constraints,
232 * but without an angle constraint. Unfortunately such incorrect use
233 * of water models can not easily be detected without checking
234 * for specific model names.
236 * The stability limit of leap-frog or velocity verlet is 4.44 steps
237 * per oscillational period.
238 * But accurate bonds distributions are lost far before that limit.
239 * To allow relatively common schemes (although not common with Gromacs)
240 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
241 * we set the note limit to 10.
243 int min_steps_warn = 5;
244 int min_steps_note = 10;
245 t_iparams *ip;
246 int molt;
247 gmx_moltype_t *moltype, *w_moltype;
248 t_atom *atom;
249 t_ilist *ilist, *ilb, *ilc, *ils;
250 int ftype;
251 int i, a1, a2, w_a1, w_a2, j;
252 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
253 gmx_bool bFound, bWater, bWarn;
254 char warn_buf[STRLEN];
256 ip = mtop->ffparams.iparams;
258 twopi2 = gmx::square(2*M_PI);
260 limit2 = gmx::square(min_steps_note*dt);
262 w_a1 = w_a2 = -1;
263 w_period2 = -1.0;
265 w_moltype = nullptr;
266 for (molt = 0; molt < mtop->nmoltype; molt++)
268 moltype = &mtop->moltype[molt];
269 atom = moltype->atoms.atom;
270 ilist = moltype->ilist;
271 ilc = &ilist[F_CONSTR];
272 ils = &ilist[F_SETTLE];
273 for (ftype = 0; ftype < F_NRE; ftype++)
275 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
277 continue;
280 ilb = &ilist[ftype];
281 for (i = 0; i < ilb->nr; i += 3)
283 fc = ip[ilb->iatoms[i]].harmonic.krA;
284 re = ip[ilb->iatoms[i]].harmonic.rA;
285 if (ftype == F_G96BONDS)
287 /* Convert squared sqaure fc to harmonic fc */
288 fc = 2*fc*re;
290 a1 = ilb->iatoms[i+1];
291 a2 = ilb->iatoms[i+2];
292 m1 = atom[a1].m;
293 m2 = atom[a2].m;
294 if (fc > 0 && m1 > 0 && m2 > 0)
296 period2 = twopi2*m1*m2/((m1 + m2)*fc);
298 else
300 period2 = GMX_FLOAT_MAX;
302 if (debug)
304 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
305 fc, m1, m2, std::sqrt(period2));
307 if (period2 < limit2)
309 bFound = FALSE;
310 for (j = 0; j < ilc->nr; j += 3)
312 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
313 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
315 bFound = TRUE;
318 for (j = 0; j < ils->nr; j += 4)
320 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
321 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
323 bFound = TRUE;
326 if (!bFound &&
327 (w_moltype == nullptr || period2 < w_period2))
329 w_moltype = moltype;
330 w_a1 = a1;
331 w_a2 = a2;
332 w_period2 = period2;
339 if (w_moltype != nullptr)
341 bWarn = (w_period2 < gmx::square(min_steps_warn*dt));
342 /* A check that would recognize most water models */
343 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
344 w_moltype->atoms.nr <= 5);
345 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"
346 "%s",
347 *w_moltype->name,
348 w_a1+1, *w_moltype->atoms.atomname[w_a1],
349 w_a2+1, *w_moltype->atoms.atomname[w_a2],
350 std::sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
351 bWater ?
352 "Maybe you asked for fexible water." :
353 "Maybe you forgot to change the constraints mdp option.");
354 if (bWarn)
356 warning(wi, warn_buf);
358 else
360 warning_note(wi, warn_buf);
365 static void check_vel(gmx_mtop_t *mtop, rvec v[])
367 gmx_mtop_atomloop_all_t aloop;
368 const t_atom *atom;
369 int a;
371 aloop = gmx_mtop_atomloop_all_init(mtop);
372 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
374 if (atom->ptype == eptShell ||
375 atom->ptype == eptBond ||
376 atom->ptype == eptVSite)
378 clear_rvec(v[a]);
383 static void check_shells_inputrec(gmx_mtop_t *mtop,
384 t_inputrec *ir,
385 warninp_t wi)
387 gmx_mtop_atomloop_all_t aloop;
388 const t_atom *atom;
389 int a, nshells = 0;
390 char warn_buf[STRLEN];
392 aloop = gmx_mtop_atomloop_all_init(mtop);
393 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
395 if (atom->ptype == eptShell ||
396 atom->ptype == eptBond)
398 nshells++;
401 if ((nshells > 0) && (ir->nstcalcenergy != 1))
403 set_warning_line(wi, "unknown", -1);
404 snprintf(warn_buf, STRLEN,
405 "There are %d shells, changing nstcalcenergy from %d to 1",
406 nshells, ir->nstcalcenergy);
407 ir->nstcalcenergy = 1;
408 warning(wi, warn_buf);
412 /* TODO Decide whether this function can be consolidated with
413 * gmx_mtop_ftype_count */
414 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
416 int nint, mb;
418 nint = 0;
419 for (mb = 0; mb < mtop->nmolblock; mb++)
421 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
424 return nint;
427 /* This routine reorders the molecule type array
428 * in the order of use in the molblocks,
429 * unused molecule types are deleted.
431 static void renumber_moltypes(gmx_mtop_t *sys,
432 int *nmolinfo, t_molinfo **molinfo)
434 int *order, norder, i;
435 int mb, mi;
436 t_molinfo *minew;
438 snew(order, *nmolinfo);
439 norder = 0;
440 for (mb = 0; mb < sys->nmolblock; mb++)
442 for (i = 0; i < norder; i++)
444 if (order[i] == sys->molblock[mb].type)
446 break;
449 if (i == norder)
451 /* This type did not occur yet, add it */
452 order[norder] = sys->molblock[mb].type;
453 /* Renumber the moltype in the topology */
454 norder++;
456 sys->molblock[mb].type = i;
459 /* We still need to reorder the molinfo structs */
460 snew(minew, norder);
461 for (mi = 0; mi < *nmolinfo; mi++)
463 for (i = 0; i < norder; i++)
465 if (order[i] == mi)
467 break;
470 if (i == norder)
472 done_mi(&(*molinfo)[mi]);
474 else
476 minew[i] = (*molinfo)[mi];
479 sfree(*molinfo);
481 *nmolinfo = norder;
482 *molinfo = minew;
485 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
487 int m;
488 gmx_moltype_t *molt;
490 mtop->nmoltype = nmi;
491 snew(mtop->moltype, nmi);
492 for (m = 0; m < nmi; m++)
494 molt = &mtop->moltype[m];
495 molt->name = mi[m].name;
496 molt->atoms = mi[m].atoms;
497 /* ilists are copied later */
498 molt->cgs = mi[m].cgs;
499 molt->excls = mi[m].excls;
503 static void
504 new_status(const char *topfile, const char *topppfile, const char *confin,
505 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
506 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
507 gpp_atomtype_t atype, gmx_mtop_t *sys,
508 int *nmi, t_molinfo **mi, t_molinfo **intermolecular_interactions,
509 t_params plist[],
510 int *comb, double *reppow, real *fudgeQQ,
511 gmx_bool bMorse,
512 warninp_t wi)
514 t_molinfo *molinfo = nullptr;
515 int nmolblock;
516 gmx_molblock_t *molblock, *molbs;
517 int mb, i, nrmols, nmismatch;
518 char buf[STRLEN];
519 gmx_bool bGB = FALSE;
520 char warn_buf[STRLEN];
522 init_mtop(sys);
524 /* Set gmx_boolean for GB */
525 if (ir->implicit_solvent)
527 bGB = TRUE;
530 /* TOPOLOGY processing */
531 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
532 plist, comb, reppow, fudgeQQ,
533 atype, &nrmols, &molinfo, intermolecular_interactions,
535 &nmolblock, &molblock, bGB,
536 wi);
538 sys->nmolblock = 0;
539 snew(sys->molblock, nmolblock);
541 sys->natoms = 0;
542 for (mb = 0; mb < nmolblock; mb++)
544 if (sys->nmolblock > 0 &&
545 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
547 /* Merge consecutive blocks with the same molecule type */
548 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
549 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
551 else if (molblock[mb].nmol > 0)
553 /* Add a new molblock to the topology */
554 molbs = &sys->molblock[sys->nmolblock];
555 *molbs = molblock[mb];
556 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
557 molbs->nposres_xA = 0;
558 molbs->nposres_xB = 0;
559 sys->natoms += molbs->nmol*molbs->natoms_mol;
560 sys->nmolblock++;
563 if (sys->nmolblock == 0)
565 gmx_fatal(FARGS, "No molecules were defined in the system");
568 renumber_moltypes(sys, &nrmols, &molinfo);
570 if (bMorse)
572 convert_harmonics(nrmols, molinfo, atype);
575 if (ir->eDisre == edrNone)
577 i = rm_interactions(F_DISRES, nrmols, molinfo);
578 if (i > 0)
580 set_warning_line(wi, "unknown", -1);
581 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
582 warning_note(wi, warn_buf);
585 if (opts->bOrire == FALSE)
587 i = rm_interactions(F_ORIRES, nrmols, molinfo);
588 if (i > 0)
590 set_warning_line(wi, "unknown", -1);
591 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
592 warning_note(wi, warn_buf);
596 /* Copy structures from msys to sys */
597 molinfo2mtop(nrmols, molinfo, sys);
599 gmx_mtop_finalize(sys);
601 /* COORDINATE file processing */
602 if (bVerbose)
604 fprintf(stderr, "processing coordinates...\n");
607 t_topology *conftop;
608 rvec *x = nullptr;
609 rvec *v = nullptr;
610 snew(conftop, 1);
611 read_tps_conf(confin, conftop, nullptr, &x, &v, state->box, FALSE);
612 state->natoms = conftop->atoms.nr;
613 if (state->natoms != sys->natoms)
615 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
616 " does not match topology (%s, %d)",
617 confin, state->natoms, topfile, sys->natoms);
619 /* It would be nice to get rid of the copies below, but we don't know
620 * a priori if the number of atoms in confin matches what we expect.
622 state->flags |= (1 << estX);
623 if (v != NULL)
625 state->flags |= (1 << estV);
627 state_change_natoms(state, state->natoms);
628 for (int i = 0; i < state->natoms; i++)
630 copy_rvec(x[i], state->x[i]);
632 sfree(x);
633 if (v != nullptr)
635 for (int i = 0; i < state->natoms; i++)
637 copy_rvec(v[i], state->v[i]);
639 sfree(v);
641 /* This call fixes the box shape for runs with pressure scaling */
642 set_box_rel(ir, state->box_rel, state->box);
644 nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
645 done_top(conftop);
646 sfree(conftop);
648 if (nmismatch)
650 sprintf(buf, "%d non-matching atom name%s\n"
651 "atom names from %s will be used\n"
652 "atom names from %s will be ignored\n",
653 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
654 warning(wi, buf);
657 /* Do more checks, mostly related to constraints */
658 if (bVerbose)
660 fprintf(stderr, "double-checking input for internal consistency...\n");
663 int bHasNormalConstraints = 0 < (nint_ftype(sys, molinfo, F_CONSTR) +
664 nint_ftype(sys, molinfo, F_CONSTRNC));
665 int bHasAnyConstraints = bHasNormalConstraints || 0 < nint_ftype(sys, molinfo, F_SETTLE);
666 double_check(ir, state->box,
667 bHasNormalConstraints,
668 bHasAnyConstraints,
669 wi);
672 if (bGenVel)
674 real *mass;
675 gmx_mtop_atomloop_all_t aloop;
676 const t_atom *atom;
678 snew(mass, state->natoms);
679 aloop = gmx_mtop_atomloop_all_init(sys);
680 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
682 mass[i] = atom->m;
685 if (opts->seed == -1)
687 opts->seed = static_cast<int>(gmx::makeRandomSeed());
688 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
690 state->flags |= (1 << estV);
691 maxwell_speed(opts->tempi, opts->seed, sys, as_rvec_array(state->v.data()));
693 stop_cm(stdout, state->natoms, mass, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()));
694 sfree(mass);
697 *nmi = nrmols;
698 *mi = molinfo;
701 static void copy_state(const char *slog, t_trxframe *fr,
702 gmx_bool bReadVel, t_state *state,
703 double *use_time)
705 int i;
707 if (fr->not_ok & FRAME_NOT_OK)
709 gmx_fatal(FARGS, "Can not start from an incomplete frame");
711 if (!fr->bX)
713 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
714 slog);
717 for (i = 0; i < state->natoms; i++)
719 copy_rvec(fr->x[i], state->x[i]);
721 if (bReadVel)
723 if (!fr->bV)
725 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
727 for (i = 0; i < state->natoms; i++)
729 copy_rvec(fr->v[i], state->v[i]);
732 if (fr->bBox)
734 copy_mat(fr->box, state->box);
737 *use_time = fr->time;
740 static void cont_status(const char *slog, const char *ener,
741 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
742 t_inputrec *ir, t_state *state,
743 gmx_mtop_t *sys,
744 const gmx_output_env_t *oenv)
745 /* If fr_time == -1 read the last frame available which is complete */
747 gmx_bool bReadVel;
748 t_trxframe fr;
749 t_trxstatus *fp;
750 int i;
751 double use_time;
753 bReadVel = (bNeedVel && !bGenVel);
755 fprintf(stderr,
756 "Reading Coordinates%s and Box size from old trajectory\n",
757 bReadVel ? ", Velocities" : "");
758 if (fr_time == -1)
760 fprintf(stderr, "Will read whole trajectory\n");
762 else
764 fprintf(stderr, "Will read till time %g\n", fr_time);
766 if (!bReadVel)
768 if (bGenVel)
770 fprintf(stderr, "Velocities generated: "
771 "ignoring velocities in input trajectory\n");
773 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
775 else
777 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
779 if (!fr.bV)
781 fprintf(stderr,
782 "\n"
783 "WARNING: Did not find a frame with velocities in file %s,\n"
784 " all velocities will be set to zero!\n\n", slog);
785 for (i = 0; i < sys->natoms; i++)
787 clear_rvec(state->v[i]);
789 close_trj(fp);
790 /* Search for a frame without velocities */
791 bReadVel = FALSE;
792 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
796 state->natoms = fr.natoms;
798 if (sys->natoms != state->natoms)
800 gmx_fatal(FARGS, "Number of atoms in Topology "
801 "is not the same as in Trajectory");
803 copy_state(slog, &fr, bReadVel, state, &use_time);
805 /* Find the appropriate frame */
806 while ((fr_time == -1 || fr.time < fr_time) &&
807 read_next_frame(oenv, fp, &fr))
809 copy_state(slog, &fr, bReadVel, state, &use_time);
812 close_trj(fp);
814 /* Set the relative box lengths for preserving the box shape.
815 * Note that this call can lead to differences in the last bit
816 * with respect to using gmx convert-tpr to create a [REF].tpx[ref] file.
818 set_box_rel(ir, state->box_rel, state->box);
820 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
821 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
823 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
825 get_enx_state(ener, use_time, &sys->groups, ir, state);
826 preserve_box_shape(ir, state->box_rel, state->boxv);
830 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
831 char *fn,
832 int rc_scaling, int ePBC,
833 rvec com,
834 warninp_t wi)
836 gmx_bool *hadAtom;
837 rvec *x, *v, *xp;
838 dvec sum;
839 double totmass;
840 t_topology *top;
841 matrix box, invbox;
842 int natoms, npbcdim = 0;
843 char warn_buf[STRLEN];
844 int a, i, ai, j, k, mb, nat_molb;
845 gmx_molblock_t *molb;
846 t_params *pr, *prfb;
847 t_atom *atom;
849 snew(top, 1);
850 read_tps_conf(fn, top, nullptr, &x, &v, box, FALSE);
851 natoms = top->atoms.nr;
852 done_top(top);
853 sfree(top);
854 if (natoms != mtop->natoms)
856 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);
857 warning(wi, warn_buf);
860 npbcdim = ePBC2npbcdim(ePBC);
861 clear_rvec(com);
862 if (rc_scaling != erscNO)
864 copy_mat(box, invbox);
865 for (j = npbcdim; j < DIM; j++)
867 clear_rvec(invbox[j]);
868 invbox[j][j] = 1;
870 gmx::invertBoxMatrix(invbox, invbox);
873 /* Copy the reference coordinates to mtop */
874 clear_dvec(sum);
875 totmass = 0;
876 a = 0;
877 snew(hadAtom, natoms);
878 for (mb = 0; mb < mtop->nmolblock; mb++)
880 molb = &mtop->molblock[mb];
881 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
882 pr = &(molinfo[molb->type].plist[F_POSRES]);
883 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
884 if (pr->nr > 0 || prfb->nr > 0)
886 atom = mtop->moltype[molb->type].atoms.atom;
887 for (i = 0; (i < pr->nr); i++)
889 ai = pr->param[i].ai();
890 if (ai >= natoms)
892 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
893 ai+1, *molinfo[molb->type].name, fn, natoms);
895 hadAtom[ai] = TRUE;
896 if (rc_scaling == erscCOM)
898 /* Determine the center of mass of the posres reference coordinates */
899 for (j = 0; j < npbcdim; j++)
901 sum[j] += atom[ai].m*x[a+ai][j];
903 totmass += atom[ai].m;
906 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
907 for (i = 0; (i < prfb->nr); i++)
909 ai = prfb->param[i].ai();
910 if (ai >= natoms)
912 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
913 ai+1, *molinfo[molb->type].name, fn, natoms);
915 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
917 /* Determine the center of mass of the posres reference coordinates */
918 for (j = 0; j < npbcdim; j++)
920 sum[j] += atom[ai].m*x[a+ai][j];
922 totmass += atom[ai].m;
925 if (!bTopB)
927 molb->nposres_xA = nat_molb;
928 snew(molb->posres_xA, molb->nposres_xA);
929 for (i = 0; i < nat_molb; i++)
931 copy_rvec(x[a+i], molb->posres_xA[i]);
934 else
936 molb->nposres_xB = nat_molb;
937 snew(molb->posres_xB, molb->nposres_xB);
938 for (i = 0; i < nat_molb; i++)
940 copy_rvec(x[a+i], molb->posres_xB[i]);
944 a += nat_molb;
946 if (rc_scaling == erscCOM)
948 if (totmass == 0)
950 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
952 for (j = 0; j < npbcdim; j++)
954 com[j] = sum[j]/totmass;
956 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]);
959 if (rc_scaling != erscNO)
961 GMX_ASSERT(npbcdim <= DIM, "Only DIM dimensions can have PBC");
963 for (mb = 0; mb < mtop->nmolblock; mb++)
965 molb = &mtop->molblock[mb];
966 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
967 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
969 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
970 for (i = 0; i < nat_molb; i++)
972 for (j = 0; j < npbcdim; j++)
974 if (rc_scaling == erscALL)
976 /* Convert from Cartesian to crystal coordinates */
977 xp[i][j] *= invbox[j][j];
978 for (k = j+1; k < npbcdim; k++)
980 xp[i][j] += invbox[k][j]*xp[i][k];
983 else if (rc_scaling == erscCOM)
985 /* Subtract the center of mass */
986 xp[i][j] -= com[j];
993 if (rc_scaling == erscCOM)
995 /* Convert the COM from Cartesian to crystal coordinates */
996 for (j = 0; j < npbcdim; j++)
998 com[j] *= invbox[j][j];
999 for (k = j+1; k < npbcdim; k++)
1001 com[j] += invbox[k][j]*com[k];
1007 sfree(x);
1008 sfree(v);
1009 sfree(hadAtom);
1012 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
1013 char *fnA, char *fnB,
1014 int rc_scaling, int ePBC,
1015 rvec com, rvec comB,
1016 warninp_t wi)
1018 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
1019 /* It is safer to simply read the b-state posres rather than trying
1020 * to be smart and copy the positions.
1022 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
1025 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
1026 t_inputrec *ir, warninp_t wi)
1028 int i;
1029 char warn_buf[STRLEN];
1031 if (ir->nwall > 0)
1033 fprintf(stderr, "Searching the wall atom type(s)\n");
1035 for (i = 0; i < ir->nwall; i++)
1037 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1038 if (ir->wall_atomtype[i] == NOTSET)
1040 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1041 warning_error(wi, warn_buf);
1046 static int nrdf_internal(t_atoms *atoms)
1048 int i, nmass, nrdf;
1050 nmass = 0;
1051 for (i = 0; i < atoms->nr; i++)
1053 /* Vsite ptype might not be set here yet, so also check the mass */
1054 if ((atoms->atom[i].ptype == eptAtom ||
1055 atoms->atom[i].ptype == eptNucleus)
1056 && atoms->atom[i].m > 0)
1058 nmass++;
1061 switch (nmass)
1063 case 0: nrdf = 0; break;
1064 case 1: nrdf = 0; break;
1065 case 2: nrdf = 1; break;
1066 default: nrdf = nmass*3 - 6; break;
1069 return nrdf;
1072 void
1073 spline1d( double dx,
1074 double * y,
1075 int n,
1076 double * u,
1077 double * y2 )
1079 int i;
1080 double p, q;
1082 y2[0] = 0.0;
1083 u[0] = 0.0;
1085 for (i = 1; i < n-1; i++)
1087 p = 0.5*y2[i-1]+2.0;
1088 y2[i] = -0.5/p;
1089 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1090 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1093 y2[n-1] = 0.0;
1095 for (i = n-2; i >= 0; i--)
1097 y2[i] = y2[i]*y2[i+1]+u[i];
1102 void
1103 interpolate1d( double xmin,
1104 double dx,
1105 double * ya,
1106 double * y2a,
1107 double x,
1108 double * y,
1109 double * y1)
1111 int ix;
1112 double a, b;
1114 ix = static_cast<int>((x-xmin)/dx);
1116 a = (xmin+(ix+1)*dx-x)/dx;
1117 b = (x-xmin-ix*dx)/dx;
1119 *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;
1120 *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];
1124 void
1125 setup_cmap (int grid_spacing,
1126 int nc,
1127 real * grid,
1128 gmx_cmap_t * cmap_grid)
1130 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1132 int i, j, k, ii, jj, kk, idx;
1133 int offset;
1134 double dx, xmin, v, v1, v2, v12;
1135 double phi, psi;
1137 snew(tmp_u, 2*grid_spacing);
1138 snew(tmp_u2, 2*grid_spacing);
1139 snew(tmp_yy, 2*grid_spacing);
1140 snew(tmp_y1, 2*grid_spacing);
1141 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1142 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1144 dx = 360.0/grid_spacing;
1145 xmin = -180.0-dx*grid_spacing/2;
1147 for (kk = 0; kk < nc; kk++)
1149 /* Compute an offset depending on which cmap we are using
1150 * Offset will be the map number multiplied with the
1151 * grid_spacing * grid_spacing * 2
1153 offset = kk * grid_spacing * grid_spacing * 2;
1155 for (i = 0; i < 2*grid_spacing; i++)
1157 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1159 for (j = 0; j < 2*grid_spacing; j++)
1161 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1162 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1166 for (i = 0; i < 2*grid_spacing; i++)
1168 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1171 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1173 ii = i-grid_spacing/2;
1174 phi = ii*dx-180.0;
1176 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1178 jj = j-grid_spacing/2;
1179 psi = jj*dx-180.0;
1181 for (k = 0; k < 2*grid_spacing; k++)
1183 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1184 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1187 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1188 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1189 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1190 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1192 idx = ii*grid_spacing+jj;
1193 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1194 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1195 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1196 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1202 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1204 int i, nelem;
1206 cmap_grid->ngrid = ngrid;
1207 cmap_grid->grid_spacing = grid_spacing;
1208 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1210 snew(cmap_grid->cmapdata, ngrid);
1212 for (i = 0; i < cmap_grid->ngrid; i++)
1214 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1219 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1221 int count, count_mol, i, mb;
1222 gmx_molblock_t *molb;
1223 t_params *plist;
1224 char buf[STRLEN];
1226 count = 0;
1227 for (mb = 0; mb < mtop->nmolblock; mb++)
1229 count_mol = 0;
1230 molb = &mtop->molblock[mb];
1231 plist = mi[molb->type].plist;
1233 for (i = 0; i < F_NRE; i++)
1235 if (i == F_SETTLE)
1237 count_mol += 3*plist[i].nr;
1239 else if (interaction_function[i].flags & IF_CONSTRAINT)
1241 count_mol += plist[i].nr;
1245 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1247 sprintf(buf,
1248 "Molecule type '%s' has %d constraints.\n"
1249 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1250 *mi[molb->type].name, count_mol,
1251 nrdf_internal(&mi[molb->type].atoms));
1252 warning(wi, buf);
1254 count += molb->nmol*count_mol;
1257 return count;
1260 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1262 int i, nmiss, natoms, mt;
1263 real q;
1264 const t_atoms *atoms;
1266 nmiss = 0;
1267 for (mt = 0; mt < sys->nmoltype; mt++)
1269 atoms = &sys->moltype[mt].atoms;
1270 natoms = atoms->nr;
1272 for (i = 0; i < natoms; i++)
1274 q = atoms->atom[i].q;
1275 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1276 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1277 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1278 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1279 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1280 q != 0)
1282 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1283 get_atomtype_name(atoms->atom[i].type, atype), q);
1284 nmiss++;
1289 if (nmiss > 0)
1291 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);
1296 static void check_gbsa_params(gpp_atomtype_t atype)
1298 int nmiss, i;
1300 /* If we are doing GBSA, check that we got the parameters we need
1301 * This checking is to see if there are GBSA paratmeters for all
1302 * atoms in the force field. To go around this for testing purposes
1303 * comment out the nerror++ counter temporarily
1305 nmiss = 0;
1306 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1308 if (get_atomtype_radius(i, atype) < 0 ||
1309 get_atomtype_vol(i, atype) < 0 ||
1310 get_atomtype_surftens(i, atype) < 0 ||
1311 get_atomtype_gb_radius(i, atype) < 0 ||
1312 get_atomtype_S_hct(i, atype) < 0)
1314 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1315 get_atomtype_name(i, atype));
1316 nmiss++;
1320 if (nmiss > 0)
1322 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);
1327 static real calc_temp(const gmx_mtop_t *mtop,
1328 const t_inputrec *ir,
1329 rvec *v)
1331 gmx_mtop_atomloop_all_t aloop;
1332 const t_atom *atom;
1333 int a;
1335 double sum_mv2 = 0;
1336 aloop = gmx_mtop_atomloop_all_init(mtop);
1337 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1339 sum_mv2 += atom->m*norm2(v[a]);
1342 double nrdf = 0;
1343 for (int g = 0; g < ir->opts.ngtc; g++)
1345 nrdf += ir->opts.nrdf[g];
1348 return sum_mv2/(nrdf*BOLTZ);
1351 static real get_max_reference_temp(const t_inputrec *ir,
1352 warninp_t wi)
1354 real ref_t;
1355 int i;
1356 gmx_bool bNoCoupl;
1358 ref_t = 0;
1359 bNoCoupl = FALSE;
1360 for (i = 0; i < ir->opts.ngtc; i++)
1362 if (ir->opts.tau_t[i] < 0)
1364 bNoCoupl = TRUE;
1366 else
1368 ref_t = std::max(ref_t, ir->opts.ref_t[i]);
1372 if (bNoCoupl)
1374 char buf[STRLEN];
1376 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.",
1377 ref_t);
1378 warning(wi, buf);
1381 return ref_t;
1384 /* Checks if there are unbound atoms in moleculetype molt.
1385 * Prints a note for each unbound atoms and a warning if any is present.
1387 static void checkForUnboundAtoms(const gmx_moltype_t *molt,
1388 gmx_bool bVerbose,
1389 warninp_t wi)
1391 const t_atoms *atoms = &molt->atoms;
1393 if (atoms->nr == 1)
1395 /* Only one atom, there can't be unbound atoms */
1396 return;
1399 std::vector<int> count(atoms->nr, 0);
1401 for (int ftype = 0; ftype < F_NRE; ftype++)
1403 if (((interaction_function[ftype].flags & IF_BOND) && ftype != F_CONNBONDS) ||
1404 (interaction_function[ftype].flags & IF_CONSTRAINT) ||
1405 ftype == F_SETTLE)
1407 const t_ilist *il = &molt->ilist[ftype];
1408 int nral = NRAL(ftype);
1410 for (int i = 0; i < il->nr; i += 1 + nral)
1412 for (int j = 0; j < nral; j++)
1414 count[il->iatoms[i + 1 + j]]++;
1420 int numDanglingAtoms = 0;
1421 for (int a = 0; a < atoms->nr; a++)
1423 if (atoms->atom[a].ptype != eptVSite &&
1424 count[a] == 0)
1426 if (bVerbose)
1428 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",
1429 a + 1, *atoms->atomname[a], *molt->name);
1431 numDanglingAtoms++;
1435 if (numDanglingAtoms > 0)
1437 char buf[STRLEN];
1438 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.",
1439 *molt->name, numDanglingAtoms);
1440 warning_note(wi, buf);
1444 /* Checks all moleculetypes for unbound atoms */
1445 static void checkForUnboundAtoms(const gmx_mtop_t *mtop,
1446 gmx_bool bVerbose,
1447 warninp_t wi)
1449 for (int mt = 0; mt < mtop->nmoltype; mt++)
1451 checkForUnboundAtoms(&mtop->moltype[mt], bVerbose, wi);
1455 /*! \brief Checks if there are decoupled modes in moleculetype \p molt.
1457 * The specific decoupled modes this routine check for are angle modes
1458 * where the two bonds are constrained and the atoms a both ends are only
1459 * involved in a single constraint; the mass of the two atoms needs to
1460 * differ by more than \p massFactorThreshold.
1462 static bool haveDecoupledModeInMol(const gmx_moltype_t *molt,
1463 const t_iparams *iparams,
1464 real massFactorThreshold)
1466 if (molt->ilist[F_CONSTR].nr == 0 && molt->ilist[F_CONSTRNC].nr == 0)
1468 return false;
1471 const t_atom * atom = molt->atoms.atom;
1473 int numFlexibleConstraints;
1474 t_blocka atomToConstraints = make_at2con(0, molt->atoms.nr,
1475 molt->ilist, iparams,
1476 FALSE,
1477 &numFlexibleConstraints);
1479 bool haveDecoupledMode = false;
1480 for (int ftype = 0; ftype < F_NRE; ftype++)
1482 if (interaction_function[ftype].flags & IF_ATYPE)
1484 const int nral = NRAL(ftype);
1485 const t_ilist *il = &molt->ilist[ftype];
1486 for (int i = 0; i < il->nr; i += 1 + nral)
1488 /* Here we check for the mass difference between the atoms
1489 * at both ends of the angle, that the atoms at the ends have
1490 * 1 contraint and the atom in the middle at least 3; we check
1491 * that the 3 atoms are linked by constraints below.
1492 * We check for at least three constraints for the middle atom,
1493 * since with only the two bonds in the angle, we have 3-atom
1494 * molecule, which has much more thermal exhange in this single
1495 * angle mode than molecules with more atoms.
1496 * Note that this check also catches molecules where more atoms
1497 * are connected to one or more atoms in the angle, but by
1498 * bond potentials instead of angles. But such cases will not
1499 * occur in "normal" molecules and it doesn't hurt running
1500 * those with higher accuracy settings as well.
1502 int a0 = il->iatoms[1 + i];
1503 int a1 = il->iatoms[1 + i + 1];
1504 int a2 = il->iatoms[1 + i + 2];
1505 if ((atom[a0].m > atom[a2].m*massFactorThreshold ||
1506 atom[a2].m > atom[a0].m*massFactorThreshold) &&
1507 atomToConstraints.index[a0 + 1] - atomToConstraints.index[a0] == 1 &&
1508 atomToConstraints.index[a2 + 1] - atomToConstraints.index[a2] == 1 &&
1509 atomToConstraints.index[a1 + 1] - atomToConstraints.index[a1] >= 3)
1511 int constraint0 = atomToConstraints.a[atomToConstraints.index[a0]];
1512 int constraint2 = atomToConstraints.a[atomToConstraints.index[a2]];
1514 bool foundAtom0 = false;
1515 bool foundAtom2 = false;
1516 for (int conIndex = atomToConstraints.index[a1]; conIndex < atomToConstraints.index[a1 + 1]; conIndex++)
1518 if (atomToConstraints.a[conIndex] == constraint0)
1520 foundAtom0 = true;
1522 if (atomToConstraints.a[conIndex] == constraint2)
1524 foundAtom2 = true;
1527 if (foundAtom0 && foundAtom2)
1529 haveDecoupledMode = true;
1536 done_blocka(&atomToConstraints);
1538 return haveDecoupledMode;
1541 /*! \brief Checks if the Verlet buffer and constraint accuracy is sufficient for decoupled dynamic modes.
1543 * When decoupled modes are present and the accuray in insufficient,
1544 * this routine issues a warning if the accuracy is insufficient.
1546 static void checkDecoupledModeAccuracy(const gmx_mtop_t *mtop,
1547 const t_inputrec *ir,
1548 warninp_t wi)
1550 /* We only have issues with decoupled modes with normal MD.
1551 * With stochastic dynamics equipartitioning is enforced strongly.
1553 if (!EI_MD(ir->eI))
1555 return;
1558 /* When atoms of very different mass are involved in an angle potential
1559 * and both bonds in the angle are constrained, the dynamic modes in such
1560 * angles have very different periods and significant energy exchange
1561 * takes several nanoseconds. Thus even a small amount of error in
1562 * different algorithms can lead to violation of equipartitioning.
1563 * The parameters below are mainly based on an all-atom chloroform model
1564 * with all bonds constrained. Equipartitioning is off by more than 1%
1565 * (need to run 5-10 ns) when the difference in mass between H and Cl
1566 * is more than a factor 13 and the accuracy is less than the thresholds
1567 * given below. This has been verified on some other molecules.
1569 * Note that the buffer and shake parameters have unit length and
1570 * energy/time, respectively, so they will "only" work correctly
1571 * for atomistic force fields using MD units.
1573 const real massFactorThreshold = 13.0;
1574 const real bufferToleranceThreshold = 1e-4;
1575 const int lincsIterationThreshold = 2;
1576 const int lincsOrderThreshold = 4;
1577 const real shakeToleranceThreshold = 0.005*ir->delta_t;
1579 bool lincsWithSufficientTolerance = (ir->eConstrAlg == econtLINCS && ir->nLincsIter >= lincsIterationThreshold && ir->nProjOrder >= lincsOrderThreshold);
1580 bool shakeWithSufficientTolerance = (ir->eConstrAlg == econtSHAKE && ir->shake_tol <= 1.1*shakeToleranceThreshold);
1581 if (ir->cutoff_scheme == ecutsVERLET &&
1582 ir->verletbuf_tol <= 1.1*bufferToleranceThreshold &&
1583 (lincsWithSufficientTolerance || shakeWithSufficientTolerance))
1585 return;
1588 bool haveDecoupledMode = false;
1589 for (int mt = 0; mt < mtop->nmoltype; mt++)
1591 if (haveDecoupledModeInMol(&mtop->moltype[mt], mtop->ffparams.iparams,
1592 massFactorThreshold))
1594 haveDecoupledMode = true;
1598 if (haveDecoupledMode)
1600 char modeMessage[STRLEN];
1601 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.",
1602 massFactorThreshold);
1603 char buf[STRLEN];
1604 if (ir->cutoff_scheme == ecutsVERLET)
1606 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",
1607 modeMessage,
1608 ei_names[eiSD1],
1609 bufferToleranceThreshold,
1610 lincsIterationThreshold, lincsOrderThreshold,
1611 shakeToleranceThreshold);
1613 else
1615 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.",
1616 modeMessage,
1617 ecutscheme_names[ecutsVERLET]);
1619 warning(wi, buf);
1623 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1624 t_inputrec *ir,
1625 real buffer_temp,
1626 matrix box,
1627 warninp_t wi)
1629 verletbuf_list_setup_t ls;
1630 real rlist_1x1;
1631 int n_nonlin_vsite;
1632 char warn_buf[STRLEN];
1634 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1636 /* Calculate the buffer size for simple atom vs atoms list */
1637 ls.cluster_size_i = 1;
1638 ls.cluster_size_j = 1;
1639 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1640 &ls, &n_nonlin_vsite, &rlist_1x1);
1642 /* Set the pair-list buffer size in ir */
1643 verletbuf_get_list_setup(FALSE, FALSE, &ls);
1644 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1645 &ls, &n_nonlin_vsite, &ir->rlist);
1647 if (n_nonlin_vsite > 0)
1649 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);
1650 warning_note(wi, warn_buf);
1653 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1654 1, 1, rlist_1x1, rlist_1x1-std::max(ir->rvdw, ir->rcoulomb));
1656 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1657 ls.cluster_size_i, ls.cluster_size_j,
1658 ir->rlist, ir->rlist-std::max(ir->rvdw, ir->rcoulomb));
1660 printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
1662 if (gmx::square(ir->rlist) >= max_cutoff2(ir->ePBC, box))
1664 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)));
1668 int gmx_grompp(int argc, char *argv[])
1670 static const char *desc[] = {
1671 "[THISMODULE] (the gromacs preprocessor)",
1672 "reads a molecular topology file, checks the validity of the",
1673 "file, expands the topology from a molecular description to an atomic",
1674 "description. The topology file contains information about",
1675 "molecule types and the number of molecules, the preprocessor",
1676 "copies each molecule as needed. ",
1677 "There is no limitation on the number of molecule types. ",
1678 "Bonds and bond-angles can be converted into constraints, separately",
1679 "for hydrogens and heavy atoms.",
1680 "Then a coordinate file is read and velocities can be generated",
1681 "from a Maxwellian distribution if requested.",
1682 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1683 "(eg. number of MD steps, time step, cut-off), and others such as",
1684 "NEMD parameters, which are corrected so that the net acceleration",
1685 "is zero.",
1686 "Eventually a binary file is produced that can serve as the sole input",
1687 "file for the MD program.[PAR]",
1689 "[THISMODULE] uses the atom names from the topology file. The atom names",
1690 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1691 "warnings when they do not match the atom names in the topology.",
1692 "Note that the atom names are irrelevant for the simulation as",
1693 "only the atom types are used for generating interaction parameters.[PAR]",
1695 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1696 "etc. The preprocessor supports the following keywords::",
1698 " #ifdef VARIABLE",
1699 " #ifndef VARIABLE",
1700 " #else",
1701 " #endif",
1702 " #define VARIABLE",
1703 " #undef VARIABLE",
1704 " #include \"filename\"",
1705 " #include <filename>",
1707 "The functioning of these statements in your topology may be modulated by",
1708 "using the following two flags in your [REF].mdp[ref] file::",
1710 " define = -DVARIABLE1 -DVARIABLE2",
1711 " include = -I/home/john/doe",
1713 "For further information a C-programming textbook may help you out.",
1714 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1715 "topology file written out so that you can verify its contents.[PAR]",
1717 "When using position restraints a file with restraint coordinates",
1718 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1719 "with respect to the conformation from the [TT]-c[tt] option.",
1720 "For free energy calculation the the coordinates for the B topology",
1721 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1722 "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 char fn[STRLEN], fnB[STRLEN];
1780 const char *mdparin;
1781 int ntype;
1782 gmx_bool bNeedVel, bGenVel;
1783 gmx_bool have_atomnumber;
1784 gmx_output_env_t *oenv;
1785 gmx_bool bVerbose = FALSE;
1786 warninp_t wi;
1787 char warn_buf[STRLEN];
1789 t_filenm fnm[] = {
1790 { efMDP, nullptr, nullptr, ffREAD },
1791 { efMDP, "-po", "mdout", ffWRITE },
1792 { efSTX, "-c", nullptr, ffREAD },
1793 { efSTX, "-r", nullptr, ffOPTRD },
1794 { efSTX, "-rb", nullptr, ffOPTRD },
1795 { efNDX, nullptr, nullptr, ffOPTRD },
1796 { efTOP, nullptr, nullptr, ffREAD },
1797 { efTOP, "-pp", "processed", ffOPTWR },
1798 { efTPR, "-o", nullptr, ffWRITE },
1799 { efTRN, "-t", nullptr, ffOPTRD },
1800 { efEDR, "-e", nullptr, ffOPTRD },
1801 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1802 { efGRO, "-imd", "imdgroup", ffOPTWR },
1803 { efTRN, "-ref", "rotref", ffOPTRW }
1805 #define NFILE asize(fnm)
1807 /* Command line options */
1808 static gmx_bool bRenum = TRUE;
1809 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1810 static int i, maxwarn = 0;
1811 static real fr_time = -1;
1812 t_pargs pa[] = {
1813 { "-v", FALSE, etBOOL, {&bVerbose},
1814 "Be loud and noisy" },
1815 { "-time", FALSE, etREAL, {&fr_time},
1816 "Take frame at or first after this time." },
1817 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1818 "Remove constant bonded interactions with virtual sites" },
1819 { "-maxwarn", FALSE, etINT, {&maxwarn},
1820 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1821 { "-zero", FALSE, etBOOL, {&bZero},
1822 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1823 { "-renum", FALSE, etBOOL, {&bRenum},
1824 "Renumber atomtypes and minimize number of atomtypes" }
1827 /* Parse the command line */
1828 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1829 asize(desc), desc, 0, nullptr, &oenv))
1831 return 0;
1834 /* Initiate some variables */
1835 gmx::MDModules mdModules;
1836 snew(opts, 1);
1837 snew(opts->include, STRLEN);
1838 snew(opts->define, STRLEN);
1840 wi = init_warning(TRUE, maxwarn);
1842 /* PARAMETER file processing */
1843 mdparin = opt2fn("-f", NFILE, fnm);
1844 set_warning_line(wi, mdparin, -1);
1845 get_ir(mdparin, opt2fn("-po", NFILE, fnm), &mdModules, opts, wi);
1846 t_inputrec *ir = mdModules.inputrec();
1848 if (bVerbose)
1850 fprintf(stderr, "checking input for internal consistency...\n");
1852 check_ir(mdparin, ir, opts, wi);
1854 if (ir->ld_seed == -1)
1856 ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
1857 fprintf(stderr, "Setting the LD random seed to %" GMX_PRId64 "\n", ir->ld_seed);
1860 if (ir->expandedvals->lmc_seed == -1)
1862 ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
1863 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1866 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1867 bGenVel = (bNeedVel && opts->bGenVel);
1868 if (bGenVel && ir->bContinuation)
1870 sprintf(warn_buf,
1871 "Generating velocities is inconsistent with attempting "
1872 "to continue a previous run. Choose only one of "
1873 "gen-vel = yes and continuation = yes.");
1874 warning_error(wi, warn_buf);
1877 snew(plist, F_NRE);
1878 init_plist(plist);
1879 snew(sys, 1);
1880 atype = init_atomtype();
1881 if (debug)
1883 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1886 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1887 if (!gmx_fexist(fn))
1889 gmx_fatal(FARGS, "%s does not exist", fn);
1892 t_state state;
1893 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1894 opts, ir, bZero, bGenVel, bVerbose, &state,
1895 atype, sys, &nmi, &mi, &intermolecular_interactions,
1896 plist, &comb, &reppow, &fudgeQQ,
1897 opts->bMorse,
1898 wi);
1900 if (debug)
1902 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1905 nvsite = 0;
1906 /* set parameters for virtual site construction (not for vsiten) */
1907 for (mt = 0; mt < sys->nmoltype; mt++)
1909 nvsite +=
1910 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1912 /* now throw away all obsolete bonds, angles and dihedrals: */
1913 /* note: constraints are ALWAYS removed */
1914 if (nvsite)
1916 for (mt = 0; mt < sys->nmoltype; mt++)
1918 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1922 if (ir->cutoff_scheme == ecutsVERLET)
1924 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1925 ecutscheme_names[ir->cutoff_scheme]);
1927 /* Remove all charge groups */
1928 gmx_mtop_remove_chargegroups(sys);
1931 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1933 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1935 sprintf(warn_buf, "Can not do %s with %s, use %s",
1936 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1937 warning_error(wi, warn_buf);
1939 if (ir->bPeriodicMols)
1941 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1942 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1943 warning_error(wi, warn_buf);
1947 if (EI_SD (ir->eI) && ir->etc != etcNO)
1949 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1952 /* If we are doing QM/MM, check that we got the atom numbers */
1953 have_atomnumber = TRUE;
1954 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1956 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1958 if (!have_atomnumber && ir->bQMMM)
1960 warning_error(wi,
1961 "\n"
1962 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1963 "field you are using does not contain atom numbers fields. This is an\n"
1964 "optional field (introduced in GROMACS 3.3) for general runs, but mandatory\n"
1965 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1966 "an integer just before the mass column in ffXXXnb.itp.\n"
1967 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1970 /* Check for errors in the input now, since they might cause problems
1971 * during processing further down.
1973 check_warning_error(wi, FARGS);
1975 if (opt2bSet("-r", NFILE, fnm))
1977 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1979 else
1981 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1983 if (opt2bSet("-rb", NFILE, fnm))
1985 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1987 else
1989 strcpy(fnB, fn);
1992 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
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;