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>
56 #include "gromacs/fileio/confio.h"
60 #include "grompp-impl.h"
61 #include "gromacs/random/random.h"
62 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
64 #include "gromacs/fileio/futil.h"
65 #include "gromacs/commandline/pargs.h"
67 #include "gromacs/gmxana/sortwater.h"
69 #include "gmx_fatal.h"
72 #include "gromacs/fileio/gmxfio.h"
73 #include "gromacs/fileio/trnio.h"
74 #include "gromacs/fileio/tpxio.h"
75 #include "gromacs/fileio/trxio.h"
76 #include "vsite_parm.h"
80 #include "gromacs/fileio/enxio.h"
82 #include "compute_io.h"
83 #include "gpp_atomtype.h"
84 #include "mtop_util.h"
86 #include "calc_verletbuf.h"
90 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
95 /* For all the molecule types */
96 for (i = 0; i < nrmols; i++)
98 n += mols[i].plist[ifunc].nr;
99 mols[i].plist[ifunc].nr = 0;
104 static int check_atom_names(const char *fn1, const char *fn2,
105 gmx_mtop_t *mtop, t_atoms *at)
107 int mb, m, i, j, nmismatch;
109 #define MAXMISMATCH 20
111 if (mtop->natoms != at->nr)
113 gmx_incons("comparing atom names");
118 for (mb = 0; mb < mtop->nmolblock; mb++)
120 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
121 for (m = 0; m < mtop->molblock[mb].nmol; m++)
123 for (j = 0; j < tat->nr; j++)
125 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
127 if (nmismatch < MAXMISMATCH)
130 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
131 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
133 else if (nmismatch == MAXMISMATCH)
135 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
147 static void check_eg_vs_cg(gmx_mtop_t *mtop)
149 int astart, mb, m, cg, j, firstj;
150 unsigned char firsteg, eg;
153 /* Go through all the charge groups and make sure all their
154 * atoms are in the same energy group.
158 for (mb = 0; mb < mtop->nmolblock; mb++)
160 molt = &mtop->moltype[mtop->molblock[mb].type];
161 for (m = 0; m < mtop->molblock[mb].nmol; m++)
163 for (cg = 0; cg < molt->cgs.nr; cg++)
165 /* Get the energy group of the first atom in this charge group */
166 firstj = astart + molt->cgs.index[cg];
167 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
168 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
170 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
173 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
174 firstj+1, astart+j+1, cg+1, *molt->name);
178 astart += molt->atoms.nr;
183 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
186 char warn_buf[STRLEN];
189 for (cg = 0; cg < cgs->nr; cg++)
191 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
194 if (maxsize > MAX_CHARGEGROUP_SIZE)
196 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
198 else if (maxsize > 10)
200 set_warning_line(wi, topfn, -1);
202 "The largest charge group contains %d atoms.\n"
203 "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"
204 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
205 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
207 warning_note(wi, warn_buf);
211 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
213 /* This check is not intended to ensure accurate integration,
214 * rather it is to signal mistakes in the mdp settings.
215 * A common mistake is to forget to turn on constraints
216 * for MD after energy minimization with flexible bonds.
217 * This check can also detect too large time steps for flexible water
218 * models, but such errors will often be masked by the constraints
219 * mdp options, which turns flexible water into water with bond constraints,
220 * but without an angle constraint. Unfortunately such incorrect use
221 * of water models can not easily be detected without checking
222 * for specific model names.
224 * The stability limit of leap-frog or velocity verlet is 4.44 steps
225 * per oscillational period.
226 * But accurate bonds distributions are lost far before that limit.
227 * To allow relatively common schemes (although not common with Gromacs)
228 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
229 * we set the note limit to 10.
231 int min_steps_warn = 5;
232 int min_steps_note = 10;
235 gmx_moltype_t *moltype, *w_moltype;
237 t_ilist *ilist, *ilb, *ilc, *ils;
239 int i, a1, a2, w_a1, w_a2, j;
240 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
241 gmx_bool bFound, bWater, bWarn;
242 char warn_buf[STRLEN];
244 ip = mtop->ffparams.iparams;
246 twopi2 = sqr(2*M_PI);
248 limit2 = sqr(min_steps_note*dt);
254 for (molt = 0; molt < mtop->nmoltype; molt++)
256 moltype = &mtop->moltype[molt];
257 atom = moltype->atoms.atom;
258 ilist = moltype->ilist;
259 ilc = &ilist[F_CONSTR];
260 ils = &ilist[F_SETTLE];
261 for (ftype = 0; ftype < F_NRE; ftype++)
263 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
269 for (i = 0; i < ilb->nr; i += 3)
271 fc = ip[ilb->iatoms[i]].harmonic.krA;
272 re = ip[ilb->iatoms[i]].harmonic.rA;
273 if (ftype == F_G96BONDS)
275 /* Convert squared sqaure fc to harmonic fc */
278 a1 = ilb->iatoms[i+1];
279 a2 = ilb->iatoms[i+2];
282 if (fc > 0 && m1 > 0 && m2 > 0)
284 period2 = twopi2*m1*m2/((m1 + m2)*fc);
288 period2 = GMX_FLOAT_MAX;
292 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
293 fc, m1, m2, sqrt(period2));
295 if (period2 < limit2)
298 for (j = 0; j < ilc->nr; j += 3)
300 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
301 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
306 for (j = 0; j < ils->nr; j += 4)
308 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
309 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
315 (w_moltype == NULL || period2 < w_period2))
327 if (w_moltype != NULL)
329 bWarn = (w_period2 < sqr(min_steps_warn*dt));
330 /* A check that would recognize most water models */
331 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
332 w_moltype->atoms.nr <= 5);
333 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"
336 w_a1+1, *w_moltype->atoms.atomname[w_a1],
337 w_a2+1, *w_moltype->atoms.atomname[w_a2],
338 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
340 "Maybe you asked for fexible water." :
341 "Maybe you forgot to change the constraints mdp option.");
344 warning(wi, warn_buf);
348 warning_note(wi, warn_buf);
353 static void check_vel(gmx_mtop_t *mtop, rvec v[])
355 gmx_mtop_atomloop_all_t aloop;
359 aloop = gmx_mtop_atomloop_all_init(mtop);
360 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
362 if (atom->ptype == eptShell ||
363 atom->ptype == eptBond ||
364 atom->ptype == eptVSite)
371 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
376 for (mb = 0; mb < mtop->nmolblock; mb++)
378 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
384 /* This routine reorders the molecule type array
385 * in the order of use in the molblocks,
386 * unused molecule types are deleted.
388 static void renumber_moltypes(gmx_mtop_t *sys,
389 int *nmolinfo, t_molinfo **molinfo)
391 int *order, norder, i;
395 snew(order, *nmolinfo);
397 for (mb = 0; mb < sys->nmolblock; mb++)
399 for (i = 0; i < norder; i++)
401 if (order[i] == sys->molblock[mb].type)
408 /* This type did not occur yet, add it */
409 order[norder] = sys->molblock[mb].type;
410 /* Renumber the moltype in the topology */
413 sys->molblock[mb].type = i;
416 /* We still need to reorder the molinfo structs */
418 for (mi = 0; mi < *nmolinfo; mi++)
420 for (i = 0; i < norder; i++)
429 done_mi(&(*molinfo)[mi]);
433 minew[i] = (*molinfo)[mi];
442 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
447 mtop->nmoltype = nmi;
448 snew(mtop->moltype, nmi);
449 for (m = 0; m < nmi; m++)
451 molt = &mtop->moltype[m];
452 molt->name = mi[m].name;
453 molt->atoms = mi[m].atoms;
454 /* ilists are copied later */
455 molt->cgs = mi[m].cgs;
456 molt->excls = mi[m].excls;
461 new_status(const char *topfile, const char *topppfile, const char *confin,
462 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
463 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
464 gpp_atomtype_t atype, gmx_mtop_t *sys,
465 int *nmi, t_molinfo **mi, t_params plist[],
466 int *comb, double *reppow, real *fudgeQQ,
470 t_molinfo *molinfo = NULL;
472 gmx_molblock_t *molblock, *molbs;
474 int mb, i, nrmols, nmismatch;
476 gmx_bool bGB = FALSE;
477 char warn_buf[STRLEN];
481 /* Set gmx_boolean for GB */
482 if (ir->implicit_solvent)
487 /* TOPOLOGY processing */
488 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
489 plist, comb, reppow, fudgeQQ,
490 atype, &nrmols, &molinfo, ir,
491 &nmolblock, &molblock, bGB,
495 snew(sys->molblock, nmolblock);
498 for (mb = 0; mb < nmolblock; mb++)
500 if (sys->nmolblock > 0 &&
501 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
503 /* Merge consecutive blocks with the same molecule type */
504 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
505 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
507 else if (molblock[mb].nmol > 0)
509 /* Add a new molblock to the topology */
510 molbs = &sys->molblock[sys->nmolblock];
511 *molbs = molblock[mb];
512 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
513 molbs->nposres_xA = 0;
514 molbs->nposres_xB = 0;
515 sys->natoms += molbs->nmol*molbs->natoms_mol;
519 if (sys->nmolblock == 0)
521 gmx_fatal(FARGS, "No molecules were defined in the system");
524 renumber_moltypes(sys, &nrmols, &molinfo);
528 convert_harmonics(nrmols, molinfo, atype);
531 if (ir->eDisre == edrNone)
533 i = rm_interactions(F_DISRES, nrmols, molinfo);
536 set_warning_line(wi, "unknown", -1);
537 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
538 warning_note(wi, warn_buf);
541 if (opts->bOrire == FALSE)
543 i = rm_interactions(F_ORIRES, nrmols, molinfo);
546 set_warning_line(wi, "unknown", -1);
547 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
548 warning_note(wi, warn_buf);
552 /* Copy structures from msys to sys */
553 molinfo2mtop(nrmols, molinfo, sys);
555 gmx_mtop_finalize(sys);
557 /* COORDINATE file processing */
560 fprintf(stderr, "processing coordinates...\n");
563 get_stx_coordnum(confin, &state->natoms);
564 if (state->natoms != sys->natoms)
566 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
567 " does not match topology (%s, %d)",
568 confin, state->natoms, topfile, sys->natoms);
572 /* make space for coordinates and velocities */
575 init_t_atoms(confat, state->natoms, FALSE);
576 init_state(state, state->natoms, 0, 0, 0, 0);
577 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
578 /* This call fixes the box shape for runs with pressure scaling */
579 set_box_rel(ir, state);
581 nmismatch = check_atom_names(topfile, confin, sys, confat);
582 free_t_atoms(confat, TRUE);
587 sprintf(buf, "%d non-matching atom name%s\n"
588 "atom names from %s will be used\n"
589 "atom names from %s will be ignored\n",
590 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
595 fprintf(stderr, "double-checking input for internal consistency...\n");
597 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
603 gmx_mtop_atomloop_all_t aloop;
607 snew(mass, state->natoms);
608 aloop = gmx_mtop_atomloop_all_init(sys);
609 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
615 if (opts->seed == -1)
617 useed = (int)gmx_rng_make_seed();
618 fprintf(stderr, "Setting gen_seed to %u\n", useed);
620 maxwell_speed(opts->tempi, useed, sys, state->v);
622 stop_cm(stdout, state->natoms, mass, state->x, state->v);
630 static void copy_state(const char *slog, t_trxframe *fr,
631 gmx_bool bReadVel, t_state *state,
636 if (fr->not_ok & FRAME_NOT_OK)
638 gmx_fatal(FARGS, "Can not start from an incomplete frame");
642 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
646 for (i = 0; i < state->natoms; i++)
648 copy_rvec(fr->x[i], state->x[i]);
654 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
656 for (i = 0; i < state->natoms; i++)
658 copy_rvec(fr->v[i], state->v[i]);
663 copy_mat(fr->box, state->box);
666 *use_time = fr->time;
669 static void cont_status(const char *slog, const char *ener,
670 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
671 t_inputrec *ir, t_state *state,
673 const output_env_t oenv)
674 /* If fr_time == -1 read the last frame available which is complete */
682 bReadVel = (bNeedVel && !bGenVel);
685 "Reading Coordinates%s and Box size from old trajectory\n",
686 bReadVel ? ", Velocities" : "");
689 fprintf(stderr, "Will read whole trajectory\n");
693 fprintf(stderr, "Will read till time %g\n", fr_time);
699 fprintf(stderr, "Velocities generated: "
700 "ignoring velocities in input trajectory\n");
702 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
706 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
712 "WARNING: Did not find a frame with velocities in file %s,\n"
713 " all velocities will be set to zero!\n\n", slog);
714 for (i = 0; i < sys->natoms; i++)
716 clear_rvec(state->v[i]);
719 /* Search for a frame without velocities */
721 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
725 state->natoms = fr.natoms;
727 if (sys->natoms != state->natoms)
729 gmx_fatal(FARGS, "Number of atoms in Topology "
730 "is not the same as in Trajectory");
732 copy_state(slog, &fr, bReadVel, state, &use_time);
734 /* Find the appropriate frame */
735 while ((fr_time == -1 || fr.time < fr_time) &&
736 read_next_frame(oenv, fp, &fr))
738 copy_state(slog, &fr, bReadVel, state, &use_time);
743 /* Set the relative box lengths for preserving the box shape.
744 * Note that this call can lead to differences in the last bit
745 * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
747 set_box_rel(ir, state);
749 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
750 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
752 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
754 get_enx_state(ener, use_time, &sys->groups, ir, state);
755 preserve_box_shape(ir, state->box_rel, state->boxv);
759 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
761 int rc_scaling, int ePBC,
765 gmx_bool bFirst = TRUE, *hadAtom;
771 int natoms, npbcdim = 0;
772 char warn_buf[STRLEN], title[STRLEN];
773 int a, i, ai, j, k, mb, nat_molb;
774 gmx_molblock_t *molb;
778 get_stx_coordnum(fn, &natoms);
779 if (natoms != mtop->natoms)
781 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);
782 warning(wi, warn_buf);
786 init_t_atoms(&dumat, natoms, FALSE);
787 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
789 npbcdim = ePBC2npbcdim(ePBC);
791 if (rc_scaling != erscNO)
793 copy_mat(box, invbox);
794 for (j = npbcdim; j < DIM; j++)
796 clear_rvec(invbox[j]);
799 m_inv_ur0(invbox, invbox);
802 /* Copy the reference coordinates to mtop */
806 snew(hadAtom, natoms);
807 for (mb = 0; mb < mtop->nmolblock; mb++)
809 molb = &mtop->molblock[mb];
810 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
811 pr = &(molinfo[molb->type].plist[F_POSRES]);
812 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
813 if (pr->nr > 0 || prfb->nr > 0)
815 atom = mtop->moltype[molb->type].atoms.atom;
816 for (i = 0; (i < pr->nr); i++)
818 ai = pr->param[i].AI;
821 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
822 ai+1, *molinfo[molb->type].name, fn, natoms);
825 if (rc_scaling == erscCOM)
827 /* Determine the center of mass of the posres reference coordinates */
828 for (j = 0; j < npbcdim; j++)
830 sum[j] += atom[ai].m*x[a+ai][j];
832 totmass += atom[ai].m;
835 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
836 for (i = 0; (i < prfb->nr); i++)
838 ai = prfb->param[i].AI;
841 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
842 ai+1, *molinfo[molb->type].name, fn, natoms);
844 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
846 /* Determine the center of mass of the posres reference coordinates */
847 for (j = 0; j < npbcdim; j++)
849 sum[j] += atom[ai].m*x[a+ai][j];
851 totmass += atom[ai].m;
856 molb->nposres_xA = nat_molb;
857 snew(molb->posres_xA, molb->nposres_xA);
858 for (i = 0; i < nat_molb; i++)
860 copy_rvec(x[a+i], molb->posres_xA[i]);
865 molb->nposres_xB = nat_molb;
866 snew(molb->posres_xB, molb->nposres_xB);
867 for (i = 0; i < nat_molb; i++)
869 copy_rvec(x[a+i], molb->posres_xB[i]);
875 if (rc_scaling == erscCOM)
879 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
881 for (j = 0; j < npbcdim; j++)
883 com[j] = sum[j]/totmass;
885 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]);
888 if (rc_scaling != erscNO)
890 for (mb = 0; mb < mtop->nmolblock; mb++)
892 molb = &mtop->molblock[mb];
893 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
894 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
896 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
897 for (i = 0; i < nat_molb; i++)
899 for (j = 0; j < npbcdim; j++)
901 if (rc_scaling == erscALL)
903 /* Convert from Cartesian to crystal coordinates */
904 xp[i][j] *= invbox[j][j];
905 for (k = j+1; k < npbcdim; k++)
907 xp[i][j] += invbox[k][j]*xp[i][k];
910 else if (rc_scaling == erscCOM)
912 /* Subtract the center of mass */
920 if (rc_scaling == erscCOM)
922 /* Convert the COM from Cartesian to crystal coordinates */
923 for (j = 0; j < npbcdim; j++)
925 com[j] *= invbox[j][j];
926 for (k = j+1; k < npbcdim; k++)
928 com[j] += invbox[k][j]*com[k];
934 free_t_atoms(&dumat, TRUE);
940 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
941 char *fnA, char *fnB,
942 int rc_scaling, int ePBC,
948 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
949 if (strcmp(fnA, fnB) != 0)
951 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
955 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
956 t_inputrec *ir, warninp_t wi)
959 char warn_buf[STRLEN];
963 fprintf(stderr, "Searching the wall atom type(s)\n");
965 for (i = 0; i < ir->nwall; i++)
967 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
968 if (ir->wall_atomtype[i] == NOTSET)
970 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
971 warning_error(wi, warn_buf);
976 static int nrdf_internal(t_atoms *atoms)
981 for (i = 0; i < atoms->nr; i++)
983 /* Vsite ptype might not be set here yet, so also check the mass */
984 if ((atoms->atom[i].ptype == eptAtom ||
985 atoms->atom[i].ptype == eptNucleus)
986 && atoms->atom[i].m > 0)
993 case 0: nrdf = 0; break;
994 case 1: nrdf = 0; break;
995 case 2: nrdf = 1; break;
996 default: nrdf = nmass*3 - 6; break;
1003 spline1d( double dx,
1015 for (i = 1; i < n-1; i++)
1017 p = 0.5*y2[i-1]+2.0;
1019 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1020 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1025 for (i = n-2; i >= 0; i--)
1027 y2[i] = y2[i]*y2[i+1]+u[i];
1033 interpolate1d( double xmin,
1046 a = (xmin+(ix+1)*dx-x)/dx;
1047 b = (x-xmin-ix*dx)/dx;
1049 *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;
1050 *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];
1055 setup_cmap (int grid_spacing,
1058 gmx_cmap_t * cmap_grid)
1060 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1062 int i, j, k, ii, jj, kk, idx;
1064 double dx, xmin, v, v1, v2, v12;
1067 snew(tmp_u, 2*grid_spacing);
1068 snew(tmp_u2, 2*grid_spacing);
1069 snew(tmp_yy, 2*grid_spacing);
1070 snew(tmp_y1, 2*grid_spacing);
1071 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1072 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1074 dx = 360.0/grid_spacing;
1075 xmin = -180.0-dx*grid_spacing/2;
1077 for (kk = 0; kk < nc; kk++)
1079 /* Compute an offset depending on which cmap we are using
1080 * Offset will be the map number multiplied with the
1081 * grid_spacing * grid_spacing * 2
1083 offset = kk * grid_spacing * grid_spacing * 2;
1085 for (i = 0; i < 2*grid_spacing; i++)
1087 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1089 for (j = 0; j < 2*grid_spacing; j++)
1091 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1092 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1096 for (i = 0; i < 2*grid_spacing; i++)
1098 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1101 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1103 ii = i-grid_spacing/2;
1106 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1108 jj = j-grid_spacing/2;
1111 for (k = 0; k < 2*grid_spacing; k++)
1113 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1114 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1117 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1118 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1119 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1120 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1122 idx = ii*grid_spacing+jj;
1123 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1124 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1125 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1126 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1132 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1136 cmap_grid->ngrid = ngrid;
1137 cmap_grid->grid_spacing = grid_spacing;
1138 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1140 snew(cmap_grid->cmapdata, ngrid);
1142 for (i = 0; i < cmap_grid->ngrid; i++)
1144 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1149 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1151 int count, count_mol, i, mb;
1152 gmx_molblock_t *molb;
1157 for (mb = 0; mb < mtop->nmolblock; mb++)
1160 molb = &mtop->molblock[mb];
1161 plist = mi[molb->type].plist;
1163 for (i = 0; i < F_NRE; i++)
1167 count_mol += 3*plist[i].nr;
1169 else if (interaction_function[i].flags & IF_CONSTRAINT)
1171 count_mol += plist[i].nr;
1175 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1178 "Molecule type '%s' has %d constraints.\n"
1179 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1180 *mi[molb->type].name, count_mol,
1181 nrdf_internal(&mi[molb->type].atoms));
1184 count += molb->nmol*count_mol;
1190 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1192 int i, nmiss, natoms, mt;
1194 const t_atoms *atoms;
1197 for (mt = 0; mt < sys->nmoltype; mt++)
1199 atoms = &sys->moltype[mt].atoms;
1202 for (i = 0; i < natoms; i++)
1204 q = atoms->atom[i].q;
1205 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1206 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1207 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1208 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1209 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1212 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1213 get_atomtype_name(atoms->atom[i].type, atype), q);
1221 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);
1226 static void check_gbsa_params(gpp_atomtype_t atype)
1230 /* If we are doing GBSA, check that we got the parameters we need
1231 * This checking is to see if there are GBSA paratmeters for all
1232 * atoms in the force field. To go around this for testing purposes
1233 * comment out the nerror++ counter temporarily
1236 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1238 if (get_atomtype_radius(i, atype) < 0 ||
1239 get_atomtype_vol(i, atype) < 0 ||
1240 get_atomtype_surftens(i, atype) < 0 ||
1241 get_atomtype_gb_radius(i, atype) < 0 ||
1242 get_atomtype_S_hct(i, atype) < 0)
1244 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1245 get_atomtype_name(i, atype));
1252 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);
1257 static real calc_temp(const gmx_mtop_t *mtop,
1258 const t_inputrec *ir,
1262 gmx_mtop_atomloop_all_t aloop;
1269 aloop = gmx_mtop_atomloop_all_init(mtop);
1270 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1272 sum_mv2 += atom->m*norm2(v[a]);
1276 for (g = 0; g < ir->opts.ngtc; g++)
1278 nrdf += ir->opts.nrdf[g];
1281 return sum_mv2/(nrdf*BOLTZ);
1284 static real get_max_reference_temp(const t_inputrec *ir,
1293 for (i = 0; i < ir->opts.ngtc; i++)
1295 if (ir->opts.tau_t[i] < 0)
1301 ref_t = max(ref_t, ir->opts.ref_t[i]);
1309 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.",
1317 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1324 verletbuf_list_setup_t ls;
1327 char warn_buf[STRLEN];
1329 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1331 /* Calculate the buffer size for simple atom vs atoms list */
1332 ls.cluster_size_i = 1;
1333 ls.cluster_size_j = 1;
1334 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1335 &ls, &n_nonlin_vsite, &rlist_1x1);
1337 /* Set the pair-list buffer size in ir */
1338 verletbuf_get_list_setup(FALSE, &ls);
1339 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1340 &ls, &n_nonlin_vsite, &ir->rlist);
1342 if (n_nonlin_vsite > 0)
1344 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);
1345 warning_note(wi, warn_buf);
1348 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1349 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1351 ir->rlistlong = ir->rlist;
1352 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1353 ls.cluster_size_i, ls.cluster_size_j,
1354 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1356 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1358 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)));
1362 int gmx_grompp(int argc, char *argv[])
1364 static const char *desc[] = {
1365 "[THISMODULE] (the gromacs preprocessor)",
1366 "reads a molecular topology file, checks the validity of the",
1367 "file, expands the topology from a molecular description to an atomic",
1368 "description. The topology file contains information about",
1369 "molecule types and the number of molecules, the preprocessor",
1370 "copies each molecule as needed. ",
1371 "There is no limitation on the number of molecule types. ",
1372 "Bonds and bond-angles can be converted into constraints, separately",
1373 "for hydrogens and heavy atoms.",
1374 "Then a coordinate file is read and velocities can be generated",
1375 "from a Maxwellian distribution if requested.",
1376 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1377 "(eg. number of MD steps, time step, cut-off), and others such as",
1378 "NEMD parameters, which are corrected so that the net acceleration",
1380 "Eventually a binary file is produced that can serve as the sole input",
1381 "file for the MD program.[PAR]",
1383 "[THISMODULE] uses the atom names from the topology file. The atom names",
1384 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1385 "warnings when they do not match the atom names in the topology.",
1386 "Note that the atom names are irrelevant for the simulation as",
1387 "only the atom types are used for generating interaction parameters.[PAR]",
1389 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1390 "etc. The preprocessor supports the following keywords:[PAR]",
1391 "#ifdef VARIABLE[BR]",
1392 "#ifndef VARIABLE[BR]",
1395 "#define VARIABLE[BR]",
1396 "#undef VARIABLE[BR]"
1397 "#include \"filename\"[BR]",
1398 "#include <filename>[PAR]",
1399 "The functioning of these statements in your topology may be modulated by",
1400 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1401 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1402 "include = -I/home/john/doe[tt][BR]",
1403 "For further information a C-programming textbook may help you out.",
1404 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1405 "topology file written out so that you can verify its contents.[PAR]",
1407 /* cpp has been unnecessary for some time, hasn't it?
1408 "If your system does not have a C-preprocessor, you can still",
1409 "use [TT]grompp[tt], but you do not have access to the features ",
1410 "from the cpp. Command line options to the C-preprocessor can be given",
1411 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1414 "When using position restraints a file with restraint coordinates",
1415 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1416 "with respect to the conformation from the [TT]-c[tt] option.",
1417 "For free energy calculation the the coordinates for the B topology",
1418 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1419 "those of the A topology.[PAR]",
1421 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1422 "The last frame with coordinates and velocities will be read,",
1423 "unless the [TT]-time[tt] option is used. Only if this information",
1424 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1425 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1426 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1427 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1430 "[THISMODULE] can be used to restart simulations (preserving",
1431 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1432 "However, for simply changing the number of run steps to extend",
1433 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1434 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1435 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1436 "like output frequency, then supplying the checkpoint file to",
1437 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1438 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1440 "By default, all bonded interactions which have constant energy due to",
1441 "virtual site constructions will be removed. If this constant energy is",
1442 "not zero, this will result in a shift in the total energy. All bonded",
1443 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1444 "all constraints for distances which will be constant anyway because",
1445 "of virtual site constructions will be removed. If any constraints remain",
1446 "which involve virtual sites, a fatal error will result.[PAR]"
1448 "To verify your run input file, please take note of all warnings",
1449 "on the screen, and correct where necessary. Do also look at the contents",
1450 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1451 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1452 "with the [TT]-debug[tt] option which will give you more information",
1453 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1454 "can see the contents of the run input file with the [gmx-dump]",
1455 "program. [gmx-check] can be used to compare the contents of two",
1456 "run input files.[PAR]"
1458 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1459 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1460 "harmless, but usually they are not. The user is advised to carefully",
1461 "interpret the output messages before attempting to bypass them with",
1468 gpp_atomtype_t atype;
1470 int natoms, nvsite, comb, mt;
1474 real max_spacing, fudgeQQ;
1476 char fn[STRLEN], fnB[STRLEN];
1477 const char *mdparin;
1479 gmx_bool bNeedVel, bGenVel;
1480 gmx_bool have_atomnumber;
1482 t_params *gb_plist = NULL;
1483 gmx_genborn_t *born = NULL;
1485 gmx_bool bVerbose = FALSE;
1487 char warn_buf[STRLEN];
1491 { efMDP, NULL, NULL, ffREAD },
1492 { efMDP, "-po", "mdout", ffWRITE },
1493 { efSTX, "-c", NULL, ffREAD },
1494 { efSTX, "-r", NULL, ffOPTRD },
1495 { efSTX, "-rb", NULL, ffOPTRD },
1496 { efNDX, NULL, NULL, ffOPTRD },
1497 { efTOP, NULL, NULL, ffREAD },
1498 { efTOP, "-pp", "processed", ffOPTWR },
1499 { efTPX, "-o", NULL, ffWRITE },
1500 { efTRN, "-t", NULL, ffOPTRD },
1501 { efEDR, "-e", NULL, ffOPTRD },
1502 { efTRN, "-ref", "rotref", ffOPTRW }
1504 #define NFILE asize(fnm)
1506 /* Command line options */
1507 static gmx_bool bRenum = TRUE;
1508 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1509 static int i, maxwarn = 0;
1510 static real fr_time = -1;
1512 { "-v", FALSE, etBOOL, {&bVerbose},
1513 "Be loud and noisy" },
1514 { "-time", FALSE, etREAL, {&fr_time},
1515 "Take frame at or first after this time." },
1516 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1517 "Remove constant bonded interactions with virtual sites" },
1518 { "-maxwarn", FALSE, etINT, {&maxwarn},
1519 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1520 { "-zero", FALSE, etBOOL, {&bZero},
1521 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1522 { "-renum", FALSE, etBOOL, {&bRenum},
1523 "Renumber atomtypes and minimize number of atomtypes" }
1526 /* Initiate some variables */
1531 /* Parse the command line */
1532 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1533 asize(desc), desc, 0, NULL, &oenv))
1538 wi = init_warning(TRUE, maxwarn);
1540 /* PARAMETER file processing */
1541 mdparin = opt2fn("-f", NFILE, fnm);
1542 set_warning_line(wi, mdparin, -1);
1543 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1547 fprintf(stderr, "checking input for internal consistency...\n");
1549 check_ir(mdparin, ir, opts, wi);
1551 if (ir->ld_seed == -1)
1553 ir->ld_seed = (int)gmx_rng_make_seed();
1554 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1557 if (ir->expandedvals->lmc_seed == -1)
1559 ir->expandedvals->lmc_seed = (int)gmx_rng_make_seed();
1560 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1563 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1564 bGenVel = (bNeedVel && opts->bGenVel);
1565 if (bGenVel && ir->bContinuation)
1568 "Generating velocities is inconsistent with attempting "
1569 "to continue a previous run. Choose only one of "
1570 "gen-vel = yes and continuation = yes.");
1571 warning_error(wi, warn_buf);
1577 atype = init_atomtype();
1580 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1583 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1584 if (!gmx_fexist(fn))
1586 gmx_fatal(FARGS, "%s does not exist", fn);
1588 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1589 opts, ir, bZero, bGenVel, bVerbose, &state,
1590 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1596 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1600 /* set parameters for virtual site construction (not for vsiten) */
1601 for (mt = 0; mt < sys->nmoltype; mt++)
1604 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1606 /* now throw away all obsolete bonds, angles and dihedrals: */
1607 /* note: constraints are ALWAYS removed */
1610 for (mt = 0; mt < sys->nmoltype; mt++)
1612 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1616 if (ir->cutoff_scheme == ecutsVERLET)
1618 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1619 ecutscheme_names[ir->cutoff_scheme]);
1621 /* Remove all charge groups */
1622 gmx_mtop_remove_chargegroups(sys);
1624 if (EVDW_PME(ir->vdwtype))
1626 gmx_fatal(FARGS, "LJ-PME not implemented together with verlet-scheme!");
1630 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1632 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1634 sprintf(warn_buf, "Can not do %s with %s, use %s",
1635 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1636 warning_error(wi, warn_buf);
1638 if (ir->bPeriodicMols)
1640 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1641 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1642 warning_error(wi, warn_buf);
1646 if (EI_SD (ir->eI) && ir->etc != etcNO)
1648 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1651 /* If we are doing QM/MM, check that we got the atom numbers */
1652 have_atomnumber = TRUE;
1653 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1655 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1657 if (!have_atomnumber && ir->bQMMM)
1661 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1662 "field you are using does not contain atom numbers fields. This is an\n"
1663 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1664 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1665 "an integer just before the mass column in ffXXXnb.itp.\n"
1666 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1671 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1673 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1675 /** \TODO check size of ex+hy width against box size */
1678 /* Check for errors in the input now, since they might cause problems
1679 * during processing further down.
1681 check_warning_error(wi, FARGS);
1683 if (opt2bSet("-r", NFILE, fnm))
1685 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1689 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1691 if (opt2bSet("-rb", NFILE, fnm))
1693 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1700 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1704 fprintf(stderr, "Reading position restraint coords from %s", fn);
1705 if (strcmp(fn, fnB) == 0)
1707 fprintf(stderr, "\n");
1711 fprintf(stderr, " and %s\n", fnB);
1714 gen_posres(sys, mi, fn, fnB,
1715 ir->refcoord_scaling, ir->ePBC,
1716 ir->posres_com, ir->posres_comB,
1720 /* If we are using CMAP, setup the pre-interpolation grid */
1721 if (plist->ncmap > 0)
1723 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1724 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1727 set_wall_atomtype(atype, opts, ir, wi);
1730 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1731 ntype = get_atomtype_ntypes(atype);
1734 if (ir->implicit_solvent != eisNO)
1736 /* Now we have renumbered the atom types, we can check the GBSA params */
1737 check_gbsa_params(atype);
1739 /* Check that all atoms that have charge and/or LJ-parameters also have
1740 * sensible GB-parameters
1742 check_gbsa_params_charged(sys, atype);
1745 /* PELA: Copy the atomtype data to the topology atomtype list */
1746 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1750 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1755 fprintf(stderr, "converting bonded parameters...\n");
1758 ntype = get_atomtype_ntypes(atype);
1759 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1763 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1766 /* set ptype to VSite for virtual sites */
1767 for (mt = 0; mt < sys->nmoltype; mt++)
1769 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1773 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1775 /* Check velocity for virtual sites and shells */
1778 check_vel(sys, state.v);
1784 for (i = 0; i < sys->nmoltype; i++)
1786 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1789 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1791 check_bonds_timestep(sys, ir->delta_t, wi);
1794 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1796 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.");
1799 check_warning_error(wi, FARGS);
1803 fprintf(stderr, "initialising group options...\n");
1805 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1807 bGenVel ? state.v : NULL,
1810 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1813 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
1817 if (EI_MD(ir->eI) && ir->etc == etcNO)
1821 buffer_temp = opts->tempi;
1825 buffer_temp = calc_temp(sys, ir, state.v);
1827 if (buffer_temp > 0)
1829 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
1830 warning_note(wi, warn_buf);
1834 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
1835 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
1836 warning_note(wi, warn_buf);
1841 buffer_temp = get_max_reference_temp(ir, wi);
1844 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
1846 /* NVE with initial T=0: we add a fixed ratio to rlist.
1847 * Since we don't actually use verletbuf_tol, we set it to -1
1848 * so it can't be misused later.
1850 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
1851 ir->verletbuf_tol = -1;
1855 /* We warn for NVE simulations with >1(.1)% drift tolerance */
1856 const real drift_tol = 0.01;
1859 /* We use 2 DOF per atom = 2kT pot+kin energy, to be on
1860 * the safe side with constraints (without constraints: 3 DOF).
1862 ener_runtime = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
1864 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
1866 ir->verletbuf_tol > 1.1*drift_tol*ener_runtime)
1868 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.",
1869 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
1870 (int)(ir->verletbuf_tol/ener_runtime*100 + 0.5),
1871 (int)(100*drift_tol + 0.5),
1872 drift_tol*ener_runtime);
1873 warning_note(wi, warn_buf);
1876 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
1881 /* Init the temperature coupling state */
1882 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1886 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1888 check_eg_vs_cg(sys);
1892 pr_symtab(debug, 0, "After index", &sys->symtab);
1895 triple_check(mdparin, ir, sys, wi);
1896 close_symtab(&sys->symtab);
1899 pr_symtab(debug, 0, "After close", &sys->symtab);
1902 /* make exclusions between QM atoms */
1905 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1907 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1911 generate_qmexcl(sys, ir, wi);
1915 if (ftp2bSet(efTRN, NFILE, fnm))
1919 fprintf(stderr, "getting data from old trajectory ...\n");
1921 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1922 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1925 if (ir->ePBC == epbcXY && ir->nwall != 2)
1927 clear_rvec(state.box[ZZ]);
1930 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1932 set_warning_line(wi, mdparin, -1);
1933 check_chargegroup_radii(sys, ir, state.x, wi);
1936 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1938 /* Calculate the optimal grid dimensions */
1939 copy_mat(state.box, box);
1940 if (ir->ePBC == epbcXY && ir->nwall == 2)
1942 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1944 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1946 /* Mark fourier_spacing as not used */
1947 ir->fourier_spacing = 0;
1949 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1951 set_warning_line(wi, mdparin, -1);
1952 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1954 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1955 &(ir->nkx), &(ir->nky), &(ir->nkz));
1958 /* MRS: eventually figure out better logic for initializing the fep
1959 values that makes declaring the lambda and declaring the state not
1960 potentially conflict if not handled correctly. */
1961 if (ir->efep != efepNO)
1963 if (EVDW_PME(ir->vdwtype))
1965 gmx_fatal(FARGS, "LJ-PME not implemented together with free energy calculations!");
1968 state.fep_state = ir->fepvals->init_fep_state;
1969 for (i = 0; i < efptNR; i++)
1971 /* init_lambda trumps state definitions*/
1972 if (ir->fepvals->init_lambda >= 0)
1974 state.lambda[i] = ir->fepvals->init_lambda;
1978 if (ir->fepvals->all_lambda[i] == NULL)
1980 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1984 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1990 if (ir->ePull != epullNO)
1992 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1997 set_reference_positions(ir->rot, state.x, state.box,
1998 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2002 /* reset_multinr(sys); */
2004 if (EEL_PME(ir->coulombtype))
2006 float ratio = pme_load_estimate(sys, ir, state.box);
2007 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2008 /* With free energy we might need to do PME both for the A and B state
2009 * charges. This will double the cost, but the optimal performance will
2010 * then probably be at a slightly larger cut-off and grid spacing.
2012 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2013 (ir->efep != efepNO && ratio > 2.0/3.0))
2016 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2017 "and for highly parallel simulations between 0.25 and 0.33,\n"
2018 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2019 if (ir->efep != efepNO)
2022 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2028 char warn_buf[STRLEN];
2029 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2030 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2033 set_warning_line(wi, mdparin, -1);
2034 warning_note(wi, warn_buf);
2038 printf("%s\n", warn_buf);
2044 fprintf(stderr, "writing run input file...\n");
2047 done_warning(wi, FARGS);
2048 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);
2050 done_atomtype(atype);
2051 done_mtop(sys, TRUE);
2052 done_inputrec_strings();