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>
53 #include "gromacs/fileio/confio.h"
60 #include "gromacs/fileio/futil.h"
63 #include "sortwater.h"
65 #include "gmx_fatal.h"
68 #include "gromacs/fileio/gmxfio.h"
69 #include "gromacs/fileio/trnio.h"
70 #include "gromacs/fileio/tpxio.h"
71 #include "gromacs/fileio/trxio.h"
72 #include "vsite_parm.h"
76 #include "gromacs/fileio/enxio.h"
78 #include "compute_io.h"
79 #include "gpp_atomtype.h"
80 #include "mtop_util.h"
82 #include "calc_verletbuf.h"
86 static int rm_interactions(int ifunc, int nrmols, t_molinfo mols[])
91 /* For all the molecule types */
92 for (i = 0; i < nrmols; i++)
94 n += mols[i].plist[ifunc].nr;
95 mols[i].plist[ifunc].nr = 0;
100 static int check_atom_names(const char *fn1, const char *fn2,
101 gmx_mtop_t *mtop, t_atoms *at)
103 int mb, m, i, j, nmismatch;
105 #define MAXMISMATCH 20
107 if (mtop->natoms != at->nr)
109 gmx_incons("comparing atom names");
114 for (mb = 0; mb < mtop->nmolblock; mb++)
116 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
117 for (m = 0; m < mtop->molblock[mb].nmol; m++)
119 for (j = 0; j < tat->nr; j++)
121 if (strcmp( *(tat->atomname[j]), *(at->atomname[i]) ) != 0)
123 if (nmismatch < MAXMISMATCH)
126 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
127 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
129 else if (nmismatch == MAXMISMATCH)
131 fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
143 static void check_eg_vs_cg(gmx_mtop_t *mtop)
145 int astart, mb, m, cg, j, firstj;
146 unsigned char firsteg, eg;
149 /* Go through all the charge groups and make sure all their
150 * atoms are in the same energy group.
154 for (mb = 0; mb < mtop->nmolblock; mb++)
156 molt = &mtop->moltype[mtop->molblock[mb].type];
157 for (m = 0; m < mtop->molblock[mb].nmol; m++)
159 for (cg = 0; cg < molt->cgs.nr; cg++)
161 /* Get the energy group of the first atom in this charge group */
162 firstj = astart + molt->cgs.index[cg];
163 firsteg = ggrpnr(&mtop->groups, egcENER, firstj);
164 for (j = molt->cgs.index[cg]+1; j < molt->cgs.index[cg+1]; j++)
166 eg = ggrpnr(&mtop->groups, egcENER, astart+j);
169 gmx_fatal(FARGS, "atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
170 firstj+1, astart+j+1, cg+1, *molt->name);
174 astart += molt->atoms.nr;
179 static void check_cg_sizes(const char *topfn, t_block *cgs, warninp_t wi)
182 char warn_buf[STRLEN];
185 for (cg = 0; cg < cgs->nr; cg++)
187 maxsize = max(maxsize, cgs->index[cg+1]-cgs->index[cg]);
190 if (maxsize > MAX_CHARGEGROUP_SIZE)
192 gmx_fatal(FARGS, "The largest charge group contains %d atoms. The maximum is %d.", maxsize, MAX_CHARGEGROUP_SIZE);
194 else if (maxsize > 10)
196 set_warning_line(wi, topfn, -1);
198 "The largest charge group contains %d atoms.\n"
199 "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"
200 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
201 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
203 warning_note(wi, warn_buf);
207 static void check_bonds_timestep(gmx_mtop_t *mtop, double dt, warninp_t wi)
209 /* This check is not intended to ensure accurate integration,
210 * rather it is to signal mistakes in the mdp settings.
211 * A common mistake is to forget to turn on constraints
212 * for MD after energy minimization with flexible bonds.
213 * This check can also detect too large time steps for flexible water
214 * models, but such errors will often be masked by the constraints
215 * mdp options, which turns flexible water into water with bond constraints,
216 * but without an angle constraint. Unfortunately such incorrect use
217 * of water models can not easily be detected without checking
218 * for specific model names.
220 * The stability limit of leap-frog or velocity verlet is 4.44 steps
221 * per oscillational period.
222 * But accurate bonds distributions are lost far before that limit.
223 * To allow relatively common schemes (although not common with Gromacs)
224 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
225 * we set the note limit to 10.
227 int min_steps_warn = 5;
228 int min_steps_note = 10;
231 gmx_moltype_t *moltype, *w_moltype;
233 t_ilist *ilist, *ilb, *ilc, *ils;
235 int i, a1, a2, w_a1, w_a2, j;
236 real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
237 gmx_bool bFound, bWater, bWarn;
238 char warn_buf[STRLEN];
240 ip = mtop->ffparams.iparams;
242 twopi2 = sqr(2*M_PI);
244 limit2 = sqr(min_steps_note*dt);
250 for (molt = 0; molt < mtop->nmoltype; molt++)
252 moltype = &mtop->moltype[molt];
253 atom = moltype->atoms.atom;
254 ilist = moltype->ilist;
255 ilc = &ilist[F_CONSTR];
256 ils = &ilist[F_SETTLE];
257 for (ftype = 0; ftype < F_NRE; ftype++)
259 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
265 for (i = 0; i < ilb->nr; i += 3)
267 fc = ip[ilb->iatoms[i]].harmonic.krA;
268 re = ip[ilb->iatoms[i]].harmonic.rA;
269 if (ftype == F_G96BONDS)
271 /* Convert squared sqaure fc to harmonic fc */
274 a1 = ilb->iatoms[i+1];
275 a2 = ilb->iatoms[i+2];
278 if (fc > 0 && m1 > 0 && m2 > 0)
280 period2 = twopi2*m1*m2/((m1 + m2)*fc);
284 period2 = GMX_FLOAT_MAX;
288 fprintf(debug, "fc %g m1 %g m2 %g period %g\n",
289 fc, m1, m2, sqrt(period2));
291 if (period2 < limit2)
294 for (j = 0; j < ilc->nr; j += 3)
296 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
297 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
302 for (j = 0; j < ils->nr; j += 4)
304 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
305 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
311 (w_moltype == NULL || period2 < w_period2))
323 if (w_moltype != NULL)
325 bWarn = (w_period2 < sqr(min_steps_warn*dt));
326 /* A check that would recognize most water models */
327 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
328 w_moltype->atoms.nr <= 5);
329 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"
332 w_a1+1, *w_moltype->atoms.atomname[w_a1],
333 w_a2+1, *w_moltype->atoms.atomname[w_a2],
334 sqrt(w_period2), bWarn ? min_steps_warn : min_steps_note, dt,
336 "Maybe you asked for fexible water." :
337 "Maybe you forgot to change the constraints mdp option.");
340 warning(wi, warn_buf);
344 warning_note(wi, warn_buf);
349 static void check_vel(gmx_mtop_t *mtop, rvec v[])
351 gmx_mtop_atomloop_all_t aloop;
355 aloop = gmx_mtop_atomloop_all_init(mtop);
356 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
358 if (atom->ptype == eptShell ||
359 atom->ptype == eptBond ||
360 atom->ptype == eptVSite)
367 static gmx_bool nint_ftype(gmx_mtop_t *mtop, t_molinfo *mi, int ftype)
372 for (mb = 0; mb < mtop->nmolblock; mb++)
374 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
380 /* This routine reorders the molecule type array
381 * in the order of use in the molblocks,
382 * unused molecule types are deleted.
384 static void renumber_moltypes(gmx_mtop_t *sys,
385 int *nmolinfo, t_molinfo **molinfo)
387 int *order, norder, i;
391 snew(order, *nmolinfo);
393 for (mb = 0; mb < sys->nmolblock; mb++)
395 for (i = 0; i < norder; i++)
397 if (order[i] == sys->molblock[mb].type)
404 /* This type did not occur yet, add it */
405 order[norder] = sys->molblock[mb].type;
406 /* Renumber the moltype in the topology */
409 sys->molblock[mb].type = i;
412 /* We still need to reorder the molinfo structs */
414 for (mi = 0; mi < *nmolinfo; mi++)
416 for (i = 0; i < norder; i++)
425 done_mi(&(*molinfo)[mi]);
429 minew[i] = (*molinfo)[mi];
438 static void molinfo2mtop(int nmi, t_molinfo *mi, gmx_mtop_t *mtop)
443 mtop->nmoltype = nmi;
444 snew(mtop->moltype, nmi);
445 for (m = 0; m < nmi; m++)
447 molt = &mtop->moltype[m];
448 molt->name = mi[m].name;
449 molt->atoms = mi[m].atoms;
450 /* ilists are copied later */
451 molt->cgs = mi[m].cgs;
452 molt->excls = mi[m].excls;
457 new_status(const char *topfile, const char *topppfile, const char *confin,
458 t_gromppopts *opts, t_inputrec *ir, gmx_bool bZero,
459 gmx_bool bGenVel, gmx_bool bVerbose, t_state *state,
460 gpp_atomtype_t atype, gmx_mtop_t *sys,
461 int *nmi, t_molinfo **mi, t_params plist[],
462 int *comb, double *reppow, real *fudgeQQ,
466 t_molinfo *molinfo = NULL;
468 gmx_molblock_t *molblock, *molbs;
470 int mb, i, nrmols, nmismatch;
472 gmx_bool bGB = FALSE;
473 char warn_buf[STRLEN];
477 /* Set gmx_boolean for GB */
478 if (ir->implicit_solvent)
483 /* TOPOLOGY processing */
484 sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
485 plist, comb, reppow, fudgeQQ,
486 atype, &nrmols, &molinfo, ir,
487 &nmolblock, &molblock, bGB,
491 snew(sys->molblock, nmolblock);
494 for (mb = 0; mb < nmolblock; mb++)
496 if (sys->nmolblock > 0 &&
497 molblock[mb].type == sys->molblock[sys->nmolblock-1].type)
499 /* Merge consecutive blocks with the same molecule type */
500 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
501 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
503 else if (molblock[mb].nmol > 0)
505 /* Add a new molblock to the topology */
506 molbs = &sys->molblock[sys->nmolblock];
507 *molbs = molblock[mb];
508 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
509 molbs->nposres_xA = 0;
510 molbs->nposres_xB = 0;
511 sys->natoms += molbs->nmol*molbs->natoms_mol;
515 if (sys->nmolblock == 0)
517 gmx_fatal(FARGS, "No molecules were defined in the system");
520 renumber_moltypes(sys, &nrmols, &molinfo);
524 convert_harmonics(nrmols, molinfo, atype);
527 if (ir->eDisre == edrNone)
529 i = rm_interactions(F_DISRES, nrmols, molinfo);
532 set_warning_line(wi, "unknown", -1);
533 sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
534 warning_note(wi, warn_buf);
537 if (opts->bOrire == FALSE)
539 i = rm_interactions(F_ORIRES, nrmols, molinfo);
542 set_warning_line(wi, "unknown", -1);
543 sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
544 warning_note(wi, warn_buf);
548 /* Copy structures from msys to sys */
549 molinfo2mtop(nrmols, molinfo, sys);
551 gmx_mtop_finalize(sys);
553 /* COORDINATE file processing */
556 fprintf(stderr, "processing coordinates...\n");
559 get_stx_coordnum(confin, &state->natoms);
560 if (state->natoms != sys->natoms)
562 gmx_fatal(FARGS, "number of coordinates in coordinate file (%s, %d)\n"
563 " does not match topology (%s, %d)",
564 confin, state->natoms, topfile, sys->natoms);
568 /* make space for coordinates and velocities */
571 init_t_atoms(confat, state->natoms, FALSE);
572 init_state(state, state->natoms, 0, 0, 0, 0);
573 read_stx_conf(confin, title, confat, state->x, state->v, NULL, state->box);
574 /* This call fixes the box shape for runs with pressure scaling */
575 set_box_rel(ir, state);
577 nmismatch = check_atom_names(topfile, confin, sys, confat);
578 free_t_atoms(confat, TRUE);
583 sprintf(buf, "%d non-matching atom name%s\n"
584 "atom names from %s will be used\n"
585 "atom names from %s will be ignored\n",
586 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
591 fprintf(stderr, "double-checking input for internal consistency...\n");
593 double_check(ir, state->box, nint_ftype(sys, molinfo, F_CONSTR), wi);
599 gmx_mtop_atomloop_all_t aloop;
602 snew(mass, state->natoms);
603 aloop = gmx_mtop_atomloop_all_init(sys);
604 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
609 if (opts->seed == -1)
611 opts->seed = make_seed();
612 fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
614 maxwell_speed(opts->tempi, opts->seed, sys, state->v);
616 stop_cm(stdout, state->natoms, mass, state->x, state->v);
624 static void copy_state(const char *slog, t_trxframe *fr,
625 gmx_bool bReadVel, t_state *state,
630 if (fr->not_ok & FRAME_NOT_OK)
632 gmx_fatal(FARGS, "Can not start from an incomplete frame");
636 gmx_fatal(FARGS, "Did not find a frame with coordinates in file %s",
640 for (i = 0; i < state->natoms; i++)
642 copy_rvec(fr->x[i], state->x[i]);
648 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
650 for (i = 0; i < state->natoms; i++)
652 copy_rvec(fr->v[i], state->v[i]);
657 copy_mat(fr->box, state->box);
660 *use_time = fr->time;
663 static void cont_status(const char *slog, const char *ener,
664 gmx_bool bNeedVel, gmx_bool bGenVel, real fr_time,
665 t_inputrec *ir, t_state *state,
667 const output_env_t oenv)
668 /* If fr_time == -1 read the last frame available which is complete */
676 bReadVel = (bNeedVel && !bGenVel);
679 "Reading Coordinates%s and Box size from old trajectory\n",
680 bReadVel ? ", Velocities" : "");
683 fprintf(stderr, "Will read whole trajectory\n");
687 fprintf(stderr, "Will read till time %g\n", fr_time);
693 fprintf(stderr, "Velocities generated: "
694 "ignoring velocities in input trajectory\n");
696 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
700 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X | TRX_NEED_V);
706 "WARNING: Did not find a frame with velocities in file %s,\n"
707 " all velocities will be set to zero!\n\n", slog);
708 for (i = 0; i < sys->natoms; i++)
710 clear_rvec(state->v[i]);
713 /* Search for a frame without velocities */
715 read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
719 state->natoms = fr.natoms;
721 if (sys->natoms != state->natoms)
723 gmx_fatal(FARGS, "Number of atoms in Topology "
724 "is not the same as in Trajectory");
726 copy_state(slog, &fr, bReadVel, state, &use_time);
728 /* Find the appropriate frame */
729 while ((fr_time == -1 || fr.time < fr_time) &&
730 read_next_frame(oenv, fp, &fr))
732 copy_state(slog, &fr, bReadVel, state, &use_time);
737 /* Set the relative box lengths for preserving the box shape.
738 * Note that this call can lead to differences in the last bit
739 * with respect to using tpbconv to create a [TT].tpx[tt] file.
741 set_box_rel(ir, state);
743 fprintf(stderr, "Using frame at t = %g ps\n", use_time);
744 fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
746 if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
748 get_enx_state(ener, use_time, &sys->groups, ir, state);
749 preserve_box_shape(ir, state->box_rel, state->boxv);
753 static void read_posres(gmx_mtop_t *mtop, t_molinfo *molinfo, gmx_bool bTopB,
755 int rc_scaling, int ePBC,
759 gmx_bool bFirst = TRUE, *hadAtom;
765 int natoms, npbcdim = 0;
766 char warn_buf[STRLEN], title[STRLEN];
767 int a, i, ai, j, k, mb, nat_molb;
768 gmx_molblock_t *molb;
772 get_stx_coordnum(fn, &natoms);
773 if (natoms != mtop->natoms)
775 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);
776 warning(wi, warn_buf);
780 init_t_atoms(&dumat, natoms, FALSE);
781 read_stx_conf(fn, title, &dumat, x, v, NULL, box);
783 npbcdim = ePBC2npbcdim(ePBC);
785 if (rc_scaling != erscNO)
787 copy_mat(box, invbox);
788 for (j = npbcdim; j < DIM; j++)
790 clear_rvec(invbox[j]);
793 m_inv_ur0(invbox, invbox);
796 /* Copy the reference coordinates to mtop */
800 snew(hadAtom, natoms);
801 for (mb = 0; mb < mtop->nmolblock; mb++)
803 molb = &mtop->molblock[mb];
804 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
805 pr = &(molinfo[molb->type].plist[F_POSRES]);
806 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
807 if (pr->nr > 0 || prfb->nr > 0)
809 atom = mtop->moltype[molb->type].atoms.atom;
810 for (i = 0; (i < pr->nr); i++)
812 ai = pr->param[i].AI;
815 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
816 ai+1, *molinfo[molb->type].name, fn, natoms);
819 if (rc_scaling == erscCOM)
821 /* Determine the center of mass of the posres reference coordinates */
822 for (j = 0; j < npbcdim; j++)
824 sum[j] += atom[ai].m*x[a+ai][j];
826 totmass += atom[ai].m;
829 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
830 for (i = 0; (i < prfb->nr); i++)
832 ai = prfb->param[i].AI;
835 gmx_fatal(FARGS, "Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
836 ai+1, *molinfo[molb->type].name, fn, natoms);
838 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE)
840 /* Determine the center of mass of the posres reference coordinates */
841 for (j = 0; j < npbcdim; j++)
843 sum[j] += atom[ai].m*x[a+ai][j];
845 totmass += atom[ai].m;
850 molb->nposres_xA = nat_molb;
851 snew(molb->posres_xA, molb->nposres_xA);
852 for (i = 0; i < nat_molb; i++)
854 copy_rvec(x[a+i], molb->posres_xA[i]);
859 molb->nposres_xB = nat_molb;
860 snew(molb->posres_xB, molb->nposres_xB);
861 for (i = 0; i < nat_molb; i++)
863 copy_rvec(x[a+i], molb->posres_xB[i]);
869 if (rc_scaling == erscCOM)
873 gmx_fatal(FARGS, "The total mass of the position restraint atoms is 0");
875 for (j = 0; j < npbcdim; j++)
877 com[j] = sum[j]/totmass;
879 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]);
882 if (rc_scaling != erscNO)
884 for (mb = 0; mb < mtop->nmolblock; mb++)
886 molb = &mtop->molblock[mb];
887 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
888 if (molb->nposres_xA > 0 || molb->nposres_xB > 0)
890 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
891 for (i = 0; i < nat_molb; i++)
893 for (j = 0; j < npbcdim; j++)
895 if (rc_scaling == erscALL)
897 /* Convert from Cartesian to crystal coordinates */
898 xp[i][j] *= invbox[j][j];
899 for (k = j+1; k < npbcdim; k++)
901 xp[i][j] += invbox[k][j]*xp[i][k];
904 else if (rc_scaling == erscCOM)
906 /* Subtract the center of mass */
914 if (rc_scaling == erscCOM)
916 /* Convert the COM from Cartesian to crystal coordinates */
917 for (j = 0; j < npbcdim; j++)
919 com[j] *= invbox[j][j];
920 for (k = j+1; k < npbcdim; k++)
922 com[j] += invbox[k][j]*com[k];
928 free_t_atoms(&dumat, TRUE);
934 static void gen_posres(gmx_mtop_t *mtop, t_molinfo *mi,
935 char *fnA, char *fnB,
936 int rc_scaling, int ePBC,
942 read_posres (mtop, mi, FALSE, fnA, rc_scaling, ePBC, com, wi);
943 if (strcmp(fnA, fnB) != 0)
945 read_posres(mtop, mi, TRUE, fnB, rc_scaling, ePBC, comB, wi);
949 static void set_wall_atomtype(gpp_atomtype_t at, t_gromppopts *opts,
950 t_inputrec *ir, warninp_t wi)
953 char warn_buf[STRLEN];
957 fprintf(stderr, "Searching the wall atom type(s)\n");
959 for (i = 0; i < ir->nwall; i++)
961 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i], at);
962 if (ir->wall_atomtype[i] == NOTSET)
964 sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
965 warning_error(wi, warn_buf);
970 static int nrdf_internal(t_atoms *atoms)
975 for (i = 0; i < atoms->nr; i++)
977 /* Vsite ptype might not be set here yet, so also check the mass */
978 if ((atoms->atom[i].ptype == eptAtom ||
979 atoms->atom[i].ptype == eptNucleus)
980 && atoms->atom[i].m > 0)
987 case 0: nrdf = 0; break;
988 case 1: nrdf = 0; break;
989 case 2: nrdf = 1; break;
990 default: nrdf = nmass*3 - 6; break;
1009 for (i = 1; i < n-1; i++)
1011 p = 0.5*y2[i-1]+2.0;
1013 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
1014 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
1019 for (i = n-2; i >= 0; i--)
1021 y2[i] = y2[i]*y2[i+1]+u[i];
1027 interpolate1d( double xmin,
1040 a = (xmin+(ix+1)*dx-x)/dx;
1041 b = (x-xmin-ix*dx)/dx;
1043 *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;
1044 *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];
1049 setup_cmap (int grid_spacing,
1052 gmx_cmap_t * cmap_grid)
1054 double *tmp_u, *tmp_u2, *tmp_yy, *tmp_y1, *tmp_t2, *tmp_grid;
1056 int i, j, k, ii, jj, kk, idx;
1058 double dx, xmin, v, v1, v2, v12;
1061 snew(tmp_u, 2*grid_spacing);
1062 snew(tmp_u2, 2*grid_spacing);
1063 snew(tmp_yy, 2*grid_spacing);
1064 snew(tmp_y1, 2*grid_spacing);
1065 snew(tmp_t2, 2*grid_spacing*2*grid_spacing);
1066 snew(tmp_grid, 2*grid_spacing*2*grid_spacing);
1068 dx = 360.0/grid_spacing;
1069 xmin = -180.0-dx*grid_spacing/2;
1071 for (kk = 0; kk < nc; kk++)
1073 /* Compute an offset depending on which cmap we are using
1074 * Offset will be the map number multiplied with the
1075 * grid_spacing * grid_spacing * 2
1077 offset = kk * grid_spacing * grid_spacing * 2;
1079 for (i = 0; i < 2*grid_spacing; i++)
1081 ii = (i+grid_spacing-grid_spacing/2)%grid_spacing;
1083 for (j = 0; j < 2*grid_spacing; j++)
1085 jj = (j+grid_spacing-grid_spacing/2)%grid_spacing;
1086 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
1090 for (i = 0; i < 2*grid_spacing; i++)
1092 spline1d(dx, &(tmp_grid[2*grid_spacing*i]), 2*grid_spacing, tmp_u, &(tmp_t2[2*grid_spacing*i]));
1095 for (i = grid_spacing/2; i < grid_spacing+grid_spacing/2; i++)
1097 ii = i-grid_spacing/2;
1100 for (j = grid_spacing/2; j < grid_spacing+grid_spacing/2; j++)
1102 jj = j-grid_spacing/2;
1105 for (k = 0; k < 2*grid_spacing; k++)
1107 interpolate1d(xmin, dx, &(tmp_grid[2*grid_spacing*k]),
1108 &(tmp_t2[2*grid_spacing*k]), psi, &tmp_yy[k], &tmp_y1[k]);
1111 spline1d(dx, tmp_yy, 2*grid_spacing, tmp_u, tmp_u2);
1112 interpolate1d(xmin, dx, tmp_yy, tmp_u2, phi, &v, &v1);
1113 spline1d(dx, tmp_y1, 2*grid_spacing, tmp_u, tmp_u2);
1114 interpolate1d(xmin, dx, tmp_y1, tmp_u2, phi, &v2, &v12);
1116 idx = ii*grid_spacing+jj;
1117 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1118 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1119 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1120 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1126 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1130 cmap_grid->ngrid = ngrid;
1131 cmap_grid->grid_spacing = grid_spacing;
1132 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1134 snew(cmap_grid->cmapdata, ngrid);
1136 for (i = 0; i < cmap_grid->ngrid; i++)
1138 snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
1143 static int count_constraints(gmx_mtop_t *mtop, t_molinfo *mi, warninp_t wi)
1145 int count, count_mol, i, mb;
1146 gmx_molblock_t *molb;
1151 for (mb = 0; mb < mtop->nmolblock; mb++)
1154 molb = &mtop->molblock[mb];
1155 plist = mi[molb->type].plist;
1157 for (i = 0; i < F_NRE; i++)
1161 count_mol += 3*plist[i].nr;
1163 else if (interaction_function[i].flags & IF_CONSTRAINT)
1165 count_mol += plist[i].nr;
1169 if (count_mol > nrdf_internal(&mi[molb->type].atoms))
1172 "Molecule type '%s' has %d constraints.\n"
1173 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1174 *mi[molb->type].name, count_mol,
1175 nrdf_internal(&mi[molb->type].atoms));
1178 count += molb->nmol*count_mol;
1184 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1186 int i, nmiss, natoms, mt;
1188 const t_atoms *atoms;
1191 for (mt = 0; mt < sys->nmoltype; mt++)
1193 atoms = &sys->moltype[mt].atoms;
1196 for (i = 0; i < natoms; i++)
1198 q = atoms->atom[i].q;
1199 if ((get_atomtype_radius(atoms->atom[i].type, atype) == 0 ||
1200 get_atomtype_vol(atoms->atom[i].type, atype) == 0 ||
1201 get_atomtype_surftens(atoms->atom[i].type, atype) == 0 ||
1202 get_atomtype_gb_radius(atoms->atom[i].type, atype) == 0 ||
1203 get_atomtype_S_hct(atoms->atom[i].type, atype) == 0) &&
1206 fprintf(stderr, "\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1207 get_atomtype_name(atoms->atom[i].type, atype), q);
1215 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);
1220 static void check_gbsa_params(gpp_atomtype_t atype)
1224 /* If we are doing GBSA, check that we got the parameters we need
1225 * This checking is to see if there are GBSA paratmeters for all
1226 * atoms in the force field. To go around this for testing purposes
1227 * comment out the nerror++ counter temporarily
1230 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1232 if (get_atomtype_radius(i, atype) < 0 ||
1233 get_atomtype_vol(i, atype) < 0 ||
1234 get_atomtype_surftens(i, atype) < 0 ||
1235 get_atomtype_gb_radius(i, atype) < 0 ||
1236 get_atomtype_S_hct(i, atype) < 0)
1238 fprintf(stderr, "\nGB parameter(s) missing or negative for atom type '%s'\n",
1239 get_atomtype_name(i, atype));
1246 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);
1251 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1254 real verletbuf_drift,
1259 verletbuf_list_setup_t ls;
1262 char warn_buf[STRLEN];
1265 for (i = 0; i < ir->opts.ngtc; i++)
1267 if (ir->opts.ref_t[i] < 0)
1269 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.");
1273 ref_T = max(ref_T, ir->opts.ref_t[i]);
1277 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n", verletbuf_drift, ref_T);
1279 for (i = 0; i < ir->opts.ngtc; i++)
1281 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1283 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.",
1284 ir->opts.nrdf[i], ir->opts.ref_t[i], ref_T);
1285 warning_note(wi, warn_buf);
1289 /* Calculate the buffer size for simple atom vs atoms list */
1290 ls.cluster_size_i = 1;
1291 ls.cluster_size_j = 1;
1292 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1293 &ls, &n_nonlin_vsite, &rlist_1x1);
1295 /* Set the pair-list buffer size in ir */
1296 verletbuf_get_list_setup(FALSE, &ls);
1297 calc_verlet_buffer_size(mtop, det(box), ir, verletbuf_drift,
1298 &ls, &n_nonlin_vsite, &ir->rlist);
1300 if (n_nonlin_vsite > 0)
1302 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);
1303 warning_note(wi, warn_buf);
1306 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1307 1, 1, rlist_1x1, rlist_1x1-max(ir->rvdw, ir->rcoulomb));
1309 ir->rlistlong = ir->rlist;
1310 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1311 ls.cluster_size_i, ls.cluster_size_j,
1312 ir->rlist, ir->rlist-max(ir->rvdw, ir->rcoulomb));
1314 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC, box))
1316 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)));
1320 int gmx_grompp(int argc, char *argv[])
1322 static const char *desc[] = {
1323 "[THISMODULE] (the gromacs preprocessor)",
1324 "reads a molecular topology file, checks the validity of the",
1325 "file, expands the topology from a molecular description to an atomic",
1326 "description. The topology file contains information about",
1327 "molecule types and the number of molecules, the preprocessor",
1328 "copies each molecule as needed. ",
1329 "There is no limitation on the number of molecule types. ",
1330 "Bonds and bond-angles can be converted into constraints, separately",
1331 "for hydrogens and heavy atoms.",
1332 "Then a coordinate file is read and velocities can be generated",
1333 "from a Maxwellian distribution if requested.",
1334 "[THISMODULE] also reads parameters for [gmx-mdrun] ",
1335 "(eg. number of MD steps, time step, cut-off), and others such as",
1336 "NEMD parameters, which are corrected so that the net acceleration",
1338 "Eventually a binary file is produced that can serve as the sole input",
1339 "file for the MD program.[PAR]",
1341 "[THISMODULE] uses the atom names from the topology file. The atom names",
1342 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1343 "warnings when they do not match the atom names in the topology.",
1344 "Note that the atom names are irrelevant for the simulation as",
1345 "only the atom types are used for generating interaction parameters.[PAR]",
1347 "[THISMODULE] uses a built-in preprocessor to resolve includes, macros, ",
1348 "etc. The preprocessor supports the following keywords:[PAR]",
1349 "#ifdef VARIABLE[BR]",
1350 "#ifndef VARIABLE[BR]",
1353 "#define VARIABLE[BR]",
1354 "#undef VARIABLE[BR]"
1355 "#include \"filename\"[BR]",
1356 "#include <filename>[PAR]",
1357 "The functioning of these statements in your topology may be modulated by",
1358 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1359 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1360 "include = -I/home/john/doe[tt][BR]",
1361 "For further information a C-programming textbook may help you out.",
1362 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1363 "topology file written out so that you can verify its contents.[PAR]",
1365 /* cpp has been unnecessary for some time, hasn't it?
1366 "If your system does not have a C-preprocessor, you can still",
1367 "use [TT]grompp[tt], but you do not have access to the features ",
1368 "from the cpp. Command line options to the C-preprocessor can be given",
1369 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1372 "When using position restraints a file with restraint coordinates",
1373 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1374 "with respect to the conformation from the [TT]-c[tt] option.",
1375 "For free energy calculation the the coordinates for the B topology",
1376 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1377 "those of the A topology.[PAR]",
1379 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1380 "The last frame with coordinates and velocities will be read,",
1381 "unless the [TT]-time[tt] option is used. Only if this information",
1382 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1383 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1384 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1385 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1388 "[THISMODULE] can be used to restart simulations (preserving",
1389 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1390 "However, for simply changing the number of run steps to extend",
1391 "a run, using [gmx-tpbconv] is more convenient than [THISMODULE].",
1392 "You then supply the old checkpoint file directly to [gmx-mdrun]",
1393 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1394 "like output frequency, then supplying the checkpoint file to",
1395 "[THISMODULE] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1396 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1398 "By default, all bonded interactions which have constant energy due to",
1399 "virtual site constructions will be removed. If this constant energy is",
1400 "not zero, this will result in a shift in the total energy. All bonded",
1401 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1402 "all constraints for distances which will be constant anyway because",
1403 "of virtual site constructions will be removed. If any constraints remain",
1404 "which involve virtual sites, a fatal error will result.[PAR]"
1406 "To verify your run input file, please take note of all warnings",
1407 "on the screen, and correct where necessary. Do also look at the contents",
1408 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1409 "the input that [THISMODULE] has read. If in doubt, you can start [THISMODULE]",
1410 "with the [TT]-debug[tt] option which will give you more information",
1411 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1412 "can see the contents of the run input file with the [gmx-dump]",
1413 "program. [gmx-check] can be used to compare the contents of two",
1414 "run input files.[PAR]"
1416 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1417 "by [THISMODULE] that otherwise halt output. In some cases, warnings are",
1418 "harmless, but usually they are not. The user is advised to carefully",
1419 "interpret the output messages before attempting to bypass them with",
1426 gpp_atomtype_t atype;
1428 int natoms, nvsite, comb, mt;
1432 real max_spacing, fudgeQQ;
1434 char fn[STRLEN], fnB[STRLEN];
1435 const char *mdparin;
1437 gmx_bool bNeedVel, bGenVel;
1438 gmx_bool have_atomnumber;
1440 t_params *gb_plist = NULL;
1441 gmx_genborn_t *born = NULL;
1443 gmx_bool bVerbose = FALSE;
1445 char warn_buf[STRLEN];
1448 { efMDP, NULL, NULL, ffREAD },
1449 { efMDP, "-po", "mdout", ffWRITE },
1450 { efSTX, "-c", NULL, ffREAD },
1451 { efSTX, "-r", NULL, ffOPTRD },
1452 { efSTX, "-rb", NULL, ffOPTRD },
1453 { efNDX, NULL, NULL, ffOPTRD },
1454 { efTOP, NULL, NULL, ffREAD },
1455 { efTOP, "-pp", "processed", ffOPTWR },
1456 { efTPX, "-o", NULL, ffWRITE },
1457 { efTRN, "-t", NULL, ffOPTRD },
1458 { efEDR, "-e", NULL, ffOPTRD },
1459 { efTRN, "-ref", "rotref", ffOPTRW }
1461 #define NFILE asize(fnm)
1463 /* Command line options */
1464 static gmx_bool bRenum = TRUE;
1465 static gmx_bool bRmVSBds = TRUE, bZero = FALSE;
1466 static int i, maxwarn = 0;
1467 static real fr_time = -1;
1469 { "-v", FALSE, etBOOL, {&bVerbose},
1470 "Be loud and noisy" },
1471 { "-time", FALSE, etREAL, {&fr_time},
1472 "Take frame at or first after this time." },
1473 { "-rmvsbds", FALSE, etBOOL, {&bRmVSBds},
1474 "Remove constant bonded interactions with virtual sites" },
1475 { "-maxwarn", FALSE, etINT, {&maxwarn},
1476 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1477 { "-zero", FALSE, etBOOL, {&bZero},
1478 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1479 { "-renum", FALSE, etBOOL, {&bRenum},
1480 "Renumber atomtypes and minimize number of atomtypes" }
1483 /* Initiate some variables */
1488 /* Parse the command line */
1489 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
1490 asize(desc), desc, 0, NULL, &oenv))
1495 wi = init_warning(TRUE, maxwarn);
1497 /* PARAMETER file processing */
1498 mdparin = opt2fn("-f", NFILE, fnm);
1499 set_warning_line(wi, mdparin, -1);
1500 get_ir(mdparin, opt2fn("-po", NFILE, fnm), ir, opts, wi);
1504 fprintf(stderr, "checking input for internal consistency...\n");
1506 check_ir(mdparin, ir, opts, wi);
1508 if (ir->ld_seed == -1)
1510 ir->ld_seed = make_seed();
1511 fprintf(stderr, "Setting the LD random seed to %d\n", ir->ld_seed);
1514 if (ir->expandedvals->lmc_seed == -1)
1516 ir->expandedvals->lmc_seed = make_seed();
1517 fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
1520 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1521 bGenVel = (bNeedVel && opts->bGenVel);
1522 if (bGenVel && ir->bContinuation)
1525 "Generating velocities is inconsistent with attempting "
1526 "to continue a previous run. Choose only one of "
1527 "gen-vel = yes and continuation = yes.");
1528 warning_error(wi, warn_buf);
1534 atype = init_atomtype();
1537 pr_symtab(debug, 0, "Just opened", &sys->symtab);
1540 strcpy(fn, ftp2fn(efTOP, NFILE, fnm));
1541 if (!gmx_fexist(fn))
1543 gmx_fatal(FARGS, "%s does not exist", fn);
1545 new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm),
1546 opts, ir, bZero, bGenVel, bVerbose, &state,
1547 atype, sys, &nmi, &mi, plist, &comb, &reppow, &fudgeQQ,
1553 pr_symtab(debug, 0, "After new_status", &sys->symtab);
1557 /* set parameters for virtual site construction (not for vsiten) */
1558 for (mt = 0; mt < sys->nmoltype; mt++)
1561 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1563 /* now throw away all obsolete bonds, angles and dihedrals: */
1564 /* note: constraints are ALWAYS removed */
1567 for (mt = 0; mt < sys->nmoltype; mt++)
1569 clean_vsite_bondeds(mi[mt].plist, sys->moltype[mt].atoms.nr, bRmVSBds);
1573 if (ir->cutoff_scheme == ecutsVERLET)
1575 fprintf(stderr, "Removing all charge groups because cutoff-scheme=%s\n",
1576 ecutscheme_names[ir->cutoff_scheme]);
1578 /* Remove all charge groups */
1579 gmx_mtop_remove_chargegroups(sys);
1582 if (count_constraints(sys, mi, wi) && (ir->eConstrAlg == econtSHAKE))
1584 if (ir->eI == eiCG || ir->eI == eiLBFGS)
1586 sprintf(warn_buf, "Can not do %s with %s, use %s",
1587 EI(ir->eI), econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1588 warning_error(wi, warn_buf);
1590 if (ir->bPeriodicMols)
1592 sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
1593 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
1594 warning_error(wi, warn_buf);
1598 if (EI_SD (ir->eI) && ir->etc != etcNO)
1600 warning_note(wi, "Temperature coupling is ignored with SD integrators.");
1603 /* If we are doing QM/MM, check that we got the atom numbers */
1604 have_atomnumber = TRUE;
1605 for (i = 0; i < get_atomtype_ntypes(atype); i++)
1607 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i, atype) >= 0);
1609 if (!have_atomnumber && ir->bQMMM)
1613 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1614 "field you are using does not contain atom numbers fields. This is an\n"
1615 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1616 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1617 "an integer just before the mass column in ffXXXnb.itp.\n"
1618 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1623 if ((ir->adress->const_wf > 1) || (ir->adress->const_wf < 0))
1625 warning_error(wi, "AdResS contant weighting function should be between 0 and 1\n\n");
1627 /** \TODO check size of ex+hy width against box size */
1630 /* Check for errors in the input now, since they might cause problems
1631 * during processing further down.
1633 check_warning_error(wi, FARGS);
1635 if (opt2bSet("-r", NFILE, fnm))
1637 sprintf(fn, "%s", opt2fn("-r", NFILE, fnm));
1641 sprintf(fn, "%s", opt2fn("-c", NFILE, fnm));
1643 if (opt2bSet("-rb", NFILE, fnm))
1645 sprintf(fnB, "%s", opt2fn("-rb", NFILE, fnm));
1652 if (nint_ftype(sys, mi, F_POSRES) > 0 || nint_ftype(sys, mi, F_FBPOSRES) > 0)
1656 fprintf(stderr, "Reading position restraint coords from %s", fn);
1657 if (strcmp(fn, fnB) == 0)
1659 fprintf(stderr, "\n");
1663 fprintf(stderr, " and %s\n", fnB);
1666 gen_posres(sys, mi, fn, fnB,
1667 ir->refcoord_scaling, ir->ePBC,
1668 ir->posres_com, ir->posres_comB,
1672 /* If we are using CMAP, setup the pre-interpolation grid */
1673 if (plist->ncmap > 0)
1675 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1676 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap, &sys->ffparams.cmap_grid);
1679 set_wall_atomtype(atype, opts, ir, wi);
1682 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1683 ntype = get_atomtype_ntypes(atype);
1686 if (ir->implicit_solvent != eisNO)
1688 /* Now we have renumbered the atom types, we can check the GBSA params */
1689 check_gbsa_params(atype);
1691 /* Check that all atoms that have charge and/or LJ-parameters also have
1692 * sensible GB-parameters
1694 check_gbsa_params_charged(sys, atype);
1697 /* PELA: Copy the atomtype data to the topology atomtype list */
1698 copy_atomtype_atomtypes(atype, &(sys->atomtypes));
1702 pr_symtab(debug, 0, "After renum_atype", &sys->symtab);
1707 fprintf(stderr, "converting bonded parameters...\n");
1710 ntype = get_atomtype_ntypes(atype);
1711 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1715 pr_symtab(debug, 0, "After convert_params", &sys->symtab);
1718 /* set ptype to VSite for virtual sites */
1719 for (mt = 0; mt < sys->nmoltype; mt++)
1721 set_vsites_ptype(FALSE, &sys->moltype[mt]);
1725 pr_symtab(debug, 0, "After virtual sites", &sys->symtab);
1727 /* Check velocity for virtual sites and shells */
1730 check_vel(sys, state.v);
1736 for (i = 0; i < sys->nmoltype; i++)
1738 check_cg_sizes(ftp2fn(efTOP, NFILE, fnm), &sys->moltype[i].cgs, wi);
1741 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1743 check_bonds_timestep(sys, ir->delta_t, wi);
1746 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1748 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.");
1751 check_warning_error(wi, FARGS);
1755 fprintf(stderr, "initialising group options...\n");
1757 do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm),
1759 bGenVel ? state.v : NULL,
1762 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1765 if (EI_DYNAMICS(ir->eI) &&
1766 !(EI_MD(ir->eI) && ir->etc == etcNO) &&
1767 inputrec2nboundeddim(ir) == 3)
1769 set_verlet_buffer(sys, ir, state.box, ir->verletbuf_drift, wi);
1773 /* Init the temperature coupling state */
1774 init_gtc_state(&state, ir->opts.ngtc, 0, ir->opts.nhchainlength); /* need to add nnhpres here? */
1778 fprintf(stderr, "Checking consistency between energy and charge groups...\n");
1780 check_eg_vs_cg(sys);
1784 pr_symtab(debug, 0, "After index", &sys->symtab);
1786 triple_check(mdparin, ir, sys, wi);
1787 close_symtab(&sys->symtab);
1790 pr_symtab(debug, 0, "After close", &sys->symtab);
1793 /* make exclusions between QM atoms */
1796 if (ir->QMMMscheme == eQMMMschemenormal && ir->ns_type == ensSIMPLE)
1798 gmx_fatal(FARGS, "electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1802 generate_qmexcl(sys, ir, wi);
1806 if (ftp2bSet(efTRN, NFILE, fnm))
1810 fprintf(stderr, "getting data from old trajectory ...\n");
1812 cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm),
1813 bNeedVel, bGenVel, fr_time, ir, &state, sys, oenv);
1816 if (ir->ePBC == epbcXY && ir->nwall != 2)
1818 clear_rvec(state.box[ZZ]);
1821 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1823 set_warning_line(wi, mdparin, -1);
1824 check_chargegroup_radii(sys, ir, state.x, wi);
1827 if (EEL_FULL(ir->coulombtype))
1829 /* Calculate the optimal grid dimensions */
1830 copy_mat(state.box, box);
1831 if (ir->ePBC == epbcXY && ir->nwall == 2)
1833 svmul(ir->wall_ewald_zfac, box[ZZ], box[ZZ]);
1835 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1837 /* Mark fourier_spacing as not used */
1838 ir->fourier_spacing = 0;
1840 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1842 set_warning_line(wi, mdparin, -1);
1843 warning_error(wi, "Some of the Fourier grid sizes are set, but all of them need to be set.");
1845 max_spacing = calc_grid(stdout, box, ir->fourier_spacing,
1846 &(ir->nkx), &(ir->nky), &(ir->nkz));
1849 /* MRS: eventually figure out better logic for initializing the fep
1850 values that makes declaring the lambda and declaring the state not
1851 potentially conflict if not handled correctly. */
1852 if (ir->efep != efepNO)
1854 state.fep_state = ir->fepvals->init_fep_state;
1855 for (i = 0; i < efptNR; i++)
1857 /* init_lambda trumps state definitions*/
1858 if (ir->fepvals->init_lambda >= 0)
1860 state.lambda[i] = ir->fepvals->init_lambda;
1864 if (ir->fepvals->all_lambda[i] == NULL)
1866 gmx_fatal(FARGS, "Values of lambda not set for a free energy calculation!");
1870 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1876 if (ir->ePull != epullNO)
1878 set_pull_init(ir, sys, state.x, state.box, state.lambda[efptMASS], oenv, opts->pull_start);
1883 set_reference_positions(ir->rot, state.x, state.box,
1884 opt2fn("-ref", NFILE, fnm), opt2bSet("-ref", NFILE, fnm),
1888 /* reset_multinr(sys); */
1890 if (EEL_PME(ir->coulombtype))
1892 float ratio = pme_load_estimate(sys, ir, state.box);
1893 fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
1894 /* With free energy we might need to do PME both for the A and B state
1895 * charges. This will double the cost, but the optimal performance will
1896 * then probably be at a slightly larger cut-off and grid spacing.
1898 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1899 (ir->efep != efepNO && ratio > 2.0/3.0))
1902 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1903 "and for highly parallel simulations between 0.25 and 0.33,\n"
1904 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1905 if (ir->efep != efepNO)
1908 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1914 char warn_buf[STRLEN];
1915 double cio = compute_io(ir, sys->natoms, &sys->groups, F_NRE, 1);
1916 sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
1919 set_warning_line(wi, mdparin, -1);
1920 warning_note(wi, warn_buf);
1924 printf("%s\n", warn_buf);
1930 fprintf(stderr, "writing run input file...\n");
1933 done_warning(wi, FARGS);
1935 write_tpx_state(ftp2fn(efTPX, NFILE, fnm), ir, &state, sys);