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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
42 #include <sys/types.h>
66 #include "sortwater.h"
68 #include "gmx_fatal.h"
74 #include "vsite_parm.h"
80 #include "compute_io.h"
81 #include "gpp_atomtype.h"
82 #include "gpp_tomorse.h"
83 #include "mtop_util.h"
85 #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 void check_shells_inputrec(gmx_mtop_t *mtop,
372 gmx_mtop_atomloop_all_t aloop;
375 char warn_buf[STRLEN];
377 aloop = gmx_mtop_atomloop_all_init(mtop);
378 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
380 if (atom->ptype == eptShell ||
381 atom->ptype == eptBond)
386 if (IR_TWINRANGE(*ir) && (nshells > 0))
388 snprintf(warn_buf, STRLEN,
389 "The combination of using shells and a twin-range cut-off is not supported");
390 warning_error(wi, warn_buf);
392 if ((nshells > 0) && (ir->nstcalcenergy != 1))
394 set_warning_line(wi, "unknown", -1);
395 snprintf(warn_buf, STRLEN,
396 "There are %d shells, changing nstcalcenergy from %d to 1",
397 nshells, ir->nstcalcenergy);
398 ir->nstcalcenergy = 1;
399 warning(wi, warn_buf);
403 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
408 for (mb = 0; mb < mtop->nmolblock; mb++)
410 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
416 /* This routine reorders the molecule type array
417 * in the order of use in the molblocks,
418 * unused molecule types are deleted.
420 static void renumber_moltypes(gmx_mtop_t *sys,
421 int *nmolinfo, t_molinfo **molinfo)
423 int *order, norder, i;
427 snew(order, *nmolinfo);
429 for (mb = 0; mb < sys->nmolblock; mb++)
431 for (i = 0; i < norder; i++)
433 if (order[i] == sys->molblock[mb].type)
440 /* This type did not occur yet, add it */
441 order[norder] = sys->molblock[mb].type;
442 /* Renumber the moltype in the topology */
445 sys->molblock[mb].type = i;
448 /* We still need to reorder the molinfo structs */
450 for (mi = 0; mi < *nmolinfo; mi++)
452 for (i = 0; i < norder; i++)
461 done_mi(&(*molinfo)[mi]);
465 minew[i] = (*molinfo)[mi];
474 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
479 mtop->nmoltype = nmi;
480 snew(mtop->moltype, nmi);
481 for (m = 0; m < nmi; m++)
483 molt = &mtop->moltype[m];
484 molt->name = mi[m].name;
485 molt->atoms = mi[m].atoms;
486 /* ilists are copied later */
487 molt->cgs = mi[m].cgs;
488 molt->excls = mi[m].excls;
493 new_status(const char *topfile, const char *topppfile, const char *confin,
494 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
495 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
496 gpp_atomtype_t atype, gmx_mtop_t *sys,
497 int *nmi, t_molinfo **mi, t_params plist[],
498 int *comb, double *reppow, real *fudgeQQ,
502 t_molinfo *molinfo = NULL;
504 gmx_molblock_t *molblock, *molbs;
506 int mb, i, nrmols, nmismatch;
508 gmx_bool bGB = FALSE;
509 char warn_buf[STRLEN];
513 /* Set gmx_boolean for GB */
514 if (ir->implicit_solvent)
519 /* TOPOLOGY processing */
520 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
521 plist, comb, reppow, fudgeQQ,
522 atype, &nrmols, &molinfo, ir,
523 &nmolblock, &molblock, bGB,
527 snew(sys->molblock, nmolblock);
530 for (mb = 0; mb < nmolblock; mb++)
532 if (sys->nmolblock > 0 &&
533 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
535 /* Merge consecutive blocks with the same molecule type */
536 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
537 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
539 else if (molblock[mb].nmol > 0)
541 /* Add a new molblock to the topology */
542 molbs = &sys->molblock[sys->nmolblock];
543 *molbs = molblock[mb];
544 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
545 molbs->nposres_xA = 0;
546 molbs->nposres_xB = 0;
547 sys->natoms += molbs->nmol*molbs->natoms_mol;
551 if (sys->nmolblock == 0)
553 gmx_fatal(FARGS, "No molecules were defined in the system");
556 renumber_moltypes(sys, &nrmols, &molinfo);
560 convert_harmonics(nrmols, molinfo, atype);
563 if (ir->eDisre == edrNone)
565 i = rm_interactions(F_DISRES, nrmols, molinfo);
568 set_warning_line(wi, "unknown", -1);
569 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
570 warning_note(wi, warn_buf);
573 if (opts->bOrire == FALSE)
575 i = rm_interactions(F_ORIRES, nrmols, molinfo);
578 set_warning_line(wi, "unknown", -1);
579 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
580 warning_note(wi, warn_buf);
584 /* Copy structures from msys to sys */
585 molinfo2mtop(nrmols, molinfo, sys);
587 gmx_mtop_finalize(sys);
589 /* COORDINATE file processing */
592 fprintf(stderr, "processing coordinates...\n");
595 get_stx_coordnum(confin, &state->natoms);
596 if (state->natoms != sys->natoms)
598 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
599 " does not match topology (%s, %d)",
600 confin, state->natoms, topfile, sys->natoms);
604 /* make space for coordinates and velocities */
607 init_t_atoms(confat, state->natoms, FALSE);
608 init_state(state, state->natoms, 0, 0, 0, 0);
609 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
610 /* This call fixes the box shape for runs with pressure scaling */
611 set_box_rel(ir, state);
613 nmismatch = check_atom_names(topfile, confin, sys, confat);
614 free_t_atoms(confat, TRUE);
619 sprintf(buf, "%d non-matching atom name%s\n"
620 "atom names from %s will be used\n"
621 "atom names from %s will be ignored\n",
622 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
627 fprintf(stderr, "double-checking input for internal consistency...\n");
629 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
635 gmx_mtop_atomloop_all_t aloop;
638 snew(mass, state->natoms);
639 aloop = gmx_mtop_atomloop_all_init(sys);
640 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
645 if (opts->seed == -1)
647 opts->seed = make_seed();
648 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
650 maxwell_speed(opts->tempi, opts->seed, sys, state->v);
652 stop_cm(stdout, state->natoms, mass, state->x, state->v);
660 static void copy_state(const char *slog, t_trxframe *fr,
661 gmx_bool bReadVel, t_state *state,
666 if (fr->not_ok & FRAME_NOT_OK)
668 gmx_fatal(FARGS, "Can not start from an incomplete frame");
672 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
676 for (i = 0; i < state->natoms; i++)
678 copy_rvec(fr->x[i], state->x[i]);
684 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
686 for (i = 0; i < state->natoms; i++)
688 copy_rvec(fr->v[i], state->v[i]);
693 copy_mat(fr->box, state->box);
696 *use_time = fr->time;
699 static void cont_status(const char *slog, const char *ener,
700 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
701 t_inputrec *ir, t_state *state,
703 const output_env_t oenv)
704 /* If fr_time == -1 read the last frame available which is complete */
712 bReadVel = (bNeedVel && !bGenVel);
715 "Reading Coordinates%s and Box size from old trajectory\n",
716 bReadVel ? ", Velocities" : "");
719 fprintf(stderr, "Will read whole trajectory\n");
723 fprintf(stderr, "Will read till time %g\n", fr_time);
729 fprintf(stderr, "Velocities generated: "
730 "ignoring velocities in input trajectory\n");
732 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
736 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
742 "WARNING: Did not find a frame with velocities in file %s,\n"
743 " all velocities will be set to zero!\n\n", slog);
744 for (i = 0; i < sys->natoms; i++)
746 clear_rvec(state->v[i]);
749 /* Search for a frame without velocities */
751 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
755 state->natoms = fr.natoms;
757 if (sys->natoms != state->natoms)
759 gmx_fatal(FARGS, "Number of atoms in Topology "
760 "is not the same as in Trajectory");
762 copy_state(slog, &fr, bReadVel, state, &use_time);
764 /* Find the appropriate frame */
765 while ((fr_time == -1 || fr.time < fr_time) &&
766 read_next_frame(oenv, fp, &fr))
768 copy_state(slog, &fr, bReadVel, state, &use_time);
773 /* Set the relative box lengths for preserving the box shape.
774 * Note that this call can lead to differences in the last bit
775 * with respect to using tpbconv to create a [TT].tpx[tt] file.
777 set_box_rel(ir, state);
779 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
780 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
782 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
784 get_enx_state(ener, use_time, &sys->groups, ir, state);
785 preserve_box_shape(ir, state->box_rel, state->boxv);
789 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
791 int rc_scaling, int ePBC,
795 gmx_bool bFirst = TRUE;
801 int natoms, npbcdim = 0;
802 char warn_buf[STRLEN], title[STRLEN];
803 int a, i, ai, j, k, mb, nat_molb;
804 gmx_molblock_t *molb;
808 get_stx_coordnum(fn, &natoms);
809 if (natoms != mtop->natoms)
811 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);
812 warning(wi, warn_buf);
816 init_t_atoms(&dumat, natoms, FALSE);
817 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
819 npbcdim = ePBC2npbcdim(ePBC);
821 if (rc_scaling != erscNO)
823 copy_mat(box, invbox);
824 for (j = npbcdim; j < DIM; j++)
826 clear_rvec(invbox[j]);
829 m_inv_ur0(invbox, invbox);
832 /* Copy the reference coordinates to mtop */
836 for (mb = 0; mb < mtop->nmolblock; mb++)
838 molb = &mtop->molblock[mb];
839 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
840 pr = &(molinfo[molb->type].plist[F_POSRES]);
843 atom = mtop->moltype[molb->type].atoms.atom;
844 for (i = 0; (i < pr->nr); i++)
846 ai = pr->param[i].AI;
849 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
850 ai+1, *molinfo[molb->type].name, fn, natoms);
852 if (rc_scaling == erscCOM)
854 /* Determine the center of mass of the posres reference coordinates */
855 for (j = 0; j < npbcdim; j++)
857 sum[j] += atom[ai].m*x[a+ai][j];
859 totmass += atom[ai].m;
864 molb->nposres_xA = nat_molb;
865 snew(molb->posres_xA, molb->nposres_xA);
866 for (i = 0; i < nat_molb; i++)
868 copy_rvec(x[a+i], molb->posres_xA[i]);
873 molb->nposres_xB = nat_molb;
874 snew(molb->posres_xB, molb->nposres_xB);
875 for (i = 0; i < nat_molb; i++)
877 copy_rvec(x[a+i], molb->posres_xB[i]);
883 if (rc_scaling == erscCOM)
887 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
889 for (j = 0; j < npbcdim; j++)
891 com[j] = sum[j]/totmass;
893 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]);
896 if (rc_scaling != erscNO)
898 for (mb = 0; mb < mtop->nmolblock; mb++)
900 molb = &mtop->molblock[mb];
901 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
902 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
904 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
905 for (i = 0; i < nat_molb; i++)
907 for (j = 0; j < npbcdim; j++)
909 if (rc_scaling == erscALL)
911 /* Convert from Cartesian to crystal coordinates */
912 xp[i][j] *= invbox[j][j];
913 for (k = j+1; k < npbcdim; k++)
915 xp[i][j] += invbox[k][j]*xp[i][k];
918 else if (rc_scaling == erscCOM)
920 /* Subtract the center of mass */
928 if (rc_scaling == erscCOM)
930 /* Convert the COM from Cartesian to crystal coordinates */
931 for (j = 0; j < npbcdim; j++)
933 com[j] *= invbox[j][j];
934 for (k = j+1; k < npbcdim; k++)
936 com[j] += invbox[k][j]*com[k];
942 free_t_atoms(&dumat, TRUE);
947 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
948 char *fnA, char *fnB,
949 int rc_scaling, int ePBC,
955 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
956 if (strcmp(fnA, fnB) != 0)
958 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
962 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
963 t_inputrec *ir, warninp_t wi)
966 char warn_buf[STRLEN];
970 fprintf(stderr, "Searching the wall atom type(s)\n");
972 for (i = 0; i < ir->nwall; i++)
974 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
975 if (ir->wall_atomtype[i] == NOTSET)
977 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
978 warning_error(wi, warn_buf);
983 static int nrdf_internal(t_atoms *atoms)
988 for (i = 0; i < atoms->nr; i++)
990 /* Vsite ptype might not be set here yet, so also check the mass */
991 if ((atoms->atom[i].ptype == eptAtom ||
992 atoms->atom[i].ptype == eptNucleus)
993 && atoms->atom[i].m > 0)
1000 case 0: nrdf = 0; break;
1001 case 1: nrdf = 0; break;
1002 case 2: nrdf = 1; break;
1003 default: nrdf = nmass*3 - 6; break;
1010 spline1d( double dx,
1022 for (i = 1; i < n-1; i++)
1024 p = 0.5*y2[i-1]+2.0;
1026 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1027 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1032 for (i = n-2; i >= 0; i--)
1034 y2[i] = y2[i]*y2[i+1]+u[i];
1040 interpolate1d( double xmin,
1053 a = (xmin+(ix+1)*dx-x)/dx;
1054 b = (x-xmin-ix*dx)/dx;
1056 *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;
1057 *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];
1062 setup_cmap (int grid_spacing,
1065 gmx_cmap_t * cmap_grid)
1067 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1069 int i, j, k, ii, jj, kk, idx;
1071 double dx, xmin, v, v1, v2, v12;
1074 snew(tmp_u, 2*grid_spacing);
1075 snew(tmp_u2, 2*grid_spacing);
1076 snew(tmp_yy, 2*grid_spacing);
1077 snew(tmp_y1, 2*grid_spacing);
1078 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1079 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1081 dx = 360.0/grid_spacing;
1082 xmin = -180.0-dx*grid_spacing/2;
1084 for (kk = 0; kk < nc; kk++)
1086 /* Compute an offset depending on which cmap we are using
1087 * Offset will be the map number multiplied with the
1088 * grid_spacing * grid_spacing * 2
1090 offset = kk * grid_spacing * grid_spacing * 2;
1092 for (i = 0; i < 2*grid_spacing; i++)
1094 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1096 for (j = 0; j < 2*grid_spacing; j++)
1098 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1099 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1103 for (i = 0; i < 2*grid_spacing; i++)
1105 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1108 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1110 ii = i-grid_spacing/2;
1113 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1115 jj = j-grid_spacing/2;
1118 for (k = 0; k < 2*grid_spacing; k++)
1120 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1121 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1124 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1125 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1126 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1127 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1129 idx = ii*grid_spacing+jj;
1130 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1131 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1132 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1133 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1139 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1143 cmap_grid->ngrid = ngrid;
1144 cmap_grid->grid_spacing = grid_spacing;
1145 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1147 snew(cmap_grid->cmapdata, ngrid);
1149 for (i = 0; i < cmap_grid->ngrid; i++)
1151 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1156 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1158 int count, count_mol, i, mb;
1159 gmx_molblock_t *molb;
1164 for (mb = 0; mb < mtop->nmolblock; mb++)
1167 molb = &mtop->molblock[mb];
1168 plist = mi[molb->type].plist;
1170 for (i = 0; i < F_NRE; i++)
1174 count_mol += 3*plist[i].nr;
1176 else if (interaction_function[i].flags & IF_CONSTRAINT)
1178 count_mol += plist[i].nr;
1182 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1185 "Molecule type '%s' has %d constraints.\n"
1186 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1187 *mi[molb->type].name, count_mol,
1188 nrdf_internal(&mi[molb->type].atoms));
1191 count += molb->nmol*count_mol;
1197 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1199 int i, nmiss, natoms, mt;
1201 const t_atoms *atoms;
1204 for (mt = 0; mt < sys->nmoltype; mt++)
1206 atoms = &sys->moltype[mt].atoms;
1209 for (i = 0; i < natoms; i++)
1211 q = atoms->atom[i].q;
1212 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1213 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1214 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1215 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1216 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1219 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1220 get_atomtype_name(atoms->atom[i].type, atype), q);
1228 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);
1233 static void check_gbsa_params(t_inputrec *ir, gpp_atomtype_t atype)
1237 /* If we are doing GBSA, check that we got the parameters we need
1238 * This checking is to see if there are GBSA paratmeters for all
1239 * atoms in the force field. To go around this for testing purposes
1240 * comment out the nerror++ counter temporarily
1243 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1245 if (get_atomtype_radius(i, atype) < 0 ||
1246 get_atomtype_vol(i, atype) < 0 ||
1247 get_atomtype_surftens(i, atype) < 0 ||
1248 get_atomtype_gb_radius(i, atype) < 0 ||
1249 get_atomtype_S_hct(i, atype) < 0)
1251 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1252 get_atomtype_name(i, atype));
1259 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);
1264 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1267 real verletbuf_drift,
1272 verletbuf_list_setup_t ls;
1275 char warn_buf[STRLEN];
1278 for (i = 0; i < ir->opts.ngtc; i++)
1280 if (ir->opts.ref_t[i] < 0)
1282 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.");
1286 ref_T = max(ref_T, ir->opts.ref_t[i]);
1290 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n", verletbuf_drift, ref_T);
1292 for (i = 0; i < ir->opts.ngtc; i++)
1294 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1296 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.",
1297 ir->opts.nrdf[i], ir->opts.ref_t[i], ref_T);
1298 warning_note(wi, warn_buf);
1302 /* Calculate the buffer size for simple atom vs atoms list */
1303 ls.cluster_size_i = 1;
1304 ls.cluster_size_j = 1;
1305 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1306 &ls, &n_nonlin_vsite, &rlist_1x1);
1308 /* Set the pair-list buffer size in ir */
1309 verletbuf_get_list_setup(FALSE, &ls);
1310 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1311 &ls, &n_nonlin_vsite, &ir->rlist);
1313 if (n_nonlin_vsite > 0)
1315 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);
1316 warning_note(wi, warn_buf);
1319 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1320 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1322 ir->rlistlong = ir->rlist;
1323 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1324 ls.cluster_size_i, ls.cluster_size_j,
1325 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1327 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1329 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)));
1333 int cmain (int argc, char *argv[])
1335 static const char *desc[] = {
1336 "The gromacs preprocessor",
1337 "reads a molecular topology file, checks the validity of the",
1338 "file, expands the topology from a molecular description to an atomic",
1339 "description. The topology file contains information about",
1340 "molecule types and the number of molecules, the preprocessor",
1341 "copies each molecule as needed. ",
1342 "There is no limitation on the number of molecule types. ",
1343 "Bonds and bond-angles can be converted into constraints, separately",
1344 "for hydrogens and heavy atoms.",
1345 "Then a coordinate file is read and velocities can be generated",
1346 "from a Maxwellian distribution if requested.",
1347 "[TT]grompp[tt] also reads parameters for the [TT]mdrun[tt] ",
1348 "(eg. number of MD steps, time step, cut-off), and others such as",
1349 "NEMD parameters, which are corrected so that the net acceleration",
1351 "Eventually a binary file is produced that can serve as the sole input",
1352 "file for the MD program.[PAR]",
1354 "[TT]grompp[tt] uses the atom names from the topology file. The atom names",
1355 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1356 "warnings when they do not match the atom names in the topology.",
1357 "Note that the atom names are irrelevant for the simulation as",
1358 "only the atom types are used for generating interaction parameters.[PAR]",
1360 "[TT]grompp[tt] uses a built-in preprocessor to resolve includes, macros, ",
1361 "etc. The preprocessor supports the following keywords:[PAR]",
1362 "#ifdef VARIABLE[BR]",
1363 "#ifndef VARIABLE[BR]",
1366 "#define VARIABLE[BR]",
1367 "#undef VARIABLE[BR]"
1368 "#include \"filename\"[BR]",
1369 "#include <filename>[PAR]",
1370 "The functioning of these statements in your topology may be modulated by",
1371 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1372 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1373 "include = -I/home/john/doe[tt][BR]",
1374 "For further information a C-programming textbook may help you out.",
1375 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1376 "topology file written out so that you can verify its contents.[PAR]",
1378 /* cpp has been unnecessary for some time, hasn't it?
1379 "If your system does not have a C-preprocessor, you can still",
1380 "use [TT]grompp[tt], but you do not have access to the features ",
1381 "from the cpp. Command line options to the C-preprocessor can be given",
1382 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1385 "When using position restraints a file with restraint coordinates",
1386 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1387 "with respect to the conformation from the [TT]-c[tt] option.",
1388 "For free energy calculation the the coordinates for the B topology",
1389 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1390 "those of the A topology.[PAR]",
1392 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1393 "The last frame with coordinates and velocities will be read,",
1394 "unless the [TT]-time[tt] option is used. Only if this information",
1395 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1396 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1397 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1398 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1401 "[TT]grompp[tt] can be used to restart simulations (preserving",
1402 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1403 "However, for simply changing the number of run steps to extend",
1404 "a run, using [TT]tpbconv[tt] is more convenient than [TT]grompp[tt].",
1405 "You then supply the old checkpoint file directly to [TT]mdrun[tt]",
1406 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1407 "like output frequency, then supplying the checkpoint file to",
1408 "[TT]grompp[tt] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1409 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1411 "By default, all bonded interactions which have constant energy due to",
1412 "virtual site constructions will be removed. If this constant energy is",
1413 "not zero, this will result in a shift in the total energy. All bonded",
1414 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1415 "all constraints for distances which will be constant anyway because",
1416 "of virtual site constructions will be removed. If any constraints remain",
1417 "which involve virtual sites, a fatal error will result.[PAR]"
1419 "To verify your run input file, please take note of all warnings",
1420 "on the screen, and correct where necessary. Do also look at the contents",
1421 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1422 "the input that [TT]grompp[tt] has read. If in doubt, you can start [TT]grompp[tt]",
1423 "with the [TT]-debug[tt] option which will give you more information",
1424 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1425 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1426 "program. [TT]gmxcheck[tt] can be used to compare the contents of two",
1427 "run input files.[PAR]"
1429 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1430 "by [TT]grompp[tt] that otherwise halt output. In some cases, warnings are",
1431 "harmless, but usually they are not. The user is advised to carefully",
1432 "interpret the output messages before attempting to bypass them with",
1439 gpp_atomtype_t atype;
1441 int natoms, nvsite, comb, mt;
1445 real max_spacing, fudgeQQ;
1447 char fn[STRLEN], fnB[STRLEN];
1448 const char *mdparin;
1450 gmx_bool bNeedVel, bGenVel;
1451 gmx_bool have_atomnumber;
1453 t_params *gb_plist = NULL;
1454 gmx_genborn_t *born = NULL;
1456 gmx_bool bVerbose = FALSE;
1458 char warn_buf[STRLEN];
1461 { efMDP, NULL, NULL, ffREAD },
1462 { efMDP, "-po", "mdout", ffWRITE },
1463 { efSTX, "-c", NULL, ffREAD },
1464 { efSTX, "-r", NULL, ffOPTRD },
1465 { efSTX, "-rb", NULL, ffOPTRD },
1466 { efNDX, NULL, NULL, ffOPTRD },
1467 { efTOP, NULL, NULL, ffREAD },
1468 { efTOP, "-pp", "processed", ffOPTWR },
1469 { efTPX, "-o", NULL, ffWRITE },
1470 { efTRN, "-t", NULL, ffOPTRD },
1471 { efEDR, "-e", NULL, ffOPTRD },
1472 { efTRN, "-ref", "rotref", ffOPTRW }
1474 #define NFILE asize(fnm)
1476 /* Command line options */
1477 static gmx_bool bRenum = TRUE;
1478 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1479 static int i, maxwarn = 0;
1480 static real fr_time = -1;
1482 { "-v", FALSE, etBOOL, {&bVerbose},
1483 "Be loud and noisy" },
1484 { "-time", FALSE, etREAL, {&fr_time},
1485 "Take frame at or first after this time." },
1486 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1487 "Remove constant bonded interactions with virtual sites" },
1488 { "-maxwarn", FALSE, etINT, {&maxwarn},
1489 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1490 { "-zero", FALSE, etBOOL, {&bZero},
1491 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1492 { "-renum", FALSE, etBOOL, {&bRenum},
1493 "Renumber atomtypes and minimize number of atomtypes" }
1496 CopyRight(stderr, argv[0]);
1498 /* Initiate some variables */
1503 /* Parse the command line */
1504 parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1505 asize(desc), desc, 0, NULL, &oenv);
1507 wi = init_warning(TRUE, maxwarn);
1509 /* PARAMETER file processing */
1510 mdparin = opt2fn("-f", NFILE, fnm);
1511 set_warning_line(wi, mdparin, -1);
1512 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1516 fprintf(stderr, "checking input for internal consistency...\n");
1518 check_ir(mdparin, ir, opts, wi);
1520 if (ir->ld_seed == -1)
1522 ir->ld_seed = make_seed();
1523 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1526 if (ir->expandedvals->lmc_seed == -1)
1528 ir->expandedvals->lmc_seed = make_seed();
1529 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1532 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1533 bGenVel = (bNeedVel && opts->bGenVel);
1534 if (bGenVel && ir->bContinuation)
1537 "Generating velocities is inconsistent with attempting "
1538 "to continue a previous run. Choose only one of "
1539 "gen-vel = yes and continuation = yes.");
1540 warning_error(wi, warn_buf);
1546 atype = init_atomtype();
1549 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1552 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1553 if (!gmx_fexist(fn))
1555 gmx_fatal(FARGS, "%s does not exist", fn);
1557 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1558 opts, ir, bZero, bGenVel, bVerbose, &state,
1559 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1565 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1569 /* set parameters for virtual site construction (not for vsiten) */
1570 for (mt = 0; mt < sys->nmoltype; mt++)
1573 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1575 /* now throw away all obsolete bonds, angles and dihedrals: */
1576 /* note: constraints are ALWAYS removed */
1579 for (mt = 0; mt < sys->nmoltype; mt++)
1581 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1585 if (nvsite && ir->eI == eiNM)
1587 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");
1590 if (ir->cutoff_scheme == ecutsVERLET)
1592 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1593 ecutscheme_names[ir->cutoff_scheme]);
1595 /* Remove all charge groups */
1596 gmx_mtop_remove_chargegroups(sys);
1599 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1601 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1603 sprintf(warn_buf, "Can not do %s with %s, use %s",
1604 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1605 warning_error(wi, warn_buf);
1607 if (ir->bPeriodicMols)
1609 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1610 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1611 warning_error(wi, warn_buf);
1615 if (EI_SD (ir->eI) && ir->etc != etcNO)
1617 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1620 /* If we are doing QM/MM, check that we got the atom numbers */
1621 have_atomnumber = TRUE;
1622 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1624 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1626 if (!have_atomnumber && ir->bQMMM)
1630 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1631 "field you are using does not contain atom numbers fields. This is an\n"
1632 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1633 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1634 "an integer just before the mass column in ffXXXnb.itp.\n"
1635 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1640 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1642 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1644 /** \TODO check size of ex+hy width against box size */
1647 /* Check for errors in the input now, since they might cause problems
1648 * during processing further down.
1650 check_warning_error(wi, FARGS);
1652 if (opt2bSet("-r", NFILE, fnm))
1654 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1658 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1660 if (opt2bSet("-rb", NFILE, fnm))
1662 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1669 if (nint_ftype(sys, mi, F_POSRES) > 0)
1673 fprintf(stderr, "Reading position restraint coords from %s", fn);
1674 if (strcmp(fn, fnB) == 0)
1676 fprintf(stderr, "\n");
1680 fprintf(stderr, " and %s\n", fnB);
1683 gen_posres(sys, mi, fn, fnB,
1684 ir->refcoord_scaling, ir->ePBC,
1685 ir->posres_com, ir->posres_comB,
1689 /* If we are using CMAP, setup the pre-interpolation grid */
1690 if (plist->ncmap > 0)
1692 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1693 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1696 set_wall_atomtype(atype, opts, ir, wi);
1699 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1700 ntype = get_atomtype_ntypes(atype);
1703 if (ir->implicit_solvent != eisNO)
1705 /* Now we have renumbered the atom types, we can check the GBSA params */
1706 check_gbsa_params(ir, atype);
1708 /* Check that all atoms that have charge and/or LJ-parameters also have
1709 * sensible GB-parameters
1711 check_gbsa_params_charged(sys, atype);
1714 /* PELA: Copy the atomtype data to the topology atomtype list */
1715 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1719 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1724 fprintf(stderr, "converting bonded parameters...\n");
1727 ntype = get_atomtype_ntypes(atype);
1728 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1732 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1735 /* set ptype to VSite for virtual sites */
1736 for (mt = 0; mt < sys->nmoltype; mt++)
1738 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1742 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1744 /* Check velocity for virtual sites and shells */
1747 check_vel(sys, state.v);
1750 /* check for shells and inpurecs */
1751 check_shells_inputrec(sys, ir, wi);
1756 for (i = 0; i < sys->nmoltype; i++)
1758 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1761 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1763 check_bonds_timestep(sys, ir->delta_t, wi);
1766 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1768 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.");
1771 check_warning_error(wi, FARGS);
1775 fprintf(stderr, "initialising group options...\n");
1777 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1779 bGenVel ? state.v : NULL,
1782 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1785 if (EI_DYNAMICS(ir->eI) &&
1786 !(EI_MD(ir->eI) && ir->etc == etcNO) &&
1787 inputrec2nboundeddim(ir) == 3)
1789 set_verlet_buffer(sys, ir, state.box, ir->verletbuf_drift, wi);
1793 /* Init the temperature coupling state */
1794 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1798 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1800 check_eg_vs_cg(sys);
1804 pr_symtab(debug, 0, "After index", &sys->symtab);
1806 triple_check(mdparin, ir, sys, wi);
1807 close_symtab(&sys->symtab);
1810 pr_symtab(debug, 0, "After close", &sys->symtab);
1813 /* make exclusions between QM atoms */
1816 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1818 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1822 generate_qmexcl(sys, ir, wi);
1826 if (ftp2bSet(efTRN, NFILE, fnm))
1830 fprintf(stderr, "getting data from old trajectory ...\n");
1832 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1833 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1836 if (ir->ePBC == epbcXY && ir->nwall != 2)
1838 clear_rvec(state.box[ZZ]);
1841 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1843 set_warning_line(wi, mdparin, -1);
1844 check_chargegroup_radii(sys, ir, state.x, wi);
1847 if (EEL_FULL(ir->coulombtype))
1849 /* Calculate the optimal grid dimensions */
1850 copy_mat(state.box, box);
1851 if (ir->ePBC == epbcXY && ir->nwall == 2)
1853 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1855 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1857 /* Mark fourier_spacing as not used */
1858 ir->fourier_spacing = 0;
1860 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1862 set_warning_line(wi, mdparin, -1);
1863 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1865 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1866 &(ir->nkx), &(ir->nky), &(ir->nkz));
1869 /* MRS: eventually figure out better logic for initializing the fep
1870 values that makes declaring the lambda and declaring the state not
1871 potentially conflict if not handled correctly. */
1872 if (ir->efep != efepNO)
1874 state.fep_state = ir->fepvals->init_fep_state;
1875 for (i = 0; i < efptNR; i++)
1877 /* init_lambda trumps state definitions*/
1878 if (ir->fepvals->init_lambda >= 0)
1880 state.lambda[i] = ir->fepvals->init_lambda;
1884 if (ir->fepvals->all_lambda[i] == NULL)
1886 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1890 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1896 if (ir->ePull != epullNO)
1898 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1903 set_reference_positions(ir->rot, sys, state.x, state.box,
1904 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1908 /* reset_multinr(sys); */
1910 if (EEL_PME(ir->coulombtype))
1912 float ratio = pme_load_estimate(sys, ir, state.box);
1913 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
1914 /* With free energy we might need to do PME both for the A and B state
1915 * charges. This will double the cost, but the optimal performance will
1916 * then probably be at a slightly larger cut-off and grid spacing.
1918 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1919 (ir->efep != efepNO && ratio > 2.0/3.0))
1922 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1923 "and for highly parallel simulations between 0.25 and 0.33,\n"
1924 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1925 if (ir->efep != efepNO)
1928 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1934 char warn_buf[STRLEN];
1935 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
1936 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
1939 set_warning_line(wi, mdparin, -1);
1940 warning_note(wi, warn_buf);
1944 printf("%s\n", warn_buf);
1950 fprintf(stderr, "writing run input file...\n");
1953 done_warning(wi, FARGS);
1955 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);