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, 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.
43 #include <sys/types.h>
57 #include "gromacs/fileio/confio.h"
61 #include "grompp-impl.h"
62 #include "gromacs/random/random.h"
63 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
65 #include "gromacs/fileio/futil.h"
66 #include "gromacs/commandline/pargs.h"
68 #include "gromacs/gmxpreprocess/sortwater.h"
70 #include "gmx_fatal.h"
73 #include "gromacs/fileio/gmxfio.h"
74 #include "gromacs/fileio/trnio.h"
75 #include "gromacs/fileio/tpxio.h"
76 #include "gromacs/fileio/trxio.h"
77 #include "vsite_parm.h"
81 #include "gromacs/fileio/enxio.h"
83 #include "compute_io.h"
84 #include "gpp_atomtype.h"
85 #include "mtop_util.h"
87 #include "calc_verletbuf.h"
89 #include "gromacs/imd/imd.h"
92 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
97 /* For all the molecule types */
98 for (i = 0; i < nrmols; i++)
100 n += mols[i].plist[ifunc].nr;
101 mols[i].plist[ifunc].nr = 0;
106 static int check_atom_names(const char *fn1, const char *fn2,
107 gmx_mtop_t *mtop, t_atoms *at)
109 int mb, m, i, j, nmismatch;
111 #define MAXMISMATCH 20
113 if (mtop->natoms != at->nr)
115 gmx_incons("comparing atom names");
120 for (mb = 0; mb < mtop->nmolblock; mb++)
122 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
123 for (m = 0; m < mtop->molblock[mb].nmol; m++)
125 for (j = 0; j < tat->nr; j++)
127 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
129 if (nmismatch < MAXMISMATCH)
132 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
133 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
135 else if (nmismatch == MAXMISMATCH)
137 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
149 static void check_eg_vs_cg(gmx_mtop_t *mtop)
151 int astart, mb, m, cg, j, firstj;
152 unsigned char firsteg, eg;
155 /* Go through all the charge groups and make sure all their
156 * atoms are in the same energy group.
160 for (mb = 0; mb < mtop->nmolblock; mb++)
162 molt = &mtop->moltype[mtop->molblock[mb].type];
163 for (m = 0; m < mtop->molblock[mb].nmol; m++)
165 for (cg = 0; cg < molt->cgs.nr; cg++)
167 /* Get the energy group of the first atom in this charge group */
168 firstj = astart + molt->cgs.index[cg];
169 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
170 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
172 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
175 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
176 firstj+1, astart+j+1, cg+1, *molt->name);
180 astart += molt->atoms.nr;
185 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
188 char warn_buf[STRLEN];
191 for (cg = 0; cg < cgs->nr; cg++)
193 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
196 if (maxsize > MAX_CHARGEGROUP_SIZE)
198 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
200 else if (maxsize > 10)
202 set_warning_line(wi, topfn, -1);
204 "The largest charge group contains %d atoms.\n"
205 "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"
206 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
207 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
209 warning_note(wi, warn_buf);
213 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
215 /* This check is not intended to ensure accurate integration,
216 * rather it is to signal mistakes in the mdp settings.
217 * A common mistake is to forget to turn on constraints
218 * for MD after energy minimization with flexible bonds.
219 * This check can also detect too large time steps for flexible water
220 * models, but such errors will often be masked by the constraints
221 * mdp options, which turns flexible water into water with bond constraints,
222 * but without an angle constraint. Unfortunately such incorrect use
223 * of water models can not easily be detected without checking
224 * for specific model names.
226 * The stability limit of leap-frog or velocity verlet is 4.44 steps
227 * per oscillational period.
228 * But accurate bonds distributions are lost far before that limit.
229 * To allow relatively common schemes (although not common with Gromacs)
230 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
231 * we set the note limit to 10.
233 int min_steps_warn = 5;
234 int min_steps_note = 10;
237 gmx_moltype_t *moltype, *w_moltype;
239 t_ilist *ilist, *ilb, *ilc, *ils;
241 int i, a1, a2, w_a1, w_a2, j;
242 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
243 gmx_bool bFound, bWater, bWarn;
244 char warn_buf[STRLEN];
246 ip = mtop->ffparams.iparams;
248 twopi2 = sqr(2*M_PI);
250 limit2 = sqr(min_steps_note*dt);
256 for (molt = 0; molt < mtop->nmoltype; molt++)
258 moltype = &mtop->moltype[molt];
259 atom = moltype->atoms.atom;
260 ilist = moltype->ilist;
261 ilc = &ilist[F_CONSTR];
262 ils = &ilist[F_SETTLE];
263 for (ftype = 0; ftype < F_NRE; ftype++)
265 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
271 for (i = 0; i < ilb->nr; i += 3)
273 fc = ip[ilb->iatoms[i]].harmonic.krA;
274 re = ip[ilb->iatoms[i]].harmonic.rA;
275 if (ftype == F_G96BONDS)
277 /* Convert squared sqaure fc to harmonic fc */
280 a1 = ilb->iatoms[i+1];
281 a2 = ilb->iatoms[i+2];
284 if (fc > 0 && m1 > 0 && m2 > 0)
286 period2 = twopi2*m1*m2/((m1 + m2)*fc);
290 period2 = GMX_FLOAT_MAX;
294 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
295 fc, m1, m2, sqrt(period2));
297 if (period2 < limit2)
300 for (j = 0; j < ilc->nr; j += 3)
302 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
303 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
308 for (j = 0; j < ils->nr; j += 4)
310 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
311 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
317 (w_moltype == NULL || period2 < w_period2))
329 if (w_moltype != NULL)
331 bWarn = (w_period2 < sqr(min_steps_warn*dt));
332 /* A check that would recognize most water models */
333 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
334 w_moltype->atoms.nr <= 5);
335 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"
338 w_a1+1, *w_moltype->atoms.atomname[w_a1],
339 w_a2+1, *w_moltype->atoms.atomname[w_a2],
340 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
342 "Maybe you asked for fexible water." :
343 "Maybe you forgot to change the constraints mdp option.");
346 warning(wi, warn_buf);
350 warning_note(wi, warn_buf);
355 static void check_vel(gmx_mtop_t *mtop, rvec v[])
357 gmx_mtop_atomloop_all_t aloop;
361 aloop = gmx_mtop_atomloop_all_init(mtop);
362 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
364 if (atom->ptype == eptShell ||
365 atom->ptype == eptBond ||
366 atom->ptype == eptVSite)
373 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
378 for (mb = 0; mb < mtop->nmolblock; mb++)
380 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
386 /* This routine reorders the molecule type array
387 * in the order of use in the molblocks,
388 * unused molecule types are deleted.
390 static void renumber_moltypes(gmx_mtop_t *sys,
391 int *nmolinfo, t_molinfo **molinfo)
393 int *order, norder, i;
397 snew(order, *nmolinfo);
399 for (mb = 0; mb < sys->nmolblock; mb++)
401 for (i = 0; i < norder; i++)
403 if (order[i] == sys->molblock[mb].type)
410 /* This type did not occur yet, add it */
411 order[norder] = sys->molblock[mb].type;
412 /* Renumber the moltype in the topology */
415 sys->molblock[mb].type = i;
418 /* We still need to reorder the molinfo structs */
420 for (mi = 0; mi < *nmolinfo; mi++)
422 for (i = 0; i < norder; i++)
431 done_mi(&(*molinfo)[mi]);
435 minew[i] = (*molinfo)[mi];
444 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
449 mtop->nmoltype = nmi;
450 snew(mtop->moltype, nmi);
451 for (m = 0; m < nmi; m++)
453 molt = &mtop->moltype[m];
454 molt->name = mi[m].name;
455 molt->atoms = mi[m].atoms;
456 /* ilists are copied later */
457 molt->cgs = mi[m].cgs;
458 molt->excls = mi[m].excls;
463 new_status(const char *topfile, const char *topppfile, const char *confin,
464 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
465 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
466 gpp_atomtype_t atype, gmx_mtop_t *sys,
467 int *nmi, t_molinfo **mi, t_params plist[],
468 int *comb, double *reppow, real *fudgeQQ,
472 t_molinfo *molinfo = NULL;
474 gmx_molblock_t *molblock, *molbs;
476 int mb, i, nrmols, nmismatch;
478 gmx_bool bGB = FALSE;
479 char warn_buf[STRLEN];
483 /* Set gmx_boolean for GB */
484 if (ir->implicit_solvent)
489 /* TOPOLOGY processing */
490 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
491 plist, comb, reppow, fudgeQQ,
492 atype, &nrmols, &molinfo, ir,
493 &nmolblock, &molblock, bGB,
497 snew(sys->molblock, nmolblock);
500 for (mb = 0; mb < nmolblock; mb++)
502 if (sys->nmolblock > 0 &&
503 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
505 /* Merge consecutive blocks with the same molecule type */
506 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
507 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
509 else if (molblock[mb].nmol > 0)
511 /* Add a new molblock to the topology */
512 molbs = &sys->molblock[sys->nmolblock];
513 *molbs = molblock[mb];
514 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
515 molbs->nposres_xA = 0;
516 molbs->nposres_xB = 0;
517 sys->natoms += molbs->nmol*molbs->natoms_mol;
521 if (sys->nmolblock == 0)
523 gmx_fatal(FARGS, "No molecules were defined in the system");
526 renumber_moltypes(sys, &nrmols, &molinfo);
530 convert_harmonics(nrmols, molinfo, atype);
533 if (ir->eDisre == edrNone)
535 i = rm_interactions(F_DISRES, nrmols, molinfo);
538 set_warning_line(wi, "unknown", -1);
539 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
540 warning_note(wi, warn_buf);
543 if (opts->bOrire == FALSE)
545 i = rm_interactions(F_ORIRES, nrmols, molinfo);
548 set_warning_line(wi, "unknown", -1);
549 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
550 warning_note(wi, warn_buf);
554 /* Copy structures from msys to sys */
555 molinfo2mtop(nrmols, molinfo, sys);
557 gmx_mtop_finalize(sys);
559 /* COORDINATE file processing */
562 fprintf(stderr, "processing coordinates...\n");
565 get_stx_coordnum(confin, &state->natoms);
566 if (state->natoms != sys->natoms)
568 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
569 " does not match topology (%s, %d)",
570 confin, state->natoms, topfile, sys->natoms);
574 /* make space for coordinates and velocities */
577 init_t_atoms(confat, state->natoms, FALSE);
578 init_state(state, state->natoms, 0, 0, 0, 0);
579 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
580 /* This call fixes the box shape for runs with pressure scaling */
581 set_box_rel(ir, state);
583 nmismatch = check_atom_names(topfile, confin, sys, confat);
584 free_t_atoms(confat, TRUE);
589 sprintf(buf, "%d non-matching atom name%s\n"
590 "atom names from %s will be used\n"
591 "atom names from %s will be ignored\n",
592 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
597 fprintf(stderr, "double-checking input for internal consistency...\n");
599 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
605 gmx_mtop_atomloop_all_t aloop;
609 snew(mass, state->natoms);
610 aloop = gmx_mtop_atomloop_all_init(sys);
611 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
617 if (opts->seed == -1)
619 useed = (int)gmx_rng_make_seed();
620 fprintf(stderr, "Setting gen_seed to %u\n", useed);
622 maxwell_speed(opts->tempi, useed, sys, state->v);
624 stop_cm(stdout, state->natoms, mass, state->x, state->v);
632 static void copy_state(const char *slog, t_trxframe *fr,
633 gmx_bool bReadVel, t_state *state,
638 if (fr->not_ok & FRAME_NOT_OK)
640 gmx_fatal(FARGS, "Can not start from an incomplete frame");
644 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
648 for (i = 0; i < state->natoms; i++)
650 copy_rvec(fr->x[i], state->x[i]);
656 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
658 for (i = 0; i < state->natoms; i++)
660 copy_rvec(fr->v[i], state->v[i]);
665 copy_mat(fr->box, state->box);
668 *use_time = fr->time;
671 static void cont_status(const char *slog, const char *ener,
672 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
673 t_inputrec *ir, t_state *state,
675 const output_env_t oenv)
676 /* If fr_time == -1 read the last frame available which is complete */
684 bReadVel = (bNeedVel && !bGenVel);
687 "Reading Coordinates%s and Box size from old trajectory\n",
688 bReadVel ? ", Velocities" : "");
691 fprintf(stderr, "Will read whole trajectory\n");
695 fprintf(stderr, "Will read till time %g\n", fr_time);
701 fprintf(stderr, "Velocities generated: "
702 "ignoring velocities in input trajectory\n");
704 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
708 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
714 "WARNING: Did not find a frame with velocities in file %s,\n"
715 " all velocities will be set to zero!\n\n", slog);
716 for (i = 0; i < sys->natoms; i++)
718 clear_rvec(state->v[i]);
721 /* Search for a frame without velocities */
723 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
727 state->natoms = fr.natoms;
729 if (sys->natoms != state->natoms)
731 gmx_fatal(FARGS, "Number of atoms in Topology "
732 "is not the same as in Trajectory");
734 copy_state(slog, &fr, bReadVel, state, &use_time);
736 /* Find the appropriate frame */
737 while ((fr_time == -1 || fr.time < fr_time) &&
738 read_next_frame(oenv, fp, &fr))
740 copy_state(slog, &fr, bReadVel, state, &use_time);
745 /* Set the relative box lengths for preserving the box shape.
746 * Note that this call can lead to differences in the last bit
747 * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
749 set_box_rel(ir, state);
751 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
752 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
754 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
756 get_enx_state(ener, use_time, &sys->groups, ir, state);
757 preserve_box_shape(ir, state->box_rel, state->boxv);
761 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
763 int rc_scaling, int ePBC,
767 gmx_bool bFirst = TRUE, *hadAtom;
773 int natoms, npbcdim = 0;
774 char warn_buf[STRLEN], title[STRLEN];
775 int a, i, ai, j, k, mb, nat_molb;
776 gmx_molblock_t *molb;
780 get_stx_coordnum(fn, &natoms);
781 if (natoms != mtop->natoms)
783 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, min(mtop->natoms, natoms), fn);
784 warning(wi, warn_buf);
788 init_t_atoms(&dumat, natoms, FALSE);
789 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
791 npbcdim = ePBC2npbcdim(ePBC);
793 if (rc_scaling != erscNO)
795 copy_mat(box, invbox);
796 for (j = npbcdim; j < DIM; j++)
798 clear_rvec(invbox[j]);
801 m_inv_ur0(invbox, invbox);
804 /* Copy the reference coordinates to mtop */
808 snew(hadAtom, natoms);
809 for (mb = 0; mb < mtop->nmolblock; mb++)
811 molb = &mtop->molblock[mb];
812 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
813 pr = &(molinfo[molb->type].plist[F_POSRES]);
814 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
815 if (pr->nr > 0 || prfb->nr > 0)
817 atom = mtop->moltype[molb->type].atoms.atom;
818 for (i = 0; (i < pr->nr); i++)
820 ai = pr->param[i].AI;
823 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
824 ai+1, *molinfo[molb->type].name, fn, natoms);
827 if (rc_scaling == erscCOM)
829 /* Determine the center of mass of the posres reference coordinates */
830 for (j = 0; j < npbcdim; j++)
832 sum[j] += atom[ai].m*x[a+ai][j];
834 totmass += atom[ai].m;
837 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
838 for (i = 0; (i < prfb->nr); i++)
840 ai = prfb->param[i].AI;
843 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
844 ai+1, *molinfo[molb->type].name, fn, natoms);
846 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
848 /* Determine the center of mass of the posres reference coordinates */
849 for (j = 0; j < npbcdim; j++)
851 sum[j] += atom[ai].m*x[a+ai][j];
853 totmass += atom[ai].m;
858 molb->nposres_xA = nat_molb;
859 snew(molb->posres_xA, molb->nposres_xA);
860 for (i = 0; i < nat_molb; i++)
862 copy_rvec(x[a+i], molb->posres_xA[i]);
867 molb->nposres_xB = nat_molb;
868 snew(molb->posres_xB, molb->nposres_xB);
869 for (i = 0; i < nat_molb; i++)
871 copy_rvec(x[a+i], molb->posres_xB[i]);
877 if (rc_scaling == erscCOM)
881 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
883 for (j = 0; j < npbcdim; j++)
885 com[j] = sum[j]/totmass;
887 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]);
890 if (rc_scaling != erscNO)
892 assert(npbcdim <= DIM);
894 for (mb = 0; mb < mtop->nmolblock; mb++)
896 molb = &mtop->molblock[mb];
897 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
898 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
900 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
901 for (i = 0; i < nat_molb; i++)
903 for (j = 0; j < npbcdim; j++)
905 if (rc_scaling == erscALL)
907 /* Convert from Cartesian to crystal coordinates */
908 xp[i][j] *= invbox[j][j];
909 for (k = j+1; k < npbcdim; k++)
911 xp[i][j] += invbox[k][j]*xp[i][k];
914 else if (rc_scaling == erscCOM)
916 /* Subtract the center of mass */
924 if (rc_scaling == erscCOM)
926 /* Convert the COM from Cartesian to crystal coordinates */
927 for (j = 0; j < npbcdim; j++)
929 com[j] *= invbox[j][j];
930 for (k = j+1; k < npbcdim; k++)
932 com[j] += invbox[k][j]*com[k];
938 free_t_atoms(&dumat, TRUE);
944 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
945 char *fnA, char *fnB,
946 int rc_scaling, int ePBC,
952 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
953 if (strcmp(fnA, fnB) != 0)
955 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
959 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
960 t_inputrec *ir, warninp_t wi)
963 char warn_buf[STRLEN];
967 fprintf(stderr, "Searching the wall atom type(s)\n");
969 for (i = 0; i < ir->nwall; i++)
971 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
972 if (ir->wall_atomtype[i] == NOTSET)
974 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
975 warning_error(wi, warn_buf);
980 static int nrdf_internal(t_atoms *atoms)
985 for (i = 0; i < atoms->nr; i++)
987 /* Vsite ptype might not be set here yet, so also check the mass */
988 if ((atoms->atom[i].ptype == eptAtom ||
989 atoms->atom[i].ptype == eptNucleus)
990 && atoms->atom[i].m > 0)
997 case 0: nrdf = 0; break;
998 case 1: nrdf = 0; break;
999 case 2: nrdf = 1; break;
1000 default: nrdf = nmass*3 - 6; break;
1007 spline1d( double dx,
1019 for (i = 1; i < n-1; i++)
1021 p = 0.5*y2[i-1]+2.0;
1023 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1024 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1029 for (i = n-2; i >= 0; i--)
1031 y2[i] = y2[i]*y2[i+1]+u[i];
1037 interpolate1d( double xmin,
1050 a = (xmin+(ix+1)*dx-x)/dx;
1051 b = (x-xmin-ix*dx)/dx;
1053 *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;
1054 *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];
1059 setup_cmap (int grid_spacing,
1062 gmx_cmap_t * cmap_grid)
1064 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1066 int i, j, k, ii, jj, kk, idx;
1068 double dx, xmin, v, v1, v2, v12;
1071 snew(tmp_u, 2*grid_spacing);
1072 snew(tmp_u2, 2*grid_spacing);
1073 snew(tmp_yy, 2*grid_spacing);
1074 snew(tmp_y1, 2*grid_spacing);
1075 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1076 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1078 dx = 360.0/grid_spacing;
1079 xmin = -180.0-dx*grid_spacing/2;
1081 for (kk = 0; kk < nc; kk++)
1083 /* Compute an offset depending on which cmap we are using
1084 * Offset will be the map number multiplied with the
1085 * grid_spacing * grid_spacing * 2
1087 offset = kk * grid_spacing * grid_spacing * 2;
1089 for (i = 0; i < 2*grid_spacing; i++)
1091 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1093 for (j = 0; j < 2*grid_spacing; j++)
1095 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1096 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1100 for (i = 0; i < 2*grid_spacing; i++)
1102 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1105 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1107 ii = i-grid_spacing/2;
1110 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1112 jj = j-grid_spacing/2;
1115 for (k = 0; k < 2*grid_spacing; k++)
1117 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1118 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1121 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1122 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1123 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1124 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1126 idx = ii*grid_spacing+jj;
1127 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1128 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1129 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1130 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1136 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1140 cmap_grid->ngrid = ngrid;
1141 cmap_grid->grid_spacing = grid_spacing;
1142 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1144 snew(cmap_grid->cmapdata, ngrid);
1146 for (i = 0; i < cmap_grid->ngrid; i++)
1148 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1153 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1155 int count, count_mol, i, mb;
1156 gmx_molblock_t *molb;
1161 for (mb = 0; mb < mtop->nmolblock; mb++)
1164 molb = &mtop->molblock[mb];
1165 plist = mi[molb->type].plist;
1167 for (i = 0; i < F_NRE; i++)
1171 count_mol += 3*plist[i].nr;
1173 else if (interaction_function[i].flags & IF_CONSTRAINT)
1175 count_mol += plist[i].nr;
1179 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1182 "Molecule type '%s' has %d constraints.\n"
1183 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1184 *mi[molb->type].name, count_mol,
1185 nrdf_internal(&mi[molb->type].atoms));
1188 count += molb->nmol*count_mol;
1194 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1196 int i, nmiss, natoms, mt;
1198 const t_atoms *atoms;
1201 for (mt = 0; mt < sys->nmoltype; mt++)
1203 atoms = &sys->moltype[mt].atoms;
1206 for (i = 0; i < natoms; i++)
1208 q = atoms->atom[i].q;
1209 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1210 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1211 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1212 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1213 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1216 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1217 get_atomtype_name(atoms->atom[i].type, atype), q);
1225 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);
1230 static void check_gbsa_params(gpp_atomtype_t atype)
1234 /* If we are doing GBSA, check that we got the parameters we need
1235 * This checking is to see if there are GBSA paratmeters for all
1236 * atoms in the force field. To go around this for testing purposes
1237 * comment out the nerror++ counter temporarily
1240 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1242 if (get_atomtype_radius(i, atype) < 0 ||
1243 get_atomtype_vol(i, atype) < 0 ||
1244 get_atomtype_surftens(i, atype) < 0 ||
1245 get_atomtype_gb_radius(i, atype) < 0 ||
1246 get_atomtype_S_hct(i, atype) < 0)
1248 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1249 get_atomtype_name(i, atype));
1256 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);
1261 static real calc_temp(const gmx_mtop_t *mtop,
1262 const t_inputrec *ir,
1266 gmx_mtop_atomloop_all_t aloop;
1273 aloop = gmx_mtop_atomloop_all_init(mtop);
1274 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1276 sum_mv2 += atom->m*norm2(v[a]);
1280 for (g = 0; g < ir->opts.ngtc; g++)
1282 nrdf += ir->opts.nrdf[g];
1285 return sum_mv2/(nrdf*BOLTZ);
1288 static real get_max_reference_temp(const t_inputrec *ir,
1297 for (i = 0; i < ir->opts.ngtc; i++)
1299 if (ir->opts.tau_t[i] < 0)
1305 ref_t = max(ref_t, ir->opts.ref_t[i]);
1313 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.",
1321 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1328 verletbuf_list_setup_t ls;
1331 char warn_buf[STRLEN];
1333 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1335 /* Calculate the buffer size for simple atom vs atoms list */
1336 ls.cluster_size_i = 1;
1337 ls.cluster_size_j = 1;
1338 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1339 &ls, &n_nonlin_vsite, &rlist_1x1);
1341 /* Set the pair-list buffer size in ir */
1342 verletbuf_get_list_setup(FALSE, &ls);
1343 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1344 &ls, &n_nonlin_vsite, &ir->rlist);
1346 if (n_nonlin_vsite > 0)
1348 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);
1349 warning_note(wi, warn_buf);
1352 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1353 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1355 ir->rlistlong = ir->rlist;
1356 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1357 ls.cluster_size_i, ls.cluster_size_j,
1358 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1360 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1362 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->rlistlong, sqrt(max_cutoff2(ir->ePBC, box)));
1366 int gmx_grompp(int argc, char *argv[])
1368 static const char *desc[] = {
1369 "[THISMODULE] (the gromacs preprocessor)",
1370 "reads a molecular topology file, checks the validity of the",
1371 "file, expands the topology from a molecular description to an atomic",
1372 "description. The topology file contains information about",
1373 "molecule types and the number of molecules, the preprocessor",
1374 "copies each molecule as needed. ",
1375 "There is no limitation on the number of molecule types. ",
1376 "Bonds and bond-angles can be converted into constraints, separately",
1377 "for hydrogens and heavy atoms.",
1378 "Then a coordinate file is read and velocities can be generated",
1379 "from a Maxwellian distribution if requested.",
1380 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1381 "(eg. number of MD steps, time step, cut-off), and others such as",
1382 "NEMD parameters, which are corrected so that the net acceleration",
1384 "Eventually a binary file is produced that can serve as the sole input",
1385 "file for the MD program.[PAR]",
1387 "[THISMODULE] uses the atom names from the topology file. The atom names",
1388 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1389 "warnings when they do not match the atom names in the topology.",
1390 "Note that the atom names are irrelevant for the simulation as",
1391 "only the atom types are used for generating interaction parameters.[PAR]",
1393 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1394 "etc. The preprocessor supports the following keywords:[PAR]",
1395 "#ifdef VARIABLE[BR]",
1396 "#ifndef VARIABLE[BR]",
1399 "#define VARIABLE[BR]",
1400 "#undef VARIABLE[BR]"
1401 "#include \"filename\"[BR]",
1402 "#include <filename>[PAR]",
1403 "The functioning of these statements in your topology may be modulated by",
1404 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1405 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1406 "include = -I/home/john/doe[tt][BR]",
1407 "For further information a C-programming textbook may help you out.",
1408 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1409 "topology file written out so that you can verify its contents.[PAR]",
1411 /* cpp has been unnecessary for some time, hasn't it?
1412 "If your system does not have a C-preprocessor, you can still",
1413 "use [TT]grompp[tt], but you do not have access to the features ",
1414 "from the cpp. Command line options to the C-preprocessor can be given",
1415 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1418 "When using position restraints a file with restraint coordinates",
1419 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1420 "with respect to the conformation from the [TT]-c[tt] option.",
1421 "For free energy calculation the the coordinates for the B topology",
1422 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1423 "those of the A topology.[PAR]",
1425 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1426 "The last frame with coordinates and velocities will be read,",
1427 "unless the [TT]-time[tt] option is used. Only if this information",
1428 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1429 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1430 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1431 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1434 "[THISMODULE] can be used to restart simulations (preserving",
1435 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1436 "However, for simply changing the number of run steps to extend",
1437 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1438 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1439 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1440 "like output frequency, then supplying the checkpoint file to",
1441 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1442 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1444 "By default, all bonded interactions which have constant energy due to",
1445 "virtual site constructions will be removed. If this constant energy is",
1446 "not zero, this will result in a shift in the total energy. All bonded",
1447 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1448 "all constraints for distances which will be constant anyway because",
1449 "of virtual site constructions will be removed. If any constraints remain",
1450 "which involve virtual sites, a fatal error will result.[PAR]"
1452 "To verify your run input file, please take note of all warnings",
1453 "on the screen, and correct where necessary. Do also look at the contents",
1454 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1455 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1456 "with the [TT]-debug[tt] option which will give you more information",
1457 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1458 "can see the contents of the run input file with the [gmx-dump]",
1459 "program. [gmx-check] can be used to compare the contents of two",
1460 "run input files.[PAR]"
1462 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1463 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1464 "harmless, but usually they are not. The user is advised to carefully",
1465 "interpret the output messages before attempting to bypass them with",
1472 gpp_atomtype_t atype;
1474 int natoms, nvsite, comb, mt;
1478 real max_spacing, fudgeQQ;
1480 char fn[STRLEN], fnB[STRLEN];
1481 const char *mdparin;
1483 gmx_bool bNeedVel, bGenVel;
1484 gmx_bool have_atomnumber;
1486 t_params *gb_plist = NULL;
1487 gmx_genborn_t *born = NULL;
1489 gmx_bool bVerbose = FALSE;
1491 char warn_buf[STRLEN];
1493 t_atoms IMDatoms; /* Atoms to be operated on interactively (IMD) */
1496 { efMDP, NULL, NULL, ffREAD },
1497 { efMDP, "-po", "mdout", ffWRITE },
1498 { efSTX, "-c", NULL, ffREAD },
1499 { efSTX, "-r", NULL, ffOPTRD },
1500 { efSTX, "-rb", NULL, ffOPTRD },
1501 { efNDX, NULL, NULL, ffOPTRD },
1502 { efTOP, NULL, NULL, ffREAD },
1503 { efTOP, "-pp", "processed", ffOPTWR },
1504 { efTPX, "-o", NULL, ffWRITE },
1505 { efTRN, "-t", NULL, ffOPTRD },
1506 { efEDR, "-e", NULL, ffOPTRD },
1507 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1508 { efGRO, "-imd", "imdgroup", ffOPTWR },
1509 { efTRN, "-ref", "rotref", ffOPTRW }
1511 #define NFILE asize(fnm)
1513 /* Command line options */
1514 static gmx_bool bRenum = TRUE;
1515 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1516 static int i, maxwarn = 0;
1517 static real fr_time = -1;
1519 { "-v", FALSE, etBOOL, {&bVerbose},
1520 "Be loud and noisy" },
1521 { "-time", FALSE, etREAL, {&fr_time},
1522 "Take frame at or first after this time." },
1523 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1524 "Remove constant bonded interactions with virtual sites" },
1525 { "-maxwarn", FALSE, etINT, {&maxwarn},
1526 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1527 { "-zero", FALSE, etBOOL, {&bZero},
1528 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1529 { "-renum", FALSE, etBOOL, {&bRenum},
1530 "Renumber atomtypes and minimize number of atomtypes" }
1533 /* Initiate some variables */
1538 /* Parse the command line */
1539 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1540 asize(desc), desc, 0, NULL, &oenv))
1545 wi = init_warning(TRUE, maxwarn);
1547 /* PARAMETER file processing */
1548 mdparin = opt2fn("-f", NFILE, fnm);
1549 set_warning_line(wi, mdparin, -1);
1550 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1554 fprintf(stderr, "checking input for internal consistency...\n");
1556 check_ir(mdparin, ir, opts, wi);
1558 if (ir->ld_seed == -1)
1560 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
1561 fprintf(stderr, "Setting the LD random seed to %"GMX_PRId64 "\n", ir->ld_seed);
1564 if (ir->expandedvals->lmc_seed == -1)
1566 ir->expandedvals->lmc_seed = (int)gmx_rng_make_seed();
1567 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1570 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1571 bGenVel = (bNeedVel && opts->bGenVel);
1572 if (bGenVel && ir->bContinuation)
1575 "Generating velocities is inconsistent with attempting "
1576 "to continue a previous run. Choose only one of "
1577 "gen-vel = yes and continuation = yes.");
1578 warning_error(wi, warn_buf);
1584 atype = init_atomtype();
1587 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1590 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1591 if (!gmx_fexist(fn))
1593 gmx_fatal(FARGS, "%s does not exist", fn);
1595 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1596 opts, ir, bZero, bGenVel, bVerbose, &state,
1597 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1603 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1607 /* set parameters for virtual site construction (not for vsiten) */
1608 for (mt = 0; mt < sys->nmoltype; mt++)
1611 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1613 /* now throw away all obsolete bonds, angles and dihedrals: */
1614 /* note: constraints are ALWAYS removed */
1617 for (mt = 0; mt < sys->nmoltype; mt++)
1619 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1623 if (ir->cutoff_scheme == ecutsVERLET)
1625 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1626 ecutscheme_names[ir->cutoff_scheme]);
1628 /* Remove all charge groups */
1629 gmx_mtop_remove_chargegroups(sys);
1632 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1634 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1636 sprintf(warn_buf, "Can not do %s with %s, use %s",
1637 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1638 warning_error(wi, warn_buf);
1640 if (ir->bPeriodicMols)
1642 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1643 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1644 warning_error(wi, warn_buf);
1648 if (EI_SD (ir->eI) && ir->etc != etcNO)
1650 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1653 /* If we are doing QM/MM, check that we got the atom numbers */
1654 have_atomnumber = TRUE;
1655 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1657 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1659 if (!have_atomnumber && ir->bQMMM)
1663 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1664 "field you are using does not contain atom numbers fields. This is an\n"
1665 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1666 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1667 "an integer just before the mass column in ffXXXnb.itp.\n"
1668 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1673 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1675 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1677 /** TODO check size of ex+hy width against box size */
1680 /* Check for errors in the input now, since they might cause problems
1681 * during processing further down.
1683 check_warning_error(wi, FARGS);
1685 if (opt2bSet("-r", NFILE, fnm))
1687 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1691 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1693 if (opt2bSet("-rb", NFILE, fnm))
1695 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1702 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1706 fprintf(stderr, "Reading position restraint coords from %s", fn);
1707 if (strcmp(fn, fnB) == 0)
1709 fprintf(stderr, "\n");
1713 fprintf(stderr, " and %s\n", fnB);
1716 gen_posres(sys, mi, fn, fnB,
1717 ir->refcoord_scaling, ir->ePBC,
1718 ir->posres_com, ir->posres_comB,
1722 /* If we are using CMAP, setup the pre-interpolation grid */
1723 if (plist->ncmap > 0)
1725 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1726 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1729 set_wall_atomtype(atype, opts, ir, wi);
1732 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1733 ntype = get_atomtype_ntypes(atype);
1736 if (ir->implicit_solvent != eisNO)
1738 /* Now we have renumbered the atom types, we can check the GBSA params */
1739 check_gbsa_params(atype);
1741 /* Check that all atoms that have charge and/or LJ-parameters also have
1742 * sensible GB-parameters
1744 check_gbsa_params_charged(sys, atype);
1747 /* PELA: Copy the atomtype data to the topology atomtype list */
1748 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1752 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1757 fprintf(stderr, "converting bonded parameters...\n");
1760 ntype = get_atomtype_ntypes(atype);
1761 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1765 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1768 /* set ptype to VSite for virtual sites */
1769 for (mt = 0; mt < sys->nmoltype; mt++)
1771 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1775 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1777 /* Check velocity for virtual sites and shells */
1780 check_vel(sys, state.v);
1786 for (i = 0; i < sys->nmoltype; i++)
1788 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1791 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1793 check_bonds_timestep(sys, ir->delta_t, wi);
1796 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1798 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.");
1801 check_warning_error(wi, FARGS);
1805 fprintf(stderr, "initialising group options...\n");
1807 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1809 bGenVel ? state.v : NULL,
1812 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1815 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
1819 if (EI_MD(ir->eI) && ir->etc == etcNO)
1823 buffer_temp = opts->tempi;
1827 buffer_temp = calc_temp(sys, ir, state.v);
1829 if (buffer_temp > 0)
1831 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
1832 warning_note(wi, warn_buf);
1836 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
1837 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
1838 warning_note(wi, warn_buf);
1843 buffer_temp = get_max_reference_temp(ir, wi);
1846 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
1848 /* NVE with initial T=0: we add a fixed ratio to rlist.
1849 * Since we don't actually use verletbuf_tol, we set it to -1
1850 * so it can't be misused later.
1852 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
1853 ir->verletbuf_tol = -1;
1857 /* We warn for NVE simulations with >1(.1)% drift tolerance */
1858 const real drift_tol = 0.01;
1861 /* We use 2 DOF per atom = 2kT pot+kin energy, to be on
1862 * the safe side with constraints (without constraints: 3 DOF).
1864 ener_runtime = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
1866 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
1868 ir->verletbuf_tol > 1.1*drift_tol*ener_runtime)
1870 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%%, you might need to set verlet-buffer-tolerance to %.1e.",
1871 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
1872 (int)(ir->verletbuf_tol/ener_runtime*100 + 0.5),
1873 (int)(100*drift_tol + 0.5),
1874 drift_tol*ener_runtime);
1875 warning_note(wi, warn_buf);
1878 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
1883 /* Init the temperature coupling state */
1884 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1888 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1890 check_eg_vs_cg(sys);
1894 pr_symtab(debug, 0, "After index", &sys->symtab);
1897 triple_check(mdparin, ir, sys, wi);
1898 close_symtab(&sys->symtab);
1901 pr_symtab(debug, 0, "After close", &sys->symtab);
1904 /* make exclusions between QM atoms */
1907 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1909 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1913 generate_qmexcl(sys, ir, wi);
1917 if (ftp2bSet(efTRN, NFILE, fnm))
1921 fprintf(stderr, "getting data from old trajectory ...\n");
1923 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1924 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1927 if (ir->ePBC == epbcXY && ir->nwall != 2)
1929 clear_rvec(state.box[ZZ]);
1932 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1934 set_warning_line(wi, mdparin, -1);
1935 check_chargegroup_radii(sys, ir, state.x, wi);
1938 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1940 /* Calculate the optimal grid dimensions */
1941 copy_mat(state.box, box);
1942 if (ir->ePBC == epbcXY && ir->nwall == 2)
1944 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1946 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1948 /* Mark fourier_spacing as not used */
1949 ir->fourier_spacing = 0;
1951 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1953 set_warning_line(wi, mdparin, -1);
1954 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1956 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1957 &(ir->nkx), &(ir->nky), &(ir->nkz));
1960 /* MRS: eventually figure out better logic for initializing the fep
1961 values that makes declaring the lambda and declaring the state not
1962 potentially conflict if not handled correctly. */
1963 if (ir->efep != efepNO)
1965 state.fep_state = ir->fepvals->init_fep_state;
1966 for (i = 0; i < efptNR; i++)
1968 /* init_lambda trumps state definitions*/
1969 if (ir->fepvals->init_lambda >= 0)
1971 state.lambda[i] = ir->fepvals->init_lambda;
1975 if (ir->fepvals->all_lambda[i] == NULL)
1977 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1981 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1987 if (ir->ePull != epullNO)
1989 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1994 set_reference_positions(ir->rot, state.x, state.box,
1995 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1999 /* reset_multinr(sys); */
2001 if (EEL_PME(ir->coulombtype))
2003 float ratio = pme_load_estimate(sys, ir, state.box);
2004 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2005 /* With free energy we might need to do PME both for the A and B state
2006 * charges. This will double the cost, but the optimal performance will
2007 * then probably be at a slightly larger cut-off and grid spacing.
2009 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2010 (ir->efep != efepNO && ratio > 2.0/3.0))
2013 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2014 "and for highly parallel simulations between 0.25 and 0.33,\n"
2015 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2016 if (ir->efep != efepNO)
2019 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2025 char warn_buf[STRLEN];
2026 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2027 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2030 set_warning_line(wi, mdparin, -1);
2031 warning_note(wi, warn_buf);
2035 printf("%s\n", warn_buf);
2041 fprintf(stderr, "writing run input file...\n");
2044 done_warning(wi, FARGS);
2045 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);
2047 /* Output IMD group, if bIMD is TRUE */
2048 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2050 done_atomtype(atype);
2051 done_mtop(sys, TRUE);
2052 done_inputrec_strings();