1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
40 #include <sys/types.h>
63 #include "sortwater.h"
65 #include "gmx_fatal.h"
71 #include "vsite_parm.h"
77 #include "compute_io.h"
78 #include "gpp_atomtype.h"
79 #include "mtop_util.h"
81 #include "calc_verletbuf.h"
85 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
90 /* For all the molecule types */
91 for (i = 0; i < nrmols; i++)
93 n += mols[i].plist[ifunc].nr;
94 mols[i].plist[ifunc].nr = 0;
99 static int check_atom_names(const char *fn1, const char *fn2,
100 gmx_mtop_t *mtop, t_atoms *at)
102 int mb, m, i, j, nmismatch;
104 #define MAXMISMATCH 20
106 if (mtop->natoms != at->nr)
108 gmx_incons("comparing atom names");
113 for (mb = 0; mb < mtop->nmolblock; mb++)
115 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
116 for (m = 0; m < mtop->molblock[mb].nmol; m++)
118 for (j = 0; j < tat->nr; j++)
120 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
122 if (nmismatch < MAXMISMATCH)
125 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
126 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
128 else if (nmismatch == MAXMISMATCH)
130 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
142 static void check_eg_vs_cg(gmx_mtop_t *mtop)
144 int astart, mb, m, cg, j, firstj;
145 unsigned char firsteg, eg;
148 /* Go through all the charge groups and make sure all their
149 * atoms are in the same energy group.
153 for (mb = 0; mb < mtop->nmolblock; mb++)
155 molt = &mtop->moltype[mtop->molblock[mb].type];
156 for (m = 0; m < mtop->molblock[mb].nmol; m++)
158 for (cg = 0; cg < molt->cgs.nr; cg++)
160 /* Get the energy group of the first atom in this charge group */
161 firstj = astart + molt->cgs.index[cg];
162 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
163 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
165 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
168 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
169 firstj+1, astart+j+1, cg+1, *molt->name);
173 astart += molt->atoms.nr;
178 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
181 char warn_buf[STRLEN];
184 for (cg = 0; cg < cgs->nr; cg++)
186 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
189 if (maxsize > MAX_CHARGEGROUP_SIZE)
191 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
193 else if (maxsize > 10)
195 set_warning_line(wi, topfn, -1);
197 "The largest charge group contains %d atoms.\n"
198 "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"
199 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
200 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
202 warning_note(wi, warn_buf);
206 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
208 /* This check is not intended to ensure accurate integration,
209 * rather it is to signal mistakes in the mdp settings.
210 * A common mistake is to forget to turn on constraints
211 * for MD after energy minimization with flexible bonds.
212 * This check can also detect too large time steps for flexible water
213 * models, but such errors will often be masked by the constraints
214 * mdp options, which turns flexible water into water with bond constraints,
215 * but without an angle constraint. Unfortunately such incorrect use
216 * of water models can not easily be detected without checking
217 * for specific model names.
219 * The stability limit of leap-frog or velocity verlet is 4.44 steps
220 * per oscillational period.
221 * But accurate bonds distributions are lost far before that limit.
222 * To allow relatively common schemes (although not common with Gromacs)
223 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
224 * we set the note limit to 10.
226 int min_steps_warn = 5;
227 int min_steps_note = 10;
230 gmx_moltype_t *moltype, *w_moltype;
232 t_ilist *ilist, *ilb, *ilc, *ils;
234 int i, a1, a2, w_a1, w_a2, j;
235 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
236 gmx_bool bFound, bWater, bWarn;
237 char warn_buf[STRLEN];
239 ip = mtop->ffparams.iparams;
241 twopi2 = sqr(2*M_PI);
243 limit2 = sqr(min_steps_note*dt);
249 for (molt = 0; molt < mtop->nmoltype; molt++)
251 moltype = &mtop->moltype[molt];
252 atom = moltype->atoms.atom;
253 ilist = moltype->ilist;
254 ilc = &ilist[F_CONSTR];
255 ils = &ilist[F_SETTLE];
256 for (ftype = 0; ftype < F_NRE; ftype++)
258 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
264 for (i = 0; i < ilb->nr; i += 3)
266 fc = ip[ilb->iatoms[i]].harmonic.krA;
267 re = ip[ilb->iatoms[i]].harmonic.rA;
268 if (ftype == F_G96BONDS)
270 /* Convert squared sqaure fc to harmonic fc */
273 a1 = ilb->iatoms[i+1];
274 a2 = ilb->iatoms[i+2];
277 if (fc > 0 && m1 > 0 && m2 > 0)
279 period2 = twopi2*m1*m2/((m1 + m2)*fc);
283 period2 = GMX_FLOAT_MAX;
287 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
288 fc, m1, m2, sqrt(period2));
290 if (period2 < limit2)
293 for (j = 0; j < ilc->nr; j += 3)
295 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
296 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
301 for (j = 0; j < ils->nr; j += 4)
303 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
304 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
310 (w_moltype == NULL || period2 < w_period2))
322 if (w_moltype != NULL)
324 bWarn = (w_period2 < sqr(min_steps_warn*dt));
325 /* A check that would recognize most water models */
326 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
327 w_moltype->atoms.nr <= 5);
328 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"
331 w_a1+1, *w_moltype->atoms.atomname[w_a1],
332 w_a2+1, *w_moltype->atoms.atomname[w_a2],
333 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
335 "Maybe you asked for fexible water." :
336 "Maybe you forgot to change the constraints mdp option.");
339 warning(wi, warn_buf);
343 warning_note(wi, warn_buf);
348 static void check_vel(gmx_mtop_t *mtop, rvec v[])
350 gmx_mtop_atomloop_all_t aloop;
354 aloop = gmx_mtop_atomloop_all_init(mtop);
355 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
357 if (atom->ptype == eptShell ||
358 atom->ptype == eptBond ||
359 atom->ptype == eptVSite)
366 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
371 for (mb = 0; mb < mtop->nmolblock; mb++)
373 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
379 /* This routine reorders the molecule type array
380 * in the order of use in the molblocks,
381 * unused molecule types are deleted.
383 static void renumber_moltypes(gmx_mtop_t *sys,
384 int *nmolinfo, t_molinfo **molinfo)
386 int *order, norder, i;
390 snew(order, *nmolinfo);
392 for (mb = 0; mb < sys->nmolblock; mb++)
394 for (i = 0; i < norder; i++)
396 if (order[i] == sys->molblock[mb].type)
403 /* This type did not occur yet, add it */
404 order[norder] = sys->molblock[mb].type;
405 /* Renumber the moltype in the topology */
408 sys->molblock[mb].type = i;
411 /* We still need to reorder the molinfo structs */
413 for (mi = 0; mi < *nmolinfo; mi++)
415 for (i = 0; i < norder; i++)
424 done_mi(&(*molinfo)[mi]);
428 minew[i] = (*molinfo)[mi];
437 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
442 mtop->nmoltype = nmi;
443 snew(mtop->moltype, nmi);
444 for (m = 0; m < nmi; m++)
446 molt = &mtop->moltype[m];
447 molt->name = mi[m].name;
448 molt->atoms = mi[m].atoms;
449 /* ilists are copied later */
450 molt->cgs = mi[m].cgs;
451 molt->excls = mi[m].excls;
456 new_status(const char *topfile, const char *topppfile, const char *confin,
457 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
458 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
459 gpp_atomtype_t atype, gmx_mtop_t *sys,
460 int *nmi, t_molinfo **mi, t_params plist[],
461 int *comb, double *reppow, real *fudgeQQ,
465 t_molinfo *molinfo = NULL;
467 gmx_molblock_t *molblock, *molbs;
469 int mb, i, nrmols, nmismatch;
471 gmx_bool bGB = FALSE;
472 char warn_buf[STRLEN];
476 /* Set gmx_boolean for GB */
477 if (ir->implicit_solvent)
482 /* TOPOLOGY processing */
483 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
484 plist, comb, reppow, fudgeQQ,
485 atype, &nrmols, &molinfo, ir,
486 &nmolblock, &molblock, bGB,
490 snew(sys->molblock, nmolblock);
493 for (mb = 0; mb < nmolblock; mb++)
495 if (sys->nmolblock > 0 &&
496 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
498 /* Merge consecutive blocks with the same molecule type */
499 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
500 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
502 else if (molblock[mb].nmol > 0)
504 /* Add a new molblock to the topology */
505 molbs = &sys->molblock[sys->nmolblock];
506 *molbs = molblock[mb];
507 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
508 molbs->nposres_xA = 0;
509 molbs->nposres_xB = 0;
510 sys->natoms += molbs->nmol*molbs->natoms_mol;
514 if (sys->nmolblock == 0)
516 gmx_fatal(FARGS, "No molecules were defined in the system");
519 renumber_moltypes(sys, &nrmols, &molinfo);
523 convert_harmonics(nrmols, molinfo, atype);
526 if (ir->eDisre == edrNone)
528 i = rm_interactions(F_DISRES, nrmols, molinfo);
531 set_warning_line(wi, "unknown", -1);
532 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
533 warning_note(wi, warn_buf);
536 if (opts->bOrire == FALSE)
538 i = rm_interactions(F_ORIRES, nrmols, molinfo);
541 set_warning_line(wi, "unknown", -1);
542 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
543 warning_note(wi, warn_buf);
547 /* Copy structures from msys to sys */
548 molinfo2mtop(nrmols, molinfo, sys);
550 gmx_mtop_finalize(sys);
552 /* COORDINATE file processing */
555 fprintf(stderr, "processing coordinates...\n");
558 get_stx_coordnum(confin, &state->natoms);
559 if (state->natoms != sys->natoms)
561 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
562 " does not match topology (%s, %d)",
563 confin, state->natoms, topfile, sys->natoms);
567 /* make space for coordinates and velocities */
570 init_t_atoms(confat, state->natoms, FALSE);
571 init_state(state, state->natoms, 0, 0, 0, 0);
572 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
573 /* This call fixes the box shape for runs with pressure scaling */
574 set_box_rel(ir, state);
576 nmismatch = check_atom_names(topfile, confin, sys, confat);
577 free_t_atoms(confat, TRUE);
582 sprintf(buf, "%d non-matching atom name%s\n"
583 "atom names from %s will be used\n"
584 "atom names from %s will be ignored\n",
585 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
590 fprintf(stderr, "double-checking input for internal consistency...\n");
592 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
598 gmx_mtop_atomloop_all_t aloop;
601 snew(mass, state->natoms);
602 aloop = gmx_mtop_atomloop_all_init(sys);
603 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
608 if (opts->seed == -1)
610 opts->seed = make_seed();
611 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
613 maxwell_speed(opts->tempi, opts->seed, sys, state->v);
615 stop_cm(stdout, state->natoms, mass, state->x, state->v);
623 static void copy_state(const char *slog, t_trxframe *fr,
624 gmx_bool bReadVel, t_state *state,
629 if (fr->not_ok & FRAME_NOT_OK)
631 gmx_fatal(FARGS, "Can not start from an incomplete frame");
635 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
639 for (i = 0; i < state->natoms; i++)
641 copy_rvec(fr->x[i], state->x[i]);
647 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
649 for (i = 0; i < state->natoms; i++)
651 copy_rvec(fr->v[i], state->v[i]);
656 copy_mat(fr->box, state->box);
659 *use_time = fr->time;
662 static void cont_status(const char *slog, const char *ener,
663 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
664 t_inputrec *ir, t_state *state,
666 const output_env_t oenv)
667 /* If fr_time == -1 read the last frame available which is complete */
675 bReadVel = (bNeedVel && !bGenVel);
678 "Reading Coordinates%s and Box size from old trajectory\n",
679 bReadVel ? ", Velocities" : "");
682 fprintf(stderr, "Will read whole trajectory\n");
686 fprintf(stderr, "Will read till time %g\n", fr_time);
692 fprintf(stderr, "Velocities generated: "
693 "ignoring velocities in input trajectory\n");
695 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
699 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
705 "WARNING: Did not find a frame with velocities in file %s,\n"
706 " all velocities will be set to zero!\n\n", slog);
707 for (i = 0; i < sys->natoms; i++)
709 clear_rvec(state->v[i]);
712 /* Search for a frame without velocities */
714 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
718 state->natoms = fr.natoms;
720 if (sys->natoms != state->natoms)
722 gmx_fatal(FARGS, "Number of atoms in Topology "
723 "is not the same as in Trajectory");
725 copy_state(slog, &fr, bReadVel, state, &use_time);
727 /* Find the appropriate frame */
728 while ((fr_time == -1 || fr.time < fr_time) &&
729 read_next_frame(oenv, fp, &fr))
731 copy_state(slog, &fr, bReadVel, state, &use_time);
736 /* Set the relative box lengths for preserving the box shape.
737 * Note that this call can lead to differences in the last bit
738 * with respect to using tpbconv to create a [TT].tpx[tt] file.
740 set_box_rel(ir, state);
742 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
743 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
745 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
747 get_enx_state(ener, use_time, &sys->groups, ir, state);
748 preserve_box_shape(ir, state->box_rel, state->boxv);
752 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
754 int rc_scaling, int ePBC,
758 gmx_bool bFirst = TRUE, *hadAtom;
764 int natoms, npbcdim = 0;
765 char warn_buf[STRLEN], title[STRLEN];
766 int a, i, ai, j, k, mb, nat_molb;
767 gmx_molblock_t *molb;
771 get_stx_coordnum(fn, &natoms);
772 if (natoms != mtop->natoms)
774 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);
775 warning(wi, warn_buf);
779 init_t_atoms(&dumat, natoms, FALSE);
780 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
782 npbcdim = ePBC2npbcdim(ePBC);
784 if (rc_scaling != erscNO)
786 copy_mat(box, invbox);
787 for (j = npbcdim; j < DIM; j++)
789 clear_rvec(invbox[j]);
792 m_inv_ur0(invbox, invbox);
795 /* Copy the reference coordinates to mtop */
799 snew(hadAtom, natoms);
800 for (mb = 0; mb < mtop->nmolblock; mb++)
802 molb = &mtop->molblock[mb];
803 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
804 pr = &(molinfo[molb->type].plist[F_POSRES]);
805 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
806 if (pr->nr > 0 || prfb->nr > 0)
808 atom = mtop->moltype[molb->type].atoms.atom;
809 for (i = 0; (i < pr->nr); i++)
811 ai = pr->param[i].AI;
814 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
815 ai+1, *molinfo[molb->type].name, fn, natoms);
818 if (rc_scaling == erscCOM)
820 /* Determine the center of mass of the posres reference coordinates */
821 for (j = 0; j < npbcdim; j++)
823 sum[j] += atom[ai].m*x[a+ai][j];
825 totmass += atom[ai].m;
828 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
829 for (i = 0; (i < prfb->nr); i++)
831 ai = prfb->param[i].AI;
834 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
835 ai+1, *molinfo[molb->type].name, fn, natoms);
837 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
839 /* Determine the center of mass of the posres reference coordinates */
840 for (j = 0; j < npbcdim; j++)
842 sum[j] += atom[ai].m*x[a+ai][j];
844 totmass += atom[ai].m;
849 molb->nposres_xA = nat_molb;
850 snew(molb->posres_xA, molb->nposres_xA);
851 for (i = 0; i < nat_molb; i++)
853 copy_rvec(x[a+i], molb->posres_xA[i]);
858 molb->nposres_xB = nat_molb;
859 snew(molb->posres_xB, molb->nposres_xB);
860 for (i = 0; i < nat_molb; i++)
862 copy_rvec(x[a+i], molb->posres_xB[i]);
868 if (rc_scaling == erscCOM)
872 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
874 for (j = 0; j < npbcdim; j++)
876 com[j] = sum[j]/totmass;
878 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]);
881 if (rc_scaling != erscNO)
883 for (mb = 0; mb < mtop->nmolblock; mb++)
885 molb = &mtop->molblock[mb];
886 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
887 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
889 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
890 for (i = 0; i < nat_molb; i++)
892 for (j = 0; j < npbcdim; j++)
894 if (rc_scaling == erscALL)
896 /* Convert from Cartesian to crystal coordinates */
897 xp[i][j] *= invbox[j][j];
898 for (k = j+1; k < npbcdim; k++)
900 xp[i][j] += invbox[k][j]*xp[i][k];
903 else if (rc_scaling == erscCOM)
905 /* Subtract the center of mass */
913 if (rc_scaling == erscCOM)
915 /* Convert the COM from Cartesian to crystal coordinates */
916 for (j = 0; j < npbcdim; j++)
918 com[j] *= invbox[j][j];
919 for (k = j+1; k < npbcdim; k++)
921 com[j] += invbox[k][j]*com[k];
927 free_t_atoms(&dumat, TRUE);
933 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
934 char *fnA, char *fnB,
935 int rc_scaling, int ePBC,
941 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
942 if (strcmp(fnA, fnB) != 0)
944 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
948 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
949 t_inputrec *ir, warninp_t wi)
952 char warn_buf[STRLEN];
956 fprintf(stderr, "Searching the wall atom type(s)\n");
958 for (i = 0; i < ir->nwall; i++)
960 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
961 if (ir->wall_atomtype[i] == NOTSET)
963 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
964 warning_error(wi, warn_buf);
969 static int nrdf_internal(t_atoms *atoms)
974 for (i = 0; i < atoms->nr; i++)
976 /* Vsite ptype might not be set here yet, so also check the mass */
977 if ((atoms->atom[i].ptype == eptAtom ||
978 atoms->atom[i].ptype == eptNucleus)
979 && atoms->atom[i].m > 0)
986 case 0: nrdf = 0; break;
987 case 1: nrdf = 0; break;
988 case 2: nrdf = 1; break;
989 default: nrdf = nmass*3 - 6; break;
1008 for (i = 1; i < n-1; i++)
1010 p = 0.5*y2[i-1]+2.0;
1012 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1013 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1018 for (i = n-2; i >= 0; i--)
1020 y2[i] = y2[i]*y2[i+1]+u[i];
1026 interpolate1d( double xmin,
1039 a = (xmin+(ix+1)*dx-x)/dx;
1040 b = (x-xmin-ix*dx)/dx;
1042 *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;
1043 *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];
1048 setup_cmap (int grid_spacing,
1051 gmx_cmap_t * cmap_grid)
1053 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1055 int i, j, k, ii, jj, kk, idx;
1057 double dx, xmin, v, v1, v2, v12;
1060 snew(tmp_u, 2*grid_spacing);
1061 snew(tmp_u2, 2*grid_spacing);
1062 snew(tmp_yy, 2*grid_spacing);
1063 snew(tmp_y1, 2*grid_spacing);
1064 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1065 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1067 dx = 360.0/grid_spacing;
1068 xmin = -180.0-dx*grid_spacing/2;
1070 for (kk = 0; kk < nc; kk++)
1072 /* Compute an offset depending on which cmap we are using
1073 * Offset will be the map number multiplied with the
1074 * grid_spacing * grid_spacing * 2
1076 offset = kk * grid_spacing * grid_spacing * 2;
1078 for (i = 0; i < 2*grid_spacing; i++)
1080 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1082 for (j = 0; j < 2*grid_spacing; j++)
1084 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1085 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1089 for (i = 0; i < 2*grid_spacing; i++)
1091 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1094 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1096 ii = i-grid_spacing/2;
1099 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1101 jj = j-grid_spacing/2;
1104 for (k = 0; k < 2*grid_spacing; k++)
1106 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1107 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1110 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1111 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1112 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1113 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1115 idx = ii*grid_spacing+jj;
1116 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1117 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1118 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1119 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1125 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1129 cmap_grid->ngrid = ngrid;
1130 cmap_grid->grid_spacing = grid_spacing;
1131 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1133 snew(cmap_grid->cmapdata, ngrid);
1135 for (i = 0; i < cmap_grid->ngrid; i++)
1137 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1142 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1144 int count, count_mol, i, mb;
1145 gmx_molblock_t *molb;
1150 for (mb = 0; mb < mtop->nmolblock; mb++)
1153 molb = &mtop->molblock[mb];
1154 plist = mi[molb->type].plist;
1156 for (i = 0; i < F_NRE; i++)
1160 count_mol += 3*plist[i].nr;
1162 else if (interaction_function[i].flags & IF_CONSTRAINT)
1164 count_mol += plist[i].nr;
1168 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1171 "Molecule type '%s' has %d constraints.\n"
1172 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1173 *mi[molb->type].name, count_mol,
1174 nrdf_internal(&mi[molb->type].atoms));
1177 count += molb->nmol*count_mol;
1183 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1185 int i, nmiss, natoms, mt;
1187 const t_atoms *atoms;
1190 for (mt = 0; mt < sys->nmoltype; mt++)
1192 atoms = &sys->moltype[mt].atoms;
1195 for (i = 0; i < natoms; i++)
1197 q = atoms->atom[i].q;
1198 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1199 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1200 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1201 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1202 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1205 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1206 get_atomtype_name(atoms->atom[i].type, atype), q);
1214 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);
1219 static void check_gbsa_params(gpp_atomtype_t atype)
1223 /* If we are doing GBSA, check that we got the parameters we need
1224 * This checking is to see if there are GBSA paratmeters for all
1225 * atoms in the force field. To go around this for testing purposes
1226 * comment out the nerror++ counter temporarily
1229 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1231 if (get_atomtype_radius(i, atype) < 0 ||
1232 get_atomtype_vol(i, atype) < 0 ||
1233 get_atomtype_surftens(i, atype) < 0 ||
1234 get_atomtype_gb_radius(i, atype) < 0 ||
1235 get_atomtype_S_hct(i, atype) < 0)
1237 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1238 get_atomtype_name(i, atype));
1245 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);
1250 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1253 real verletbuf_drift,
1258 verletbuf_list_setup_t ls;
1261 char warn_buf[STRLEN];
1264 for (i = 0; i < ir->opts.ngtc; i++)
1266 if (ir->opts.ref_t[i] < 0)
1268 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.");
1272 ref_T = max(ref_T, ir->opts.ref_t[i]);
1276 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n", verletbuf_drift, ref_T);
1278 for (i = 0; i < ir->opts.ngtc; i++)
1280 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1282 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.",
1283 ir->opts.nrdf[i], ir->opts.ref_t[i], ref_T);
1284 warning_note(wi, warn_buf);
1288 /* Calculate the buffer size for simple atom vs atoms list */
1289 ls.cluster_size_i = 1;
1290 ls.cluster_size_j = 1;
1291 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1292 &ls, &n_nonlin_vsite, &rlist_1x1);
1294 /* Set the pair-list buffer size in ir */
1295 verletbuf_get_list_setup(FALSE, &ls);
1296 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1297 &ls, &n_nonlin_vsite, &ir->rlist);
1299 if (n_nonlin_vsite > 0)
1301 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);
1302 warning_note(wi, warn_buf);
1305 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1306 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1308 ir->rlistlong = ir->rlist;
1309 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1310 ls.cluster_size_i, ls.cluster_size_j,
1311 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1313 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1315 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)));
1319 int gmx_grompp(int argc, char *argv[])
1321 static const char *desc[] = {
1322 "The gromacs preprocessor",
1323 "reads a molecular topology file, checks the validity of the",
1324 "file, expands the topology from a molecular description to an atomic",
1325 "description. The topology file contains information about",
1326 "molecule types and the number of molecules, the preprocessor",
1327 "copies each molecule as needed. ",
1328 "There is no limitation on the number of molecule types. ",
1329 "Bonds and bond-angles can be converted into constraints, separately",
1330 "for hydrogens and heavy atoms.",
1331 "Then a coordinate file is read and velocities can be generated",
1332 "from a Maxwellian distribution if requested.",
1333 "[TT]grompp[tt] also reads parameters for the [TT]mdrun[tt] ",
1334 "(eg. number of MD steps, time step, cut-off), and others such as",
1335 "NEMD parameters, which are corrected so that the net acceleration",
1337 "Eventually a binary file is produced that can serve as the sole input",
1338 "file for the MD program.[PAR]",
1340 "[TT]grompp[tt] uses the atom names from the topology file. The atom names",
1341 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1342 "warnings when they do not match the atom names in the topology.",
1343 "Note that the atom names are irrelevant for the simulation as",
1344 "only the atom types are used for generating interaction parameters.[PAR]",
1346 "[TT]grompp[tt] uses a built-in preprocessor to resolve includes, macros, ",
1347 "etc. The preprocessor supports the following keywords:[PAR]",
1348 "#ifdef VARIABLE[BR]",
1349 "#ifndef VARIABLE[BR]",
1352 "#define VARIABLE[BR]",
1353 "#undef VARIABLE[BR]"
1354 "#include \"filename\"[BR]",
1355 "#include <filename>[PAR]",
1356 "The functioning of these statements in your topology may be modulated by",
1357 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1358 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1359 "include = -I/home/john/doe[tt][BR]",
1360 "For further information a C-programming textbook may help you out.",
1361 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1362 "topology file written out so that you can verify its contents.[PAR]",
1364 /* cpp has been unnecessary for some time, hasn't it?
1365 "If your system does not have a C-preprocessor, you can still",
1366 "use [TT]grompp[tt], but you do not have access to the features ",
1367 "from the cpp. Command line options to the C-preprocessor can be given",
1368 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1371 "When using position restraints a file with restraint coordinates",
1372 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1373 "with respect to the conformation from the [TT]-c[tt] option.",
1374 "For free energy calculation the the coordinates for the B topology",
1375 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1376 "those of the A topology.[PAR]",
1378 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1379 "The last frame with coordinates and velocities will be read,",
1380 "unless the [TT]-time[tt] option is used. Only if this information",
1381 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1382 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1383 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1384 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1387 "[TT]grompp[tt] can be used to restart simulations (preserving",
1388 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1389 "However, for simply changing the number of run steps to extend",
1390 "a run, using [TT]tpbconv[tt] is more convenient than [TT]grompp[tt].",
1391 "You then supply the old checkpoint file directly to [TT]mdrun[tt]",
1392 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1393 "like output frequency, then supplying the checkpoint file to",
1394 "[TT]grompp[tt] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1395 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1397 "By default, all bonded interactions which have constant energy due to",
1398 "virtual site constructions will be removed. If this constant energy is",
1399 "not zero, this will result in a shift in the total energy. All bonded",
1400 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1401 "all constraints for distances which will be constant anyway because",
1402 "of virtual site constructions will be removed. If any constraints remain",
1403 "which involve virtual sites, a fatal error will result.[PAR]"
1405 "To verify your run input file, please take note of all warnings",
1406 "on the screen, and correct where necessary. Do also look at the contents",
1407 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1408 "the input that [TT]grompp[tt] has read. If in doubt, you can start [TT]grompp[tt]",
1409 "with the [TT]-debug[tt] option which will give you more information",
1410 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1411 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1412 "program. [TT]gmxcheck[tt] can be used to compare the contents of two",
1413 "run input files.[PAR]"
1415 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1416 "by [TT]grompp[tt] that otherwise halt output. In some cases, warnings are",
1417 "harmless, but usually they are not. The user is advised to carefully",
1418 "interpret the output messages before attempting to bypass them with",
1425 gpp_atomtype_t atype;
1427 int natoms, nvsite, comb, mt;
1431 real max_spacing, fudgeQQ;
1433 char fn[STRLEN], fnB[STRLEN];
1434 const char *mdparin;
1436 gmx_bool bNeedVel, bGenVel;
1437 gmx_bool have_atomnumber;
1439 t_params *gb_plist = NULL;
1440 gmx_genborn_t *born = NULL;
1442 gmx_bool bVerbose = FALSE;
1444 char warn_buf[STRLEN];
1447 { efMDP, NULL, NULL, ffREAD },
1448 { efMDP, "-po", "mdout", ffWRITE },
1449 { efSTX, "-c", NULL, ffREAD },
1450 { efSTX, "-r", NULL, ffOPTRD },
1451 { efSTX, "-rb", NULL, ffOPTRD },
1452 { efNDX, NULL, NULL, ffOPTRD },
1453 { efTOP, NULL, NULL, ffREAD },
1454 { efTOP, "-pp", "processed", ffOPTWR },
1455 { efTPX, "-o", NULL, ffWRITE },
1456 { efTRN, "-t", NULL, ffOPTRD },
1457 { efEDR, "-e", NULL, ffOPTRD },
1458 { efTRN, "-ref", "rotref", ffOPTRW }
1460 #define NFILE asize(fnm)
1462 /* Command line options */
1463 static gmx_bool bRenum = TRUE;
1464 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1465 static int i, maxwarn = 0;
1466 static real fr_time = -1;
1468 { "-v", FALSE, etBOOL, {&bVerbose},
1469 "Be loud and noisy" },
1470 { "-time", FALSE, etREAL, {&fr_time},
1471 "Take frame at or first after this time." },
1472 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1473 "Remove constant bonded interactions with virtual sites" },
1474 { "-maxwarn", FALSE, etINT, {&maxwarn},
1475 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1476 { "-zero", FALSE, etBOOL, {&bZero},
1477 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1478 { "-renum", FALSE, etBOOL, {&bRenum},
1479 "Renumber atomtypes and minimize number of atomtypes" }
1482 /* Initiate some variables */
1487 /* Parse the command line */
1488 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1489 asize(desc), desc, 0, NULL, &oenv))
1494 wi = init_warning(TRUE, maxwarn);
1496 /* PARAMETER file processing */
1497 mdparin = opt2fn("-f", NFILE, fnm);
1498 set_warning_line(wi, mdparin, -1);
1499 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1503 fprintf(stderr, "checking input for internal consistency...\n");
1505 check_ir(mdparin, ir, opts, wi);
1507 if (ir->ld_seed == -1)
1509 ir->ld_seed = make_seed();
1510 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1513 if (ir->expandedvals->lmc_seed == -1)
1515 ir->expandedvals->lmc_seed = make_seed();
1516 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1519 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1520 bGenVel = (bNeedVel && opts->bGenVel);
1521 if (bGenVel && ir->bContinuation)
1524 "Generating velocities is inconsistent with attempting "
1525 "to continue a previous run. Choose only one of "
1526 "gen-vel = yes and continuation = yes.");
1527 warning_error(wi, warn_buf);
1533 atype = init_atomtype();
1536 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1539 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1540 if (!gmx_fexist(fn))
1542 gmx_fatal(FARGS, "%s does not exist", fn);
1544 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1545 opts, ir, bZero, bGenVel, bVerbose, &state,
1546 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1552 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1556 /* set parameters for virtual site construction (not for vsiten) */
1557 for (mt = 0; mt < sys->nmoltype; mt++)
1560 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1562 /* now throw away all obsolete bonds, angles and dihedrals: */
1563 /* note: constraints are ALWAYS removed */
1566 for (mt = 0; mt < sys->nmoltype; mt++)
1568 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1572 if (ir->cutoff_scheme == ecutsVERLET)
1574 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1575 ecutscheme_names[ir->cutoff_scheme]);
1577 /* Remove all charge groups */
1578 gmx_mtop_remove_chargegroups(sys);
1581 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1583 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1585 sprintf(warn_buf, "Can not do %s with %s, use %s",
1586 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1587 warning_error(wi, warn_buf);
1589 if (ir->bPeriodicMols)
1591 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1592 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1593 warning_error(wi, warn_buf);
1597 if (EI_SD (ir->eI) && ir->etc != etcNO)
1599 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1602 /* If we are doing QM/MM, check that we got the atom numbers */
1603 have_atomnumber = TRUE;
1604 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1606 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1608 if (!have_atomnumber && ir->bQMMM)
1612 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1613 "field you are using does not contain atom numbers fields. This is an\n"
1614 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1615 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1616 "an integer just before the mass column in ffXXXnb.itp.\n"
1617 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1622 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1624 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1626 /** \TODO check size of ex+hy width against box size */
1629 /* Check for errors in the input now, since they might cause problems
1630 * during processing further down.
1632 check_warning_error(wi, FARGS);
1634 if (opt2bSet("-r", NFILE, fnm))
1636 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1640 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1642 if (opt2bSet("-rb", NFILE, fnm))
1644 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1651 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1655 fprintf(stderr, "Reading position restraint coords from %s", fn);
1656 if (strcmp(fn, fnB) == 0)
1658 fprintf(stderr, "\n");
1662 fprintf(stderr, " and %s\n", fnB);
1665 gen_posres(sys, mi, fn, fnB,
1666 ir->refcoord_scaling, ir->ePBC,
1667 ir->posres_com, ir->posres_comB,
1671 /* If we are using CMAP, setup the pre-interpolation grid */
1672 if (plist->ncmap > 0)
1674 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1675 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1678 set_wall_atomtype(atype, opts, ir, wi);
1681 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1682 ntype = get_atomtype_ntypes(atype);
1685 if (ir->implicit_solvent != eisNO)
1687 /* Now we have renumbered the atom types, we can check the GBSA params */
1688 check_gbsa_params(atype);
1690 /* Check that all atoms that have charge and/or LJ-parameters also have
1691 * sensible GB-parameters
1693 check_gbsa_params_charged(sys, atype);
1696 /* PELA: Copy the atomtype data to the topology atomtype list */
1697 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1701 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1706 fprintf(stderr, "converting bonded parameters...\n");
1709 ntype = get_atomtype_ntypes(atype);
1710 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1714 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1717 /* set ptype to VSite for virtual sites */
1718 for (mt = 0; mt < sys->nmoltype; mt++)
1720 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1724 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1726 /* Check velocity for virtual sites and shells */
1729 check_vel(sys, state.v);
1735 for (i = 0; i < sys->nmoltype; i++)
1737 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1740 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1742 check_bonds_timestep(sys, ir->delta_t, wi);
1745 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1747 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.");
1750 check_warning_error(wi, FARGS);
1754 fprintf(stderr, "initialising group options...\n");
1756 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1758 bGenVel ? state.v : NULL,
1761 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1764 if (EI_DYNAMICS(ir->eI) &&
1765 !(EI_MD(ir->eI) && ir->etc == etcNO) &&
1766 inputrec2nboundeddim(ir) == 3)
1768 set_verlet_buffer(sys, ir, state.box, ir->verletbuf_drift, wi);
1772 /* Init the temperature coupling state */
1773 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1777 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1779 check_eg_vs_cg(sys);
1783 pr_symtab(debug, 0, "After index", &sys->symtab);
1785 triple_check(mdparin, ir, sys, wi);
1786 close_symtab(&sys->symtab);
1789 pr_symtab(debug, 0, "After close", &sys->symtab);
1792 /* make exclusions between QM atoms */
1795 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1797 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1801 generate_qmexcl(sys, ir, wi);
1805 if (ftp2bSet(efTRN, NFILE, fnm))
1809 fprintf(stderr, "getting data from old trajectory ...\n");
1811 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1812 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1815 if (ir->ePBC == epbcXY && ir->nwall != 2)
1817 clear_rvec(state.box[ZZ]);
1820 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1822 set_warning_line(wi, mdparin, -1);
1823 check_chargegroup_radii(sys, ir, state.x, wi);
1826 if (EEL_FULL(ir->coulombtype))
1828 /* Calculate the optimal grid dimensions */
1829 copy_mat(state.box, box);
1830 if (ir->ePBC == epbcXY && ir->nwall == 2)
1832 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1834 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1836 /* Mark fourier_spacing as not used */
1837 ir->fourier_spacing = 0;
1839 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1841 set_warning_line(wi, mdparin, -1);
1842 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1844 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1845 &(ir->nkx), &(ir->nky), &(ir->nkz));
1848 /* MRS: eventually figure out better logic for initializing the fep
1849 values that makes declaring the lambda and declaring the state not
1850 potentially conflict if not handled correctly. */
1851 if (ir->efep != efepNO)
1853 state.fep_state = ir->fepvals->init_fep_state;
1854 for (i = 0; i < efptNR; i++)
1856 /* init_lambda trumps state definitions*/
1857 if (ir->fepvals->init_lambda >= 0)
1859 state.lambda[i] = ir->fepvals->init_lambda;
1863 if (ir->fepvals->all_lambda[i] == NULL)
1865 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1869 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1875 if (ir->ePull != epullNO)
1877 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1882 set_reference_positions(ir->rot, state.x, state.box,
1883 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1887 /* reset_multinr(sys); */
1889 if (EEL_PME(ir->coulombtype))
1891 float ratio = pme_load_estimate(sys, ir, state.box);
1892 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
1893 /* With free energy we might need to do PME both for the A and B state
1894 * charges. This will double the cost, but the optimal performance will
1895 * then probably be at a slightly larger cut-off and grid spacing.
1897 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1898 (ir->efep != efepNO && ratio > 2.0/3.0))
1901 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1902 "and for highly parallel simulations between 0.25 and 0.33,\n"
1903 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1904 if (ir->efep != efepNO)
1907 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1913 char warn_buf[STRLEN];
1914 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
1915 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
1918 set_warning_line(wi, mdparin, -1);
1919 warning_note(wi, warn_buf);
1923 printf("%s\n", warn_buf);
1929 fprintf(stderr, "writing run input file...\n");
1932 done_warning(wi, FARGS);
1934 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);