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.
41 #include <sys/types.h>
48 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/fileio/confio.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "grompp-impl.h"
56 #include "gromacs/gmxpreprocess/gen_maxwell_velocities.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/legacyheaders/splitter.h"
60 #include "gromacs/gmxpreprocess/sortwater.h"
62 #include "gromacs/legacyheaders/warninp.h"
63 #include "gromacs/fileio/gmxfio.h"
64 #include "gromacs/fileio/trnio.h"
65 #include "gromacs/fileio/tpxio.h"
66 #include "gromacs/fileio/trxio.h"
67 #include "vsite_parm.h"
68 #include "gromacs/legacyheaders/txtdump.h"
69 #include "gromacs/legacyheaders/calcgrid.h"
71 #include "gromacs/fileio/enxio.h"
72 #include "gromacs/legacyheaders/perf_est.h"
73 #include "compute_io.h"
74 #include "gpp_atomtype.h"
75 #include "gromacs/topology/mtop_util.h"
76 #include "gromacs/legacyheaders/genborn.h"
77 #include "calc_verletbuf.h"
79 #include "gromacs/imd/imd.h"
80 #include "gromacs/utility/cstringutil.h"
82 #include "gromacs/commandline/pargs.h"
83 #include "gromacs/pbcutil/pbc.h"
84 #include "gromacs/random/random.h"
85 #include "gromacs/topology/symtab.h"
86 #include "gromacs/topology/topology.h"
87 #include "gromacs/utility/fatalerror.h"
88 #include "gromacs/utility/smalloc.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 void check_shells_inputrec(gmx_mtop_t *mtop,
375 gmx_mtop_atomloop_all_t aloop;
378 char warn_buf[STRLEN];
380 aloop = gmx_mtop_atomloop_all_init(mtop);
381 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
383 if (atom->ptype == eptShell ||
384 atom->ptype == eptBond)
389 if (IR_TWINRANGE(*ir) && (nshells > 0))
391 snprintf(warn_buf, STRLEN,
392 "The combination of using shells and a twin-range cut-off is not supported");
393 warning_error(wi, warn_buf);
395 if ((nshells > 0) && (ir->nstcalcenergy != 1))
397 set_warning_line(wi, "unknown", -1);
398 snprintf(warn_buf, STRLEN,
399 "There are %d shells, changing nstcalcenergy from %d to 1",
400 nshells, ir->nstcalcenergy);
401 ir->nstcalcenergy = 1;
402 warning(wi, warn_buf);
406 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
411 for (mb = 0; mb < mtop->nmolblock; mb++)
413 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
419 /* This routine reorders the molecule type array
420 * in the order of use in the molblocks,
421 * unused molecule types are deleted.
423 static void renumber_moltypes(gmx_mtop_t *sys,
424 int *nmolinfo, t_molinfo **molinfo)
426 int *order, norder, i;
430 snew(order, *nmolinfo);
432 for (mb = 0; mb < sys->nmolblock; mb++)
434 for (i = 0; i < norder; i++)
436 if (order[i] == sys->molblock[mb].type)
443 /* This type did not occur yet, add it */
444 order[norder] = sys->molblock[mb].type;
445 /* Renumber the moltype in the topology */
448 sys->molblock[mb].type = i;
451 /* We still need to reorder the molinfo structs */
453 for (mi = 0; mi < *nmolinfo; mi++)
455 for (i = 0; i < norder; i++)
464 done_mi(&(*molinfo)[mi]);
468 minew[i] = (*molinfo)[mi];
477 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
482 mtop->nmoltype = nmi;
483 snew(mtop->moltype, nmi);
484 for (m = 0; m < nmi; m++)
486 molt = &mtop->moltype[m];
487 molt->name = mi[m].name;
488 molt->atoms = mi[m].atoms;
489 /* ilists are copied later */
490 molt->cgs = mi[m].cgs;
491 molt->excls = mi[m].excls;
496 new_status(const char *topfile, const char *topppfile, const char *confin,
497 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
498 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
499 gpp_atomtype_t atype, gmx_mtop_t *sys,
500 int *nmi, t_molinfo **mi, t_params plist[],
501 int *comb, double *reppow, real *fudgeQQ,
505 t_molinfo *molinfo = NULL;
507 gmx_molblock_t *molblock, *molbs;
509 int mb, i, nrmols, nmismatch;
511 gmx_bool bGB = FALSE;
512 char warn_buf[STRLEN];
516 /* Set gmx_boolean for GB */
517 if (ir->implicit_solvent)
522 /* TOPOLOGY processing */
523 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
524 plist, comb, reppow, fudgeQQ,
525 atype, &nrmols, &molinfo, ir,
526 &nmolblock, &molblock, bGB,
530 snew(sys->molblock, nmolblock);
533 for (mb = 0; mb < nmolblock; mb++)
535 if (sys->nmolblock > 0 &&
536 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
538 /* Merge consecutive blocks with the same molecule type */
539 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
540 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
542 else if (molblock[mb].nmol > 0)
544 /* Add a new molblock to the topology */
545 molbs = &sys->molblock[sys->nmolblock];
546 *molbs = molblock[mb];
547 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
548 molbs->nposres_xA = 0;
549 molbs->nposres_xB = 0;
550 sys->natoms += molbs->nmol*molbs->natoms_mol;
554 if (sys->nmolblock == 0)
556 gmx_fatal(FARGS, "No molecules were defined in the system");
559 renumber_moltypes(sys, &nrmols, &molinfo);
563 convert_harmonics(nrmols, molinfo, atype);
566 if (ir->eDisre == edrNone)
568 i = rm_interactions(F_DISRES, nrmols, molinfo);
571 set_warning_line(wi, "unknown", -1);
572 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
573 warning_note(wi, warn_buf);
576 if (opts->bOrire == FALSE)
578 i = rm_interactions(F_ORIRES, nrmols, molinfo);
581 set_warning_line(wi, "unknown", -1);
582 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
583 warning_note(wi, warn_buf);
587 /* Copy structures from msys to sys */
588 molinfo2mtop(nrmols, molinfo, sys);
590 gmx_mtop_finalize(sys);
592 /* COORDINATE file processing */
595 fprintf(stderr, "processing coordinates...\n");
598 get_stx_coordnum(confin, &state->natoms);
599 if (state->natoms != sys->natoms)
601 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
602 " does not match topology (%s, %d)",
603 confin, state->natoms, topfile, sys->natoms);
607 /* make space for coordinates and velocities */
610 init_t_atoms(confat, state->natoms, FALSE);
611 init_state(state, state->natoms, 0, 0, 0, 0);
612 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
613 /* This call fixes the box shape for runs with pressure scaling */
614 set_box_rel(ir, state);
616 nmismatch = check_atom_names(topfile, confin, sys, confat);
617 free_t_atoms(confat, TRUE);
622 sprintf(buf, "%d non-matching atom name%s\n"
623 "atom names from %s will be used\n"
624 "atom names from %s will be ignored\n",
625 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
630 fprintf(stderr, "double-checking input for internal consistency...\n");
632 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
638 gmx_mtop_atomloop_all_t aloop;
642 snew(mass, state->natoms);
643 aloop = gmx_mtop_atomloop_all_init(sys);
644 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
650 if (opts->seed == -1)
652 useed = (int)gmx_rng_make_seed();
653 fprintf(stderr, "Setting gen_seed to %u\n", useed);
655 maxwell_speed(opts->tempi, useed, sys, state->v);
657 stop_cm(stdout, state->natoms, mass, state->x, state->v);
665 static void copy_state(const char *slog, t_trxframe *fr,
666 gmx_bool bReadVel, t_state *state,
671 if (fr->not_ok & FRAME_NOT_OK)
673 gmx_fatal(FARGS, "Can not start from an incomplete frame");
677 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
681 for (i = 0; i < state->natoms; i++)
683 copy_rvec(fr->x[i], state->x[i]);
689 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
691 for (i = 0; i < state->natoms; i++)
693 copy_rvec(fr->v[i], state->v[i]);
698 copy_mat(fr->box, state->box);
701 *use_time = fr->time;
704 static void cont_status(const char *slog, const char *ener,
705 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
706 t_inputrec *ir, t_state *state,
708 const output_env_t oenv)
709 /* If fr_time == -1 read the last frame available which is complete */
717 bReadVel = (bNeedVel && !bGenVel);
720 "Reading Coordinates%s and Box size from old trajectory\n",
721 bReadVel ? ", Velocities" : "");
724 fprintf(stderr, "Will read whole trajectory\n");
728 fprintf(stderr, "Will read till time %g\n", fr_time);
734 fprintf(stderr, "Velocities generated: "
735 "ignoring velocities in input trajectory\n");
737 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
741 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
747 "WARNING: Did not find a frame with velocities in file %s,\n"
748 " all velocities will be set to zero!\n\n", slog);
749 for (i = 0; i < sys->natoms; i++)
751 clear_rvec(state->v[i]);
754 /* Search for a frame without velocities */
756 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
760 state->natoms = fr.natoms;
762 if (sys->natoms != state->natoms)
764 gmx_fatal(FARGS, "Number of atoms in Topology "
765 "is not the same as in Trajectory");
767 copy_state(slog, &fr, bReadVel, state, &use_time);
769 /* Find the appropriate frame */
770 while ((fr_time == -1 || fr.time < fr_time) &&
771 read_next_frame(oenv, fp, &fr))
773 copy_state(slog, &fr, bReadVel, state, &use_time);
778 /* Set the relative box lengths for preserving the box shape.
779 * Note that this call can lead to differences in the last bit
780 * with respect to using gmx convert-tpr to create a [TT].tpx[tt] file.
782 set_box_rel(ir, state);
784 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
785 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
787 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
789 get_enx_state(ener, use_time, &sys->groups, ir, state);
790 preserve_box_shape(ir, state->box_rel, state->boxv);
794 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
796 int rc_scaling, int ePBC,
800 gmx_bool bFirst = TRUE, *hadAtom;
806 int natoms, npbcdim = 0;
807 char warn_buf[STRLEN], title[STRLEN];
808 int a, i, ai, j, k, mb, nat_molb;
809 gmx_molblock_t *molb;
813 get_stx_coordnum(fn, &natoms);
814 if (natoms != mtop->natoms)
816 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);
817 warning(wi, warn_buf);
821 init_t_atoms(&dumat, natoms, FALSE);
822 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
824 npbcdim = ePBC2npbcdim(ePBC);
826 if (rc_scaling != erscNO)
828 copy_mat(box, invbox);
829 for (j = npbcdim; j < DIM; j++)
831 clear_rvec(invbox[j]);
834 m_inv_ur0(invbox, invbox);
837 /* Copy the reference coordinates to mtop */
841 snew(hadAtom, natoms);
842 for (mb = 0; mb < mtop->nmolblock; mb++)
844 molb = &mtop->molblock[mb];
845 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
846 pr = &(molinfo[molb->type].plist[F_POSRES]);
847 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
848 if (pr->nr > 0 || prfb->nr > 0)
850 atom = mtop->moltype[molb->type].atoms.atom;
851 for (i = 0; (i < pr->nr); i++)
853 ai = pr->param[i].AI;
856 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
857 ai+1, *molinfo[molb->type].name, fn, natoms);
860 if (rc_scaling == erscCOM)
862 /* Determine the center of mass of the posres reference coordinates */
863 for (j = 0; j < npbcdim; j++)
865 sum[j] += atom[ai].m*x[a+ai][j];
867 totmass += atom[ai].m;
870 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
871 for (i = 0; (i < prfb->nr); i++)
873 ai = prfb->param[i].AI;
876 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
877 ai+1, *molinfo[molb->type].name, fn, natoms);
879 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
881 /* Determine the center of mass of the posres reference coordinates */
882 for (j = 0; j < npbcdim; j++)
884 sum[j] += atom[ai].m*x[a+ai][j];
886 totmass += atom[ai].m;
891 molb->nposres_xA = nat_molb;
892 snew(molb->posres_xA, molb->nposres_xA);
893 for (i = 0; i < nat_molb; i++)
895 copy_rvec(x[a+i], molb->posres_xA[i]);
900 molb->nposres_xB = nat_molb;
901 snew(molb->posres_xB, molb->nposres_xB);
902 for (i = 0; i < nat_molb; i++)
904 copy_rvec(x[a+i], molb->posres_xB[i]);
910 if (rc_scaling == erscCOM)
914 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
916 for (j = 0; j < npbcdim; j++)
918 com[j] = sum[j]/totmass;
920 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]);
923 if (rc_scaling != erscNO)
925 assert(npbcdim <= DIM);
927 for (mb = 0; mb < mtop->nmolblock; mb++)
929 molb = &mtop->molblock[mb];
930 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
931 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
933 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
934 for (i = 0; i < nat_molb; i++)
936 for (j = 0; j < npbcdim; j++)
938 if (rc_scaling == erscALL)
940 /* Convert from Cartesian to crystal coordinates */
941 xp[i][j] *= invbox[j][j];
942 for (k = j+1; k < npbcdim; k++)
944 xp[i][j] += invbox[k][j]*xp[i][k];
947 else if (rc_scaling == erscCOM)
949 /* Subtract the center of mass */
957 if (rc_scaling == erscCOM)
959 /* Convert the COM from Cartesian to crystal coordinates */
960 for (j = 0; j < npbcdim; j++)
962 com[j] *= invbox[j][j];
963 for (k = j+1; k < npbcdim; k++)
965 com[j] += invbox[k][j]*com[k];
971 free_t_atoms(&dumat, TRUE);
977 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
978 char *fnA, char *fnB,
979 int rc_scaling, int ePBC,
985 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
986 /* It is safer to simply read the b-state posres rather than trying
987 * to be smart and copy the positions.
989 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
992 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
993 t_inputrec *ir, warninp_t wi)
996 char warn_buf[STRLEN];
1000 fprintf(stderr, "Searching the wall atom type(s)\n");
1002 for (i = 0; i < ir->nwall; i++)
1004 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
1005 if (ir->wall_atomtype[i] == NOTSET)
1007 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
1008 warning_error(wi, warn_buf);
1013 static int nrdf_internal(t_atoms *atoms)
1018 for (i = 0; i < atoms->nr; i++)
1020 /* Vsite ptype might not be set here yet, so also check the mass */
1021 if ((atoms->atom[i].ptype == eptAtom ||
1022 atoms->atom[i].ptype == eptNucleus)
1023 && atoms->atom[i].m > 0)
1030 case 0: nrdf = 0; break;
1031 case 1: nrdf = 0; break;
1032 case 2: nrdf = 1; break;
1033 default: nrdf = nmass*3 - 6; break;
1040 spline1d( double dx,
1052 for (i = 1; i < n-1; i++)
1054 p = 0.5*y2[i-1]+2.0;
1056 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1057 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1062 for (i = n-2; i >= 0; i--)
1064 y2[i] = y2[i]*y2[i+1]+u[i];
1070 interpolate1d( double xmin,
1083 a = (xmin+(ix+1)*dx-x)/dx;
1084 b = (x-xmin-ix*dx)/dx;
1086 *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;
1087 *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];
1092 setup_cmap (int grid_spacing,
1095 gmx_cmap_t * cmap_grid)
1097 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1099 int i, j, k, ii, jj, kk, idx;
1101 double dx, xmin, v, v1, v2, v12;
1104 snew(tmp_u, 2*grid_spacing);
1105 snew(tmp_u2, 2*grid_spacing);
1106 snew(tmp_yy, 2*grid_spacing);
1107 snew(tmp_y1, 2*grid_spacing);
1108 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1109 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1111 dx = 360.0/grid_spacing;
1112 xmin = -180.0-dx*grid_spacing/2;
1114 for (kk = 0; kk < nc; kk++)
1116 /* Compute an offset depending on which cmap we are using
1117 * Offset will be the map number multiplied with the
1118 * grid_spacing * grid_spacing * 2
1120 offset = kk * grid_spacing * grid_spacing * 2;
1122 for (i = 0; i < 2*grid_spacing; i++)
1124 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1126 for (j = 0; j < 2*grid_spacing; j++)
1128 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1129 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1133 for (i = 0; i < 2*grid_spacing; i++)
1135 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1138 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1140 ii = i-grid_spacing/2;
1143 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1145 jj = j-grid_spacing/2;
1148 for (k = 0; k < 2*grid_spacing; k++)
1150 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1151 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1154 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1155 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1156 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1157 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1159 idx = ii*grid_spacing+jj;
1160 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1161 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1162 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1163 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1169 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1173 cmap_grid->ngrid = ngrid;
1174 cmap_grid->grid_spacing = grid_spacing;
1175 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1177 snew(cmap_grid->cmapdata, ngrid);
1179 for (i = 0; i < cmap_grid->ngrid; i++)
1181 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1186 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1188 int count, count_mol, i, mb;
1189 gmx_molblock_t *molb;
1194 for (mb = 0; mb < mtop->nmolblock; mb++)
1197 molb = &mtop->molblock[mb];
1198 plist = mi[molb->type].plist;
1200 for (i = 0; i < F_NRE; i++)
1204 count_mol += 3*plist[i].nr;
1206 else if (interaction_function[i].flags & IF_CONSTRAINT)
1208 count_mol += plist[i].nr;
1212 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1215 "Molecule type '%s' has %d constraints.\n"
1216 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1217 *mi[molb->type].name, count_mol,
1218 nrdf_internal(&mi[molb->type].atoms));
1221 count += molb->nmol*count_mol;
1227 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1229 int i, nmiss, natoms, mt;
1231 const t_atoms *atoms;
1234 for (mt = 0; mt < sys->nmoltype; mt++)
1236 atoms = &sys->moltype[mt].atoms;
1239 for (i = 0; i < natoms; i++)
1241 q = atoms->atom[i].q;
1242 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1243 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1244 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1245 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1246 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1249 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1250 get_atomtype_name(atoms->atom[i].type, atype), q);
1258 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);
1263 static void check_gbsa_params(gpp_atomtype_t atype)
1267 /* If we are doing GBSA, check that we got the parameters we need
1268 * This checking is to see if there are GBSA paratmeters for all
1269 * atoms in the force field. To go around this for testing purposes
1270 * comment out the nerror++ counter temporarily
1273 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1275 if (get_atomtype_radius(i, atype) < 0 ||
1276 get_atomtype_vol(i, atype) < 0 ||
1277 get_atomtype_surftens(i, atype) < 0 ||
1278 get_atomtype_gb_radius(i, atype) < 0 ||
1279 get_atomtype_S_hct(i, atype) < 0)
1281 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1282 get_atomtype_name(i, atype));
1289 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);
1294 static real calc_temp(const gmx_mtop_t *mtop,
1295 const t_inputrec *ir,
1299 gmx_mtop_atomloop_all_t aloop;
1306 aloop = gmx_mtop_atomloop_all_init(mtop);
1307 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
1309 sum_mv2 += atom->m*norm2(v[a]);
1313 for (g = 0; g < ir->opts.ngtc; g++)
1315 nrdf += ir->opts.nrdf[g];
1318 return sum_mv2/(nrdf*BOLTZ);
1321 static real get_max_reference_temp(const t_inputrec *ir,
1330 for (i = 0; i < ir->opts.ngtc; i++)
1332 if (ir->opts.tau_t[i] < 0)
1338 ref_t = max(ref_t, ir->opts.ref_t[i]);
1346 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.",
1354 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1361 verletbuf_list_setup_t ls;
1364 char warn_buf[STRLEN];
1366 printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol, buffer_temp);
1368 /* Calculate the buffer size for simple atom vs atoms list */
1369 ls.cluster_size_i = 1;
1370 ls.cluster_size_j = 1;
1371 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1372 &ls, &n_nonlin_vsite, &rlist_1x1);
1374 /* Set the pair-list buffer size in ir */
1375 verletbuf_get_list_setup(FALSE, &ls);
1376 calc_verlet_buffer_size(mtop, det(box), ir, buffer_temp,
1377 &ls, &n_nonlin_vsite, &ir->rlist);
1379 if (n_nonlin_vsite > 0)
1381 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);
1382 warning_note(wi, warn_buf);
1385 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1386 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1388 ir->rlistlong = ir->rlist;
1389 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1390 ls.cluster_size_i, ls.cluster_size_j,
1391 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1393 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1395 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)));
1399 int gmx_grompp(int argc, char *argv[])
1401 static const char *desc[] = {
1402 "[THISMODULE] (the gromacs preprocessor)",
1403 "reads a molecular topology file, checks the validity of the",
1404 "file, expands the topology from a molecular description to an atomic",
1405 "description. The topology file contains information about",
1406 "molecule types and the number of molecules, the preprocessor",
1407 "copies each molecule as needed. ",
1408 "There is no limitation on the number of molecule types. ",
1409 "Bonds and bond-angles can be converted into constraints, separately",
1410 "for hydrogens and heavy atoms.",
1411 "Then a coordinate file is read and velocities can be generated",
1412 "from a Maxwellian distribution if requested.",
1413 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1414 "(eg. number of MD steps, time step, cut-off), and others such as",
1415 "NEMD parameters, which are corrected so that the net acceleration",
1417 "Eventually a binary file is produced that can serve as the sole input",
1418 "file for the MD program.[PAR]",
1420 "[THISMODULE] uses the atom names from the topology file. The atom names",
1421 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1422 "warnings when they do not match the atom names in the topology.",
1423 "Note that the atom names are irrelevant for the simulation as",
1424 "only the atom types are used for generating interaction parameters.[PAR]",
1426 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1427 "etc. The preprocessor supports the following keywords:[PAR]",
1428 "#ifdef VARIABLE[BR]",
1429 "#ifndef VARIABLE[BR]",
1432 "#define VARIABLE[BR]",
1433 "#undef VARIABLE[BR]"
1434 "#include \"filename\"[BR]",
1435 "#include <filename>[PAR]",
1436 "The functioning of these statements in your topology may be modulated by",
1437 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1438 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1439 "include = -I/home/john/doe[tt][BR]",
1440 "For further information a C-programming textbook may help you out.",
1441 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1442 "topology file written out so that you can verify its contents.[PAR]",
1444 /* cpp has been unnecessary for some time, hasn't it?
1445 "If your system does not have a C-preprocessor, you can still",
1446 "use [TT]grompp[tt], but you do not have access to the features ",
1447 "from the cpp. Command line options to the C-preprocessor can be given",
1448 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1451 "When using position restraints a file with restraint coordinates",
1452 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1453 "with respect to the conformation from the [TT]-c[tt] option.",
1454 "For free energy calculation the the coordinates for the B topology",
1455 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1456 "those of the A topology.[PAR]",
1458 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1459 "The last frame with coordinates and velocities will be read,",
1460 "unless the [TT]-time[tt] option is used. Only if this information",
1461 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1462 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1463 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1464 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1467 "[THISMODULE] can be used to restart simulations (preserving",
1468 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1469 "However, for simply changing the number of run steps to extend",
1470 "a run, using [gmx-convert-tpr] is more convenient than [THISMODULE].",
1471 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1472 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1473 "like output frequency, then supplying the checkpoint file to",
1474 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1475 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1477 "By default, all bonded interactions which have constant energy due to",
1478 "virtual site constructions will be removed. If this constant energy is",
1479 "not zero, this will result in a shift in the total energy. All bonded",
1480 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1481 "all constraints for distances which will be constant anyway because",
1482 "of virtual site constructions will be removed. If any constraints remain",
1483 "which involve virtual sites, a fatal error will result.[PAR]"
1485 "To verify your run input file, please take note of all warnings",
1486 "on the screen, and correct where necessary. Do also look at the contents",
1487 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1488 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1489 "with the [TT]-debug[tt] option which will give you more information",
1490 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1491 "can see the contents of the run input file with the [gmx-dump]",
1492 "program. [gmx-check] can be used to compare the contents of two",
1493 "run input files.[PAR]"
1495 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1496 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1497 "harmless, but usually they are not. The user is advised to carefully",
1498 "interpret the output messages before attempting to bypass them with",
1505 gpp_atomtype_t atype;
1507 int natoms, nvsite, comb, mt;
1511 real max_spacing, fudgeQQ;
1513 char fn[STRLEN], fnB[STRLEN];
1514 const char *mdparin;
1516 gmx_bool bNeedVel, bGenVel;
1517 gmx_bool have_atomnumber;
1519 t_params *gb_plist = NULL;
1520 gmx_genborn_t *born = NULL;
1522 gmx_bool bVerbose = FALSE;
1524 char warn_buf[STRLEN];
1526 t_atoms IMDatoms; /* Atoms to be operated on interactively (IMD) */
1529 { efMDP, NULL, NULL, ffREAD },
1530 { efMDP, "-po", "mdout", ffWRITE },
1531 { efSTX, "-c", NULL, ffREAD },
1532 { efSTX, "-r", NULL, ffOPTRD },
1533 { efSTX, "-rb", NULL, ffOPTRD },
1534 { efNDX, NULL, NULL, ffOPTRD },
1535 { efTOP, NULL, NULL, ffREAD },
1536 { efTOP, "-pp", "processed", ffOPTWR },
1537 { efTPR, "-o", NULL, ffWRITE },
1538 { efTRN, "-t", NULL, ffOPTRD },
1539 { efEDR, "-e", NULL, ffOPTRD },
1540 /* This group is needed by the VMD viewer as the start configuration for IMD sessions: */
1541 { efGRO, "-imd", "imdgroup", ffOPTWR },
1542 { efTRN, "-ref", "rotref", ffOPTRW }
1544 #define NFILE asize(fnm)
1546 /* Command line options */
1547 static gmx_bool bRenum = TRUE;
1548 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1549 static int i, maxwarn = 0;
1550 static real fr_time = -1;
1552 { "-v", FALSE, etBOOL, {&bVerbose},
1553 "Be loud and noisy" },
1554 { "-time", FALSE, etREAL, {&fr_time},
1555 "Take frame at or first after this time." },
1556 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1557 "Remove constant bonded interactions with virtual sites" },
1558 { "-maxwarn", FALSE, etINT, {&maxwarn},
1559 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1560 { "-zero", FALSE, etBOOL, {&bZero},
1561 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1562 { "-renum", FALSE, etBOOL, {&bRenum},
1563 "Renumber atomtypes and minimize number of atomtypes" }
1566 /* Initiate some variables */
1571 /* Parse the command line */
1572 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1573 asize(desc), desc, 0, NULL, &oenv))
1578 wi = init_warning(TRUE, maxwarn);
1580 /* PARAMETER file processing */
1581 mdparin = opt2fn("-f", NFILE, fnm);
1582 set_warning_line(wi, mdparin, -1);
1583 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1587 fprintf(stderr, "checking input for internal consistency...\n");
1589 check_ir(mdparin, ir, opts, wi);
1591 if (ir->ld_seed == -1)
1593 ir->ld_seed = (gmx_int64_t)gmx_rng_make_seed();
1594 fprintf(stderr, "Setting the LD random seed to %"GMX_PRId64 "\n", ir->ld_seed);
1597 if (ir->expandedvals->lmc_seed == -1)
1599 ir->expandedvals->lmc_seed = (int)gmx_rng_make_seed();
1600 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1603 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1604 bGenVel = (bNeedVel && opts->bGenVel);
1605 if (bGenVel && ir->bContinuation)
1608 "Generating velocities is inconsistent with attempting "
1609 "to continue a previous run. Choose only one of "
1610 "gen-vel = yes and continuation = yes.");
1611 warning_error(wi, warn_buf);
1617 atype = init_atomtype();
1620 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1623 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1624 if (!gmx_fexist(fn))
1626 gmx_fatal(FARGS, "%s does not exist", fn);
1628 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1629 opts, ir, bZero, bGenVel, bVerbose, &state,
1630 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1636 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1640 /* set parameters for virtual site construction (not for vsiten) */
1641 for (mt = 0; mt < sys->nmoltype; mt++)
1644 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1646 /* now throw away all obsolete bonds, angles and dihedrals: */
1647 /* note: constraints are ALWAYS removed */
1650 for (mt = 0; mt < sys->nmoltype; mt++)
1652 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1656 if (nvsite && ir->eI == eiNM)
1658 gmx_fatal(FARGS, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
1661 if (ir->cutoff_scheme == ecutsVERLET)
1663 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1664 ecutscheme_names[ir->cutoff_scheme]);
1666 /* Remove all charge groups */
1667 gmx_mtop_remove_chargegroups(sys);
1670 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1672 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1674 sprintf(warn_buf, "Can not do %s with %s, use %s",
1675 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1676 warning_error(wi, warn_buf);
1678 if (ir->bPeriodicMols)
1680 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1681 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1682 warning_error(wi, warn_buf);
1686 if (EI_SD (ir->eI) && ir->etc != etcNO)
1688 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1691 /* If we are doing QM/MM, check that we got the atom numbers */
1692 have_atomnumber = TRUE;
1693 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1695 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1697 if (!have_atomnumber && ir->bQMMM)
1701 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1702 "field you are using does not contain atom numbers fields. This is an\n"
1703 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1704 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1705 "an integer just before the mass column in ffXXXnb.itp.\n"
1706 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1711 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1713 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1715 /** TODO check size of ex+hy width against box size */
1718 /* Check for errors in the input now, since they might cause problems
1719 * during processing further down.
1721 check_warning_error(wi, FARGS);
1723 if (opt2bSet("-r", NFILE, fnm))
1725 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1729 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1731 if (opt2bSet("-rb", NFILE, fnm))
1733 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1740 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1744 fprintf(stderr, "Reading position restraint coords from %s", fn);
1745 if (strcmp(fn, fnB) == 0)
1747 fprintf(stderr, "\n");
1751 fprintf(stderr, " and %s\n", fnB);
1754 gen_posres(sys, mi, fn, fnB,
1755 ir->refcoord_scaling, ir->ePBC,
1756 ir->posres_com, ir->posres_comB,
1760 /* If we are using CMAP, setup the pre-interpolation grid */
1761 if (plist[F_CMAP].ncmap > 0)
1763 init_cmap_grid(&sys->ffparams.cmap_grid, plist[F_CMAP].nc, plist[F_CMAP].grid_spacing);
1764 setup_cmap(plist[F_CMAP].grid_spacing, plist[F_CMAP].nc, plist[F_CMAP].cmap, &sys->ffparams.cmap_grid);
1767 set_wall_atomtype(atype, opts, ir, wi);
1770 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1771 ntype = get_atomtype_ntypes(atype);
1774 if (ir->implicit_solvent != eisNO)
1776 /* Now we have renumbered the atom types, we can check the GBSA params */
1777 check_gbsa_params(atype);
1779 /* Check that all atoms that have charge and/or LJ-parameters also have
1780 * sensible GB-parameters
1782 check_gbsa_params_charged(sys, atype);
1785 /* PELA: Copy the atomtype data to the topology atomtype list */
1786 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1790 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1795 fprintf(stderr, "converting bonded parameters...\n");
1798 ntype = get_atomtype_ntypes(atype);
1799 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1803 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1806 /* set ptype to VSite for virtual sites */
1807 for (mt = 0; mt < sys->nmoltype; mt++)
1809 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1813 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1815 /* Check velocity for virtual sites and shells */
1818 check_vel(sys, state.v);
1821 /* check for shells and inpurecs */
1822 check_shells_inputrec(sys, ir, wi);
1827 for (i = 0; i < sys->nmoltype; i++)
1829 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1832 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1834 check_bonds_timestep(sys, ir->delta_t, wi);
1837 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1839 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.");
1842 check_warning_error(wi, FARGS);
1846 fprintf(stderr, "initialising group options...\n");
1848 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1850 bGenVel ? state.v : NULL,
1853 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_tol > 0 &&
1856 if (EI_DYNAMICS(ir->eI) && inputrec2nboundeddim(ir) == 3)
1860 if (EI_MD(ir->eI) && ir->etc == etcNO)
1864 buffer_temp = opts->tempi;
1868 buffer_temp = calc_temp(sys, ir, state.v);
1870 if (buffer_temp > 0)
1872 sprintf(warn_buf, "NVE simulation: will use the initial temperature of %.3f K for determining the Verlet buffer size", buffer_temp);
1873 warning_note(wi, warn_buf);
1877 sprintf(warn_buf, "NVE simulation with an initial temperature of zero: will use a Verlet buffer of %d%%. Check your energy drift!",
1878 (int)(verlet_buffer_ratio_NVE_T0*100 + 0.5));
1879 warning_note(wi, warn_buf);
1884 buffer_temp = get_max_reference_temp(ir, wi);
1887 if (EI_MD(ir->eI) && ir->etc == etcNO && buffer_temp == 0)
1889 /* NVE with initial T=0: we add a fixed ratio to rlist.
1890 * Since we don't actually use verletbuf_tol, we set it to -1
1891 * so it can't be misused later.
1893 ir->rlist *= 1.0 + verlet_buffer_ratio_NVE_T0;
1894 ir->verletbuf_tol = -1;
1898 /* We warn for NVE simulations with >1(.1)% drift tolerance */
1899 const real drift_tol = 0.01;
1902 /* We use 2 DOF per atom = 2kT pot+kin energy, to be on
1903 * the safe side with constraints (without constraints: 3 DOF).
1905 ener_runtime = 2*BOLTZ*buffer_temp/(ir->nsteps*ir->delta_t);
1907 if (EI_MD(ir->eI) && ir->etc == etcNO && ir->nstlist > 1 &&
1909 ir->verletbuf_tol > 1.1*drift_tol*ener_runtime)
1911 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.",
1912 ir->verletbuf_tol, ir->nsteps*ir->delta_t,
1913 (int)(ir->verletbuf_tol/ener_runtime*100 + 0.5),
1914 (int)(100*drift_tol + 0.5),
1915 drift_tol*ener_runtime);
1916 warning_note(wi, warn_buf);
1919 set_verlet_buffer(sys, ir, buffer_temp, state.box, wi);
1924 /* Init the temperature coupling state */
1925 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1929 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1931 check_eg_vs_cg(sys);
1935 pr_symtab(debug, 0, "After index", &sys->symtab);
1938 triple_check(mdparin, ir, sys, wi);
1939 close_symtab(&sys->symtab);
1942 pr_symtab(debug, 0, "After close", &sys->symtab);
1945 /* make exclusions between QM atoms */
1948 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1950 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1954 generate_qmexcl(sys, ir, wi);
1958 if (ftp2bSet(efTRN, NFILE, fnm))
1962 fprintf(stderr, "getting data from old trajectory ...\n");
1964 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1965 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1968 if (ir->ePBC == epbcXY && ir->nwall != 2)
1970 clear_rvec(state.box[ZZ]);
1973 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1975 set_warning_line(wi, mdparin, -1);
1976 check_chargegroup_radii(sys, ir, state.x, wi);
1979 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1981 /* Calculate the optimal grid dimensions */
1982 copy_mat(state.box, box);
1983 if (ir->ePBC == epbcXY && ir->nwall == 2)
1985 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1987 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1989 /* Mark fourier_spacing as not used */
1990 ir->fourier_spacing = 0;
1992 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1994 set_warning_line(wi, mdparin, -1);
1995 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1997 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1998 &(ir->nkx), &(ir->nky), &(ir->nkz));
2001 /* MRS: eventually figure out better logic for initializing the fep
2002 values that makes declaring the lambda and declaring the state not
2003 potentially conflict if not handled correctly. */
2004 if (ir->efep != efepNO)
2006 state.fep_state = ir->fepvals->init_fep_state;
2007 for (i = 0; i < efptNR; i++)
2009 /* init_lambda trumps state definitions*/
2010 if (ir->fepvals->init_lambda >= 0)
2012 state.lambda[i] = ir->fepvals->init_lambda;
2016 if (ir->fepvals->all_lambda[i] == NULL)
2018 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
2022 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
2028 if (ir->ePull != epullNO)
2030 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
2035 set_reference_positions(ir->rot, state.x, state.box,
2036 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
2040 /* reset_multinr(sys); */
2042 if (EEL_PME(ir->coulombtype))
2044 float ratio = pme_load_estimate(sys, ir, state.box);
2045 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
2046 /* With free energy we might need to do PME both for the A and B state
2047 * charges. This will double the cost, but the optimal performance will
2048 * then probably be at a slightly larger cut-off and grid spacing.
2050 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
2051 (ir->efep != efepNO && ratio > 2.0/3.0))
2054 "The optimal PME mesh load for parallel simulations is below 0.5\n"
2055 "and for highly parallel simulations between 0.25 and 0.33,\n"
2056 "for higher performance, increase the cut-off and the PME grid spacing.\n");
2057 if (ir->efep != efepNO)
2060 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
2066 char warn_buf[STRLEN];
2067 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
2068 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
2071 set_warning_line(wi, mdparin, -1);
2072 warning_note(wi, warn_buf);
2076 printf("%s\n", warn_buf);
2082 fprintf(stderr, "writing run input file...\n");
2085 done_warning(wi, FARGS);
2086 write_tpx_state(ftp2fn(efTPR, NFILE, fnm), ir, &state, sys);
2088 /* Output IMD group, if bIMD is TRUE */
2089 write_IMDgroup_to_file(ir->bIMD, ir, &state, sys, NFILE, fnm);
2091 done_atomtype(atype);
2092 done_mtop(sys, TRUE);
2093 done_inputrec_strings();