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, 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>
54 #include "gromacs/fileio/confio.h"
61 #include "gromacs/fileio/futil.h"
64 #include "sortwater.h"
66 #include "gmx_fatal.h"
69 #include "gromacs/fileio/gmxfio.h"
70 #include "gromacs/fileio/trnio.h"
71 #include "gromacs/fileio/tpxio.h"
72 #include "gromacs/fileio/trxio.h"
73 #include "vsite_parm.h"
77 #include "gromacs/fileio/enxio.h"
79 #include "compute_io.h"
80 #include "gpp_atomtype.h"
81 #include "mtop_util.h"
83 #include "calc_verletbuf.h"
87 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
92 /* For all the molecule types */
93 for (i = 0; i < nrmols; i++)
95 n += mols[i].plist[ifunc].nr;
96 mols[i].plist[ifunc].nr = 0;
101 static int check_atom_names(const char *fn1, const char *fn2,
102 gmx_mtop_t *mtop, t_atoms *at)
104 int mb, m, i, j, nmismatch;
106 #define MAXMISMATCH 20
108 if (mtop->natoms != at->nr)
110 gmx_incons("comparing atom names");
115 for (mb = 0; mb < mtop->nmolblock; mb++)
117 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
118 for (m = 0; m < mtop->molblock[mb].nmol; m++)
120 for (j = 0; j < tat->nr; j++)
122 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
124 if (nmismatch < MAXMISMATCH)
127 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
128 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
130 else if (nmismatch == MAXMISMATCH)
132 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
144 static void check_eg_vs_cg(gmx_mtop_t *mtop)
146 int astart, mb, m, cg, j, firstj;
147 unsigned char firsteg, eg;
150 /* Go through all the charge groups and make sure all their
151 * atoms are in the same energy group.
155 for (mb = 0; mb < mtop->nmolblock; mb++)
157 molt = &mtop->moltype[mtop->molblock[mb].type];
158 for (m = 0; m < mtop->molblock[mb].nmol; m++)
160 for (cg = 0; cg < molt->cgs.nr; cg++)
162 /* Get the energy group of the first atom in this charge group */
163 firstj = astart + molt->cgs.index[cg];
164 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
165 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
167 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
170 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
171 firstj+1, astart+j+1, cg+1, *molt->name);
175 astart += molt->atoms.nr;
180 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
183 char warn_buf[STRLEN];
186 for (cg = 0; cg < cgs->nr; cg++)
188 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
191 if (maxsize > MAX_CHARGEGROUP_SIZE)
193 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
195 else if (maxsize > 10)
197 set_warning_line(wi, topfn, -1);
199 "The largest charge group contains %d atoms.\n"
200 "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"
201 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
202 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
204 warning_note(wi, warn_buf);
208 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
210 /* This check is not intended to ensure accurate integration,
211 * rather it is to signal mistakes in the mdp settings.
212 * A common mistake is to forget to turn on constraints
213 * for MD after energy minimization with flexible bonds.
214 * This check can also detect too large time steps for flexible water
215 * models, but such errors will often be masked by the constraints
216 * mdp options, which turns flexible water into water with bond constraints,
217 * but without an angle constraint. Unfortunately such incorrect use
218 * of water models can not easily be detected without checking
219 * for specific model names.
221 * The stability limit of leap-frog or velocity verlet is 4.44 steps
222 * per oscillational period.
223 * But accurate bonds distributions are lost far before that limit.
224 * To allow relatively common schemes (although not common with Gromacs)
225 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
226 * we set the note limit to 10.
228 int min_steps_warn = 5;
229 int min_steps_note = 10;
232 gmx_moltype_t *moltype, *w_moltype;
234 t_ilist *ilist, *ilb, *ilc, *ils;
236 int i, a1, a2, w_a1, w_a2, j;
237 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
238 gmx_bool bFound, bWater, bWarn;
239 char warn_buf[STRLEN];
241 ip = mtop->ffparams.iparams;
243 twopi2 = sqr(2*M_PI);
245 limit2 = sqr(min_steps_note*dt);
251 for (molt = 0; molt < mtop->nmoltype; molt++)
253 moltype = &mtop->moltype[molt];
254 atom = moltype->atoms.atom;
255 ilist = moltype->ilist;
256 ilc = &ilist[F_CONSTR];
257 ils = &ilist[F_SETTLE];
258 for (ftype = 0; ftype < F_NRE; ftype++)
260 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
266 for (i = 0; i < ilb->nr; i += 3)
268 fc = ip[ilb->iatoms[i]].harmonic.krA;
269 re = ip[ilb->iatoms[i]].harmonic.rA;
270 if (ftype == F_G96BONDS)
272 /* Convert squared sqaure fc to harmonic fc */
275 a1 = ilb->iatoms[i+1];
276 a2 = ilb->iatoms[i+2];
279 if (fc > 0 && m1 > 0 && m2 > 0)
281 period2 = twopi2*m1*m2/((m1 + m2)*fc);
285 period2 = GMX_FLOAT_MAX;
289 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
290 fc, m1, m2, sqrt(period2));
292 if (period2 < limit2)
295 for (j = 0; j < ilc->nr; j += 3)
297 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
298 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
303 for (j = 0; j < ils->nr; j += 4)
305 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
306 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
312 (w_moltype == NULL || period2 < w_period2))
324 if (w_moltype != NULL)
326 bWarn = (w_period2 < sqr(min_steps_warn*dt));
327 /* A check that would recognize most water models */
328 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
329 w_moltype->atoms.nr <= 5);
330 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"
333 w_a1+1, *w_moltype->atoms.atomname[w_a1],
334 w_a2+1, *w_moltype->atoms.atomname[w_a2],
335 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
337 "Maybe you asked for fexible water." :
338 "Maybe you forgot to change the constraints mdp option.");
341 warning(wi, warn_buf);
345 warning_note(wi, warn_buf);
350 static void check_vel(gmx_mtop_t *mtop, rvec v[])
352 gmx_mtop_atomloop_all_t aloop;
356 aloop = gmx_mtop_atomloop_all_init(mtop);
357 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
359 if (atom->ptype == eptShell ||
360 atom->ptype == eptBond ||
361 atom->ptype == eptVSite)
368 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
373 for (mb = 0; mb < mtop->nmolblock; mb++)
375 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
381 /* This routine reorders the molecule type array
382 * in the order of use in the molblocks,
383 * unused molecule types are deleted.
385 static void renumber_moltypes(gmx_mtop_t *sys,
386 int *nmolinfo, t_molinfo **molinfo)
388 int *order, norder, i;
392 snew(order, *nmolinfo);
394 for (mb = 0; mb < sys->nmolblock; mb++)
396 for (i = 0; i < norder; i++)
398 if (order[i] == sys->molblock[mb].type)
405 /* This type did not occur yet, add it */
406 order[norder] = sys->molblock[mb].type;
407 /* Renumber the moltype in the topology */
410 sys->molblock[mb].type = i;
413 /* We still need to reorder the molinfo structs */
415 for (mi = 0; mi < *nmolinfo; mi++)
417 for (i = 0; i < norder; i++)
426 done_mi(&(*molinfo)[mi]);
430 minew[i] = (*molinfo)[mi];
439 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
444 mtop->nmoltype = nmi;
445 snew(mtop->moltype, nmi);
446 for (m = 0; m < nmi; m++)
448 molt = &mtop->moltype[m];
449 molt->name = mi[m].name;
450 molt->atoms = mi[m].atoms;
451 /* ilists are copied later */
452 molt->cgs = mi[m].cgs;
453 molt->excls = mi[m].excls;
458 new_status(const char *topfile, const char *topppfile, const char *confin,
459 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
460 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
461 gpp_atomtype_t atype, gmx_mtop_t *sys,
462 int *nmi, t_molinfo **mi, t_params plist[],
463 int *comb, double *reppow, real *fudgeQQ,
467 t_molinfo *molinfo = NULL;
469 gmx_molblock_t *molblock, *molbs;
471 int mb, i, nrmols, nmismatch;
473 gmx_bool bGB = FALSE;
474 char warn_buf[STRLEN];
478 /* Set gmx_boolean for GB */
479 if (ir->implicit_solvent)
484 /* TOPOLOGY processing */
485 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
486 plist, comb, reppow, fudgeQQ,
487 atype, &nrmols, &molinfo, ir,
488 &nmolblock, &molblock, bGB,
492 snew(sys->molblock, nmolblock);
495 for (mb = 0; mb < nmolblock; mb++)
497 if (sys->nmolblock > 0 &&
498 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
500 /* Merge consecutive blocks with the same molecule type */
501 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
502 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
504 else if (molblock[mb].nmol > 0)
506 /* Add a new molblock to the topology */
507 molbs = &sys->molblock[sys->nmolblock];
508 *molbs = molblock[mb];
509 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
510 molbs->nposres_xA = 0;
511 molbs->nposres_xB = 0;
512 sys->natoms += molbs->nmol*molbs->natoms_mol;
516 if (sys->nmolblock == 0)
518 gmx_fatal(FARGS, "No molecules were defined in the system");
521 renumber_moltypes(sys, &nrmols, &molinfo);
525 convert_harmonics(nrmols, molinfo, atype);
528 if (ir->eDisre == edrNone)
530 i = rm_interactions(F_DISRES, nrmols, molinfo);
533 set_warning_line(wi, "unknown", -1);
534 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
535 warning_note(wi, warn_buf);
538 if (opts->bOrire == FALSE)
540 i = rm_interactions(F_ORIRES, nrmols, molinfo);
543 set_warning_line(wi, "unknown", -1);
544 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
545 warning_note(wi, warn_buf);
549 /* Copy structures from msys to sys */
550 molinfo2mtop(nrmols, molinfo, sys);
552 gmx_mtop_finalize(sys);
554 /* COORDINATE file processing */
557 fprintf(stderr, "processing coordinates...\n");
560 get_stx_coordnum(confin, &state->natoms);
561 if (state->natoms != sys->natoms)
563 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
564 " does not match topology (%s, %d)",
565 confin, state->natoms, topfile, sys->natoms);
569 /* make space for coordinates and velocities */
572 init_t_atoms(confat, state->natoms, FALSE);
573 init_state(state, state->natoms, 0, 0, 0, 0);
574 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
575 /* This call fixes the box shape for runs with pressure scaling */
576 set_box_rel(ir, state);
578 nmismatch = check_atom_names(topfile, confin, sys, confat);
579 free_t_atoms(confat, TRUE);
584 sprintf(buf, "%d non-matching atom name%s\n"
585 "atom names from %s will be used\n"
586 "atom names from %s will be ignored\n",
587 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
592 fprintf(stderr, "double-checking input for internal consistency...\n");
594 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
600 gmx_mtop_atomloop_all_t aloop;
603 snew(mass, state->natoms);
604 aloop = gmx_mtop_atomloop_all_init(sys);
605 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
610 if (opts->seed == -1)
612 opts->seed = make_seed();
613 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
615 maxwell_speed(opts->tempi, opts->seed, sys, state->v);
617 stop_cm(stdout, state->natoms, mass, state->x, state->v);
625 static void copy_state(const char *slog, t_trxframe *fr,
626 gmx_bool bReadVel, t_state *state,
631 if (fr->not_ok & FRAME_NOT_OK)
633 gmx_fatal(FARGS, "Can not start from an incomplete frame");
637 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
641 for (i = 0; i < state->natoms; i++)
643 copy_rvec(fr->x[i], state->x[i]);
649 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
651 for (i = 0; i < state->natoms; i++)
653 copy_rvec(fr->v[i], state->v[i]);
658 copy_mat(fr->box, state->box);
661 *use_time = fr->time;
664 static void cont_status(const char *slog, const char *ener,
665 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
666 t_inputrec *ir, t_state *state,
668 const output_env_t oenv)
669 /* If fr_time == -1 read the last frame available which is complete */
677 bReadVel = (bNeedVel && !bGenVel);
680 "Reading Coordinates%s and Box size from old trajectory\n",
681 bReadVel ? ", Velocities" : "");
684 fprintf(stderr, "Will read whole trajectory\n");
688 fprintf(stderr, "Will read till time %g\n", fr_time);
694 fprintf(stderr, "Velocities generated: "
695 "ignoring velocities in input trajectory\n");
697 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
701 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
707 "WARNING: Did not find a frame with velocities in file %s,\n"
708 " all velocities will be set to zero!\n\n", slog);
709 for (i = 0; i < sys->natoms; i++)
711 clear_rvec(state->v[i]);
714 /* Search for a frame without velocities */
716 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
720 state->natoms = fr.natoms;
722 if (sys->natoms != state->natoms)
724 gmx_fatal(FARGS, "Number of atoms in Topology "
725 "is not the same as in Trajectory");
727 copy_state(slog, &fr, bReadVel, state, &use_time);
729 /* Find the appropriate frame */
730 while ((fr_time == -1 || fr.time < fr_time) &&
731 read_next_frame(oenv, fp, &fr))
733 copy_state(slog, &fr, bReadVel, state, &use_time);
738 /* Set the relative box lengths for preserving the box shape.
739 * Note that this call can lead to differences in the last bit
740 * with respect to using tpbconv to create a [TT].tpx[tt] file.
742 set_box_rel(ir, state);
744 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
745 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
747 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
749 get_enx_state(ener, use_time, &sys->groups, ir, state);
750 preserve_box_shape(ir, state->box_rel, state->boxv);
754 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
756 int rc_scaling, int ePBC,
760 gmx_bool bFirst = TRUE, *hadAtom;
766 int natoms, npbcdim = 0;
767 char warn_buf[STRLEN], title[STRLEN];
768 int a, i, ai, j, k, mb, nat_molb;
769 gmx_molblock_t *molb;
773 get_stx_coordnum(fn, &natoms);
774 if (natoms != mtop->natoms)
776 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);
777 warning(wi, warn_buf);
781 init_t_atoms(&dumat, natoms, FALSE);
782 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
784 npbcdim = ePBC2npbcdim(ePBC);
786 if (rc_scaling != erscNO)
788 copy_mat(box, invbox);
789 for (j = npbcdim; j < DIM; j++)
791 clear_rvec(invbox[j]);
794 m_inv_ur0(invbox, invbox);
797 /* Copy the reference coordinates to mtop */
801 snew(hadAtom, natoms);
802 for (mb = 0; mb < mtop->nmolblock; mb++)
804 molb = &mtop->molblock[mb];
805 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
806 pr = &(molinfo[molb->type].plist[F_POSRES]);
807 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
808 if (pr->nr > 0 || prfb->nr > 0)
810 atom = mtop->moltype[molb->type].atoms.atom;
811 for (i = 0; (i < pr->nr); i++)
813 ai = pr->param[i].AI;
816 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
817 ai+1, *molinfo[molb->type].name, fn, natoms);
820 if (rc_scaling == erscCOM)
822 /* Determine the center of mass of the posres reference coordinates */
823 for (j = 0; j < npbcdim; j++)
825 sum[j] += atom[ai].m*x[a+ai][j];
827 totmass += atom[ai].m;
830 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
831 for (i = 0; (i < prfb->nr); i++)
833 ai = prfb->param[i].AI;
836 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
837 ai+1, *molinfo[molb->type].name, fn, natoms);
839 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
841 /* Determine the center of mass of the posres reference coordinates */
842 for (j = 0; j < npbcdim; j++)
844 sum[j] += atom[ai].m*x[a+ai][j];
846 totmass += atom[ai].m;
851 molb->nposres_xA = nat_molb;
852 snew(molb->posres_xA, molb->nposres_xA);
853 for (i = 0; i < nat_molb; i++)
855 copy_rvec(x[a+i], molb->posres_xA[i]);
860 molb->nposres_xB = nat_molb;
861 snew(molb->posres_xB, molb->nposres_xB);
862 for (i = 0; i < nat_molb; i++)
864 copy_rvec(x[a+i], molb->posres_xB[i]);
870 if (rc_scaling == erscCOM)
874 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
876 for (j = 0; j < npbcdim; j++)
878 com[j] = sum[j]/totmass;
880 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]);
883 if (rc_scaling != erscNO)
885 for (mb = 0; mb < mtop->nmolblock; mb++)
887 molb = &mtop->molblock[mb];
888 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
889 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
891 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
892 for (i = 0; i < nat_molb; i++)
894 for (j = 0; j < npbcdim; j++)
896 if (rc_scaling == erscALL)
898 /* Convert from Cartesian to crystal coordinates */
899 xp[i][j] *= invbox[j][j];
900 for (k = j+1; k < npbcdim; k++)
902 xp[i][j] += invbox[k][j]*xp[i][k];
905 else if (rc_scaling == erscCOM)
907 /* Subtract the center of mass */
915 if (rc_scaling == erscCOM)
917 /* Convert the COM from Cartesian to crystal coordinates */
918 for (j = 0; j < npbcdim; j++)
920 com[j] *= invbox[j][j];
921 for (k = j+1; k < npbcdim; k++)
923 com[j] += invbox[k][j]*com[k];
929 free_t_atoms(&dumat, TRUE);
935 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
936 char *fnA, char *fnB,
937 int rc_scaling, int ePBC,
943 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
944 if (strcmp(fnA, fnB) != 0)
946 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
950 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
951 t_inputrec *ir, warninp_t wi)
954 char warn_buf[STRLEN];
958 fprintf(stderr, "Searching the wall atom type(s)\n");
960 for (i = 0; i < ir->nwall; i++)
962 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
963 if (ir->wall_atomtype[i] == NOTSET)
965 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
966 warning_error(wi, warn_buf);
971 static int nrdf_internal(t_atoms *atoms)
976 for (i = 0; i < atoms->nr; i++)
978 /* Vsite ptype might not be set here yet, so also check the mass */
979 if ((atoms->atom[i].ptype == eptAtom ||
980 atoms->atom[i].ptype == eptNucleus)
981 && atoms->atom[i].m > 0)
988 case 0: nrdf = 0; break;
989 case 1: nrdf = 0; break;
990 case 2: nrdf = 1; break;
991 default: nrdf = nmass*3 - 6; break;
1010 for (i = 1; i < n-1; i++)
1012 p = 0.5*y2[i-1]+2.0;
1014 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1015 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1020 for (i = n-2; i >= 0; i--)
1022 y2[i] = y2[i]*y2[i+1]+u[i];
1028 interpolate1d( double xmin,
1041 a = (xmin+(ix+1)*dx-x)/dx;
1042 b = (x-xmin-ix*dx)/dx;
1044 *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;
1045 *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];
1050 setup_cmap (int grid_spacing,
1053 gmx_cmap_t * cmap_grid)
1055 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1057 int i, j, k, ii, jj, kk, idx;
1059 double dx, xmin, v, v1, v2, v12;
1062 snew(tmp_u, 2*grid_spacing);
1063 snew(tmp_u2, 2*grid_spacing);
1064 snew(tmp_yy, 2*grid_spacing);
1065 snew(tmp_y1, 2*grid_spacing);
1066 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1067 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1069 dx = 360.0/grid_spacing;
1070 xmin = -180.0-dx*grid_spacing/2;
1072 for (kk = 0; kk < nc; kk++)
1074 /* Compute an offset depending on which cmap we are using
1075 * Offset will be the map number multiplied with the
1076 * grid_spacing * grid_spacing * 2
1078 offset = kk * grid_spacing * grid_spacing * 2;
1080 for (i = 0; i < 2*grid_spacing; i++)
1082 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1084 for (j = 0; j < 2*grid_spacing; j++)
1086 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1087 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1091 for (i = 0; i < 2*grid_spacing; i++)
1093 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1096 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1098 ii = i-grid_spacing/2;
1101 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1103 jj = j-grid_spacing/2;
1106 for (k = 0; k < 2*grid_spacing; k++)
1108 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1109 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1112 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1113 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1114 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1115 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1117 idx = ii*grid_spacing+jj;
1118 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1119 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1120 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1121 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1127 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1131 cmap_grid->ngrid = ngrid;
1132 cmap_grid->grid_spacing = grid_spacing;
1133 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1135 snew(cmap_grid->cmapdata, ngrid);
1137 for (i = 0; i < cmap_grid->ngrid; i++)
1139 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1144 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1146 int count, count_mol, i, mb;
1147 gmx_molblock_t *molb;
1152 for (mb = 0; mb < mtop->nmolblock; mb++)
1155 molb = &mtop->molblock[mb];
1156 plist = mi[molb->type].plist;
1158 for (i = 0; i < F_NRE; i++)
1162 count_mol += 3*plist[i].nr;
1164 else if (interaction_function[i].flags & IF_CONSTRAINT)
1166 count_mol += plist[i].nr;
1170 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1173 "Molecule type '%s' has %d constraints.\n"
1174 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1175 *mi[molb->type].name, count_mol,
1176 nrdf_internal(&mi[molb->type].atoms));
1179 count += molb->nmol*count_mol;
1185 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1187 int i, nmiss, natoms, mt;
1189 const t_atoms *atoms;
1192 for (mt = 0; mt < sys->nmoltype; mt++)
1194 atoms = &sys->moltype[mt].atoms;
1197 for (i = 0; i < natoms; i++)
1199 q = atoms->atom[i].q;
1200 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1201 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1202 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1203 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1204 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1207 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1208 get_atomtype_name(atoms->atom[i].type, atype), q);
1216 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);
1221 static void check_gbsa_params(gpp_atomtype_t atype)
1225 /* If we are doing GBSA, check that we got the parameters we need
1226 * This checking is to see if there are GBSA paratmeters for all
1227 * atoms in the force field. To go around this for testing purposes
1228 * comment out the nerror++ counter temporarily
1231 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1233 if (get_atomtype_radius(i, atype) < 0 ||
1234 get_atomtype_vol(i, atype) < 0 ||
1235 get_atomtype_surftens(i, atype) < 0 ||
1236 get_atomtype_gb_radius(i, atype) < 0 ||
1237 get_atomtype_S_hct(i, atype) < 0)
1239 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1240 get_atomtype_name(i, atype));
1247 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);
1252 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1255 real verletbuf_drift,
1260 verletbuf_list_setup_t ls;
1263 char warn_buf[STRLEN];
1266 for (i = 0; i < ir->opts.ngtc; i++)
1268 if (ir->opts.ref_t[i] < 0)
1270 warning(wi, "Some atom groups do not use temperature coupling. This cannot be accounted for in the energy drift estimation for the Verlet buffer size. The energy drift and the Verlet buffer might be underestimated.");
1274 ref_T = max(ref_T, ir->opts.ref_t[i]);
1278 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n", verletbuf_drift, ref_T);
1280 for (i = 0; i < ir->opts.ngtc; i++)
1282 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1284 sprintf(warn_buf, "ref_T for group of %.1f DOFs is %g K, which is smaller than the maximum of %g K used for the buffer size calculation. The buffer size might be on the conservative (large) side.",
1285 ir->opts.nrdf[i], ir->opts.ref_t[i], ref_T);
1286 warning_note(wi, warn_buf);
1290 /* Calculate the buffer size for simple atom vs atoms list */
1291 ls.cluster_size_i = 1;
1292 ls.cluster_size_j = 1;
1293 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1294 &ls, &n_nonlin_vsite, &rlist_1x1);
1296 /* Set the pair-list buffer size in ir */
1297 verletbuf_get_list_setup(FALSE, &ls);
1298 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1299 &ls, &n_nonlin_vsite, &ir->rlist);
1301 if (n_nonlin_vsite > 0)
1303 sprintf(warn_buf, "There are %d non-linear virtual site constructions. Their contribution to the energy drift is approximated. In most cases this does not affect the energy drift significantly.", n_nonlin_vsite);
1304 warning_note(wi, warn_buf);
1307 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1308 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1310 ir->rlistlong = ir->rlist;
1311 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1312 ls.cluster_size_i, ls.cluster_size_j,
1313 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1315 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1317 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-drift.", ir->rlistlong, sqrt(max_cutoff2(ir->ePBC, box)));
1321 int gmx_grompp(int argc, char *argv[])
1323 static const char *desc[] = {
1324 "[THISMODULE] (the gromacs preprocessor)",
1325 "reads a molecular topology file, checks the validity of the",
1326 "file, expands the topology from a molecular description to an atomic",
1327 "description. The topology file contains information about",
1328 "molecule types and the number of molecules, the preprocessor",
1329 "copies each molecule as needed. ",
1330 "There is no limitation on the number of molecule types. ",
1331 "Bonds and bond-angles can be converted into constraints, separately",
1332 "for hydrogens and heavy atoms.",
1333 "Then a coordinate file is read and velocities can be generated",
1334 "from a Maxwellian distribution if requested.",
1335 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1336 "(eg. number of MD steps, time step, cut-off), and others such as",
1337 "NEMD parameters, which are corrected so that the net acceleration",
1339 "Eventually a binary file is produced that can serve as the sole input",
1340 "file for the MD program.[PAR]",
1342 "[THISMODULE] uses the atom names from the topology file. The atom names",
1343 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1344 "warnings when they do not match the atom names in the topology.",
1345 "Note that the atom names are irrelevant for the simulation as",
1346 "only the atom types are used for generating interaction parameters.[PAR]",
1348 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1349 "etc. The preprocessor supports the following keywords:[PAR]",
1350 "#ifdef VARIABLE[BR]",
1351 "#ifndef VARIABLE[BR]",
1354 "#define VARIABLE[BR]",
1355 "#undef VARIABLE[BR]"
1356 "#include \"filename\"[BR]",
1357 "#include <filename>[PAR]",
1358 "The functioning of these statements in your topology may be modulated by",
1359 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1360 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1361 "include = -I/home/john/doe[tt][BR]",
1362 "For further information a C-programming textbook may help you out.",
1363 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1364 "topology file written out so that you can verify its contents.[PAR]",
1366 /* cpp has been unnecessary for some time, hasn't it?
1367 "If your system does not have a C-preprocessor, you can still",
1368 "use [TT]grompp[tt], but you do not have access to the features ",
1369 "from the cpp. Command line options to the C-preprocessor can be given",
1370 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1373 "When using position restraints a file with restraint coordinates",
1374 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1375 "with respect to the conformation from the [TT]-c[tt] option.",
1376 "For free energy calculation the the coordinates for the B topology",
1377 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1378 "those of the A topology.[PAR]",
1380 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1381 "The last frame with coordinates and velocities will be read,",
1382 "unless the [TT]-time[tt] option is used. Only if this information",
1383 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1384 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1385 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1386 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1389 "[THISMODULE] can be used to restart simulations (preserving",
1390 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1391 "However, for simply changing the number of run steps to extend",
1392 "a run, using [gmx-tpbconv] is more convenient than [THISMODULE].",
1393 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1394 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1395 "like output frequency, then supplying the checkpoint file to",
1396 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1397 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1399 "By default, all bonded interactions which have constant energy due to",
1400 "virtual site constructions will be removed. If this constant energy is",
1401 "not zero, this will result in a shift in the total energy. All bonded",
1402 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1403 "all constraints for distances which will be constant anyway because",
1404 "of virtual site constructions will be removed. If any constraints remain",
1405 "which involve virtual sites, a fatal error will result.[PAR]"
1407 "To verify your run input file, please take note of all warnings",
1408 "on the screen, and correct where necessary. Do also look at the contents",
1409 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1410 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1411 "with the [TT]-debug[tt] option which will give you more information",
1412 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1413 "can see the contents of the run input file with the [gmx-dump]",
1414 "program. [gmx-check] can be used to compare the contents of two",
1415 "run input files.[PAR]"
1417 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1418 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1419 "harmless, but usually they are not. The user is advised to carefully",
1420 "interpret the output messages before attempting to bypass them with",
1427 gpp_atomtype_t atype;
1429 int natoms, nvsite, comb, mt;
1433 real max_spacing, fudgeQQ;
1435 char fn[STRLEN], fnB[STRLEN];
1436 const char *mdparin;
1438 gmx_bool bNeedVel, bGenVel;
1439 gmx_bool have_atomnumber;
1441 t_params *gb_plist = NULL;
1442 gmx_genborn_t *born = NULL;
1444 gmx_bool bVerbose = FALSE;
1446 char warn_buf[STRLEN];
1449 { efMDP, NULL, NULL, ffREAD },
1450 { efMDP, "-po", "mdout", ffWRITE },
1451 { efSTX, "-c", NULL, ffREAD },
1452 { efSTX, "-r", NULL, ffOPTRD },
1453 { efSTX, "-rb", NULL, ffOPTRD },
1454 { efNDX, NULL, NULL, ffOPTRD },
1455 { efTOP, NULL, NULL, ffREAD },
1456 { efTOP, "-pp", "processed", ffOPTWR },
1457 { efTPX, "-o", NULL, ffWRITE },
1458 { efTRN, "-t", NULL, ffOPTRD },
1459 { efEDR, "-e", NULL, ffOPTRD },
1460 { efTRN, "-ref", "rotref", ffOPTRW }
1462 #define NFILE asize(fnm)
1464 /* Command line options */
1465 static gmx_bool bRenum = TRUE;
1466 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1467 static int i, maxwarn = 0;
1468 static real fr_time = -1;
1470 { "-v", FALSE, etBOOL, {&bVerbose},
1471 "Be loud and noisy" },
1472 { "-time", FALSE, etREAL, {&fr_time},
1473 "Take frame at or first after this time." },
1474 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1475 "Remove constant bonded interactions with virtual sites" },
1476 { "-maxwarn", FALSE, etINT, {&maxwarn},
1477 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1478 { "-zero", FALSE, etBOOL, {&bZero},
1479 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1480 { "-renum", FALSE, etBOOL, {&bRenum},
1481 "Renumber atomtypes and minimize number of atomtypes" }
1484 /* Initiate some variables */
1489 /* Parse the command line */
1490 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1491 asize(desc), desc, 0, NULL, &oenv))
1496 wi = init_warning(TRUE, maxwarn);
1498 /* PARAMETER file processing */
1499 mdparin = opt2fn("-f", NFILE, fnm);
1500 set_warning_line(wi, mdparin, -1);
1501 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1505 fprintf(stderr, "checking input for internal consistency...\n");
1507 check_ir(mdparin, ir, opts, wi);
1509 if (ir->ld_seed == -1)
1511 ir->ld_seed = make_seed();
1512 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1515 if (ir->expandedvals->lmc_seed == -1)
1517 ir->expandedvals->lmc_seed = make_seed();
1518 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1521 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1522 bGenVel = (bNeedVel && opts->bGenVel);
1523 if (bGenVel && ir->bContinuation)
1526 "Generating velocities is inconsistent with attempting "
1527 "to continue a previous run. Choose only one of "
1528 "gen-vel = yes and continuation = yes.");
1529 warning_error(wi, warn_buf);
1535 atype = init_atomtype();
1538 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1541 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1542 if (!gmx_fexist(fn))
1544 gmx_fatal(FARGS, "%s does not exist", fn);
1546 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1547 opts, ir, bZero, bGenVel, bVerbose, &state,
1548 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1554 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1558 /* set parameters for virtual site construction (not for vsiten) */
1559 for (mt = 0; mt < sys->nmoltype; mt++)
1562 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1564 /* now throw away all obsolete bonds, angles and dihedrals: */
1565 /* note: constraints are ALWAYS removed */
1568 for (mt = 0; mt < sys->nmoltype; mt++)
1570 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1574 if (ir->cutoff_scheme == ecutsVERLET)
1576 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1577 ecutscheme_names[ir->cutoff_scheme]);
1579 /* Remove all charge groups */
1580 gmx_mtop_remove_chargegroups(sys);
1582 if (EVDW_PME(ir->vdwtype))
1584 gmx_fatal(FARGS, "LJ-PME not implemented together with verlet-scheme!");
1588 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1590 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1592 sprintf(warn_buf, "Can not do %s with %s, use %s",
1593 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1594 warning_error(wi, warn_buf);
1596 if (ir->bPeriodicMols)
1598 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1599 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1600 warning_error(wi, warn_buf);
1604 if (EI_SD (ir->eI) && ir->etc != etcNO)
1606 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1609 /* If we are doing QM/MM, check that we got the atom numbers */
1610 have_atomnumber = TRUE;
1611 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1613 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1615 if (!have_atomnumber && ir->bQMMM)
1619 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1620 "field you are using does not contain atom numbers fields. This is an\n"
1621 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1622 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1623 "an integer just before the mass column in ffXXXnb.itp.\n"
1624 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1629 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1631 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1633 /** \TODO check size of ex+hy width against box size */
1636 /* Check for errors in the input now, since they might cause problems
1637 * during processing further down.
1639 check_warning_error(wi, FARGS);
1641 if (opt2bSet("-r", NFILE, fnm))
1643 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1647 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1649 if (opt2bSet("-rb", NFILE, fnm))
1651 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1658 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1662 fprintf(stderr, "Reading position restraint coords from %s", fn);
1663 if (strcmp(fn, fnB) == 0)
1665 fprintf(stderr, "\n");
1669 fprintf(stderr, " and %s\n", fnB);
1672 gen_posres(sys, mi, fn, fnB,
1673 ir->refcoord_scaling, ir->ePBC,
1674 ir->posres_com, ir->posres_comB,
1678 /* If we are using CMAP, setup the pre-interpolation grid */
1679 if (plist->ncmap > 0)
1681 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1682 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1685 set_wall_atomtype(atype, opts, ir, wi);
1688 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1689 ntype = get_atomtype_ntypes(atype);
1692 if (ir->implicit_solvent != eisNO)
1694 /* Now we have renumbered the atom types, we can check the GBSA params */
1695 check_gbsa_params(atype);
1697 /* Check that all atoms that have charge and/or LJ-parameters also have
1698 * sensible GB-parameters
1700 check_gbsa_params_charged(sys, atype);
1703 /* PELA: Copy the atomtype data to the topology atomtype list */
1704 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1708 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1713 fprintf(stderr, "converting bonded parameters...\n");
1716 ntype = get_atomtype_ntypes(atype);
1717 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1721 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1724 /* set ptype to VSite for virtual sites */
1725 for (mt = 0; mt < sys->nmoltype; mt++)
1727 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1731 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1733 /* Check velocity for virtual sites and shells */
1736 check_vel(sys, state.v);
1742 for (i = 0; i < sys->nmoltype; i++)
1744 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1747 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1749 check_bonds_timestep(sys, ir->delta_t, wi);
1752 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1754 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.");
1757 check_warning_error(wi, FARGS);
1761 fprintf(stderr, "initialising group options...\n");
1763 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1765 bGenVel ? state.v : NULL,
1768 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1771 if (EI_DYNAMICS(ir->eI) &&
1772 !(EI_MD(ir->eI) && ir->etc == etcNO) &&
1773 inputrec2nboundeddim(ir) == 3)
1775 set_verlet_buffer(sys, ir, state.box, ir->verletbuf_drift, wi);
1779 /* Init the temperature coupling state */
1780 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1784 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1786 check_eg_vs_cg(sys);
1790 pr_symtab(debug, 0, "After index", &sys->symtab);
1793 triple_check(mdparin, ir, sys, wi);
1794 close_symtab(&sys->symtab);
1797 pr_symtab(debug, 0, "After close", &sys->symtab);
1800 /* make exclusions between QM atoms */
1803 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1805 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1809 generate_qmexcl(sys, ir, wi);
1813 if (ftp2bSet(efTRN, NFILE, fnm))
1817 fprintf(stderr, "getting data from old trajectory ...\n");
1819 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1820 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1823 if (ir->ePBC == epbcXY && ir->nwall != 2)
1825 clear_rvec(state.box[ZZ]);
1828 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1830 set_warning_line(wi, mdparin, -1);
1831 check_chargegroup_radii(sys, ir, state.x, wi);
1834 if (EEL_FULL(ir->coulombtype) || EVDW_PME(ir->vdwtype))
1836 /* Calculate the optimal grid dimensions */
1837 copy_mat(state.box, box);
1838 if (ir->ePBC == epbcXY && ir->nwall == 2)
1840 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1842 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1844 /* Mark fourier_spacing as not used */
1845 ir->fourier_spacing = 0;
1847 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1849 set_warning_line(wi, mdparin, -1);
1850 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1852 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1853 &(ir->nkx), &(ir->nky), &(ir->nkz));
1856 /* MRS: eventually figure out better logic for initializing the fep
1857 values that makes declaring the lambda and declaring the state not
1858 potentially conflict if not handled correctly. */
1859 if (ir->efep != efepNO)
1861 if (EVDW_PME(ir->vdwtype))
1863 gmx_fatal(FARGS, "LJ-PME not implemented together with free energy calculations!");
1866 state.fep_state = ir->fepvals->init_fep_state;
1867 for (i = 0; i < efptNR; i++)
1869 /* init_lambda trumps state definitions*/
1870 if (ir->fepvals->init_lambda >= 0)
1872 state.lambda[i] = ir->fepvals->init_lambda;
1876 if (ir->fepvals->all_lambda[i] == NULL)
1878 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1882 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1888 if (ir->ePull != epullNO)
1890 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1895 set_reference_positions(ir->rot, state.x, state.box,
1896 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1900 /* reset_multinr(sys); */
1902 if (EEL_PME(ir->coulombtype))
1904 float ratio = pme_load_estimate(sys, ir, state.box);
1905 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
1906 /* With free energy we might need to do PME both for the A and B state
1907 * charges. This will double the cost, but the optimal performance will
1908 * then probably be at a slightly larger cut-off and grid spacing.
1910 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1911 (ir->efep != efepNO && ratio > 2.0/3.0))
1914 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1915 "and for highly parallel simulations between 0.25 and 0.33,\n"
1916 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1917 if (ir->efep != efepNO)
1920 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1926 char warn_buf[STRLEN];
1927 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
1928 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
1931 set_warning_line(wi, mdparin, -1);
1932 warning_note(wi, warn_buf);
1936 printf("%s\n", warn_buf);
1942 fprintf(stderr, "writing run input file...\n");
1945 done_warning(wi, FARGS);
1947 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);