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>
64 #include "sortwater.h"
66 #include "gmx_fatal.h"
72 #include "vsite_parm.h"
78 #include "compute_io.h"
79 #include "gpp_atomtype.h"
80 #include "gpp_tomorse.h"
81 #include "mtop_util.h"
83 #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++) {
92 n += mols[i].plist[ifunc].nr;
93 mols[i].plist[ifunc].nr=0;
98 static int check_atom_names(const char *fn1, const char *fn2,
99 gmx_mtop_t *mtop, t_atoms *at)
101 int mb,m,i,j,nmismatch;
103 #define MAXMISMATCH 20
105 if (mtop->natoms != at->nr)
106 gmx_incons("comparing atom names");
110 for(mb=0; mb<mtop->nmolblock; mb++) {
111 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
112 for(m=0; m<mtop->molblock[mb].nmol; m++) {
113 for(j=0; j < tat->nr; j++) {
114 if (strcmp( *(tat->atomname[j]) , *(at->atomname[i]) ) != 0) {
115 if (nmismatch < MAXMISMATCH) {
117 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
118 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
119 } else if (nmismatch == MAXMISMATCH) {
120 fprintf(stderr,"(more than %d non-matching atom names)\n",MAXMISMATCH);
132 static void check_eg_vs_cg(gmx_mtop_t *mtop)
134 int astart,mb,m,cg,j,firstj;
135 unsigned char firsteg,eg;
138 /* Go through all the charge groups and make sure all their
139 * atoms are in the same energy group.
143 for(mb=0; mb<mtop->nmolblock; mb++) {
144 molt = &mtop->moltype[mtop->molblock[mb].type];
145 for(m=0; m<mtop->molblock[mb].nmol; m++) {
146 for(cg=0; cg<molt->cgs.nr;cg++) {
147 /* Get the energy group of the first atom in this charge group */
148 firstj = astart + molt->cgs.index[cg];
149 firsteg = ggrpnr(&mtop->groups,egcENER,firstj);
150 for(j=molt->cgs.index[cg]+1;j<molt->cgs.index[cg+1];j++) {
151 eg = ggrpnr(&mtop->groups,egcENER,astart+j);
153 gmx_fatal(FARGS,"atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
154 firstj+1,astart+j+1,cg+1,*molt->name);
158 astart += molt->atoms.nr;
163 static void check_cg_sizes(const char *topfn,t_block *cgs,warninp_t wi)
166 char warn_buf[STRLEN];
169 for(cg=0; cg<cgs->nr; cg++)
171 maxsize = max(maxsize,cgs->index[cg+1]-cgs->index[cg]);
174 if (maxsize > MAX_CHARGEGROUP_SIZE)
176 gmx_fatal(FARGS,"The largest charge group contains %d atoms. The maximum is %d.",maxsize,MAX_CHARGEGROUP_SIZE);
178 else if (maxsize > 10)
180 set_warning_line(wi,topfn,-1);
182 "The largest charge group contains %d atoms.\n"
183 "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"
184 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
185 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
187 warning_note(wi,warn_buf);
191 static void check_bonds_timestep(gmx_mtop_t *mtop,double dt,warninp_t wi)
193 /* This check is not intended to ensure accurate integration,
194 * rather it is to signal mistakes in the mdp settings.
195 * A common mistake is to forget to turn on constraints
196 * for MD after energy minimization with flexible bonds.
197 * This check can also detect too large time steps for flexible water
198 * models, but such errors will often be masked by the constraints
199 * mdp options, which turns flexible water into water with bond constraints,
200 * but without an angle constraint. Unfortunately such incorrect use
201 * of water models can not easily be detected without checking
202 * for specific model names.
204 * The stability limit of leap-frog or velocity verlet is 4.44 steps
205 * per oscillational period.
206 * But accurate bonds distributions are lost far before that limit.
207 * To allow relatively common schemes (although not common with Gromacs)
208 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
209 * we set the note limit to 10.
211 int min_steps_warn=5;
212 int min_steps_note=10;
215 gmx_moltype_t *moltype,*w_moltype;
217 t_ilist *ilist,*ilb,*ilc,*ils;
219 int i,a1,a2,w_a1,w_a2,j;
220 real twopi2,limit2,fc,re,m1,m2,period2,w_period2;
221 gmx_bool bFound,bWater,bWarn;
222 char warn_buf[STRLEN];
224 ip = mtop->ffparams.iparams;
226 twopi2 = sqr(2*M_PI);
228 limit2 = sqr(min_steps_note*dt);
234 for(molt=0; molt<mtop->nmoltype; molt++)
236 moltype = &mtop->moltype[molt];
237 atom = moltype->atoms.atom;
238 ilist = moltype->ilist;
239 ilc = &ilist[F_CONSTR];
240 ils = &ilist[F_SETTLE];
241 for(ftype=0; ftype<F_NRE; ftype++)
243 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
249 for(i=0; i<ilb->nr; i+=3)
251 fc = ip[ilb->iatoms[i]].harmonic.krA;
252 re = ip[ilb->iatoms[i]].harmonic.rA;
253 if (ftype == F_G96BONDS)
255 /* Convert squared sqaure fc to harmonic fc */
258 a1 = ilb->iatoms[i+1];
259 a2 = ilb->iatoms[i+2];
262 if (fc > 0 && m1 > 0 && m2 > 0)
264 period2 = twopi2*m1*m2/((m1 + m2)*fc);
268 period2 = GMX_FLOAT_MAX;
272 fprintf(debug,"fc %g m1 %g m2 %g period %g\n",
273 fc,m1,m2,sqrt(period2));
275 if (period2 < limit2)
278 for(j=0; j<ilc->nr; j+=3)
280 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
281 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
286 for(j=0; j<ils->nr; j+=4)
288 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
289 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
295 (w_moltype == NULL || period2 < w_period2))
307 if (w_moltype != NULL)
309 bWarn = (w_period2 < sqr(min_steps_warn*dt));
310 /* A check that would recognize most water models */
311 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
312 w_moltype->atoms.nr <= 5);
313 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"
316 w_a1+1,*w_moltype->atoms.atomname[w_a1],
317 w_a2+1,*w_moltype->atoms.atomname[w_a2],
318 sqrt(w_period2),bWarn ? min_steps_warn : min_steps_note,dt,
320 "Maybe you asked for fexible water." :
321 "Maybe you forgot to change the constraints mdp option.");
324 warning(wi,warn_buf);
328 warning_note(wi,warn_buf);
333 static void check_vel(gmx_mtop_t *mtop,rvec v[])
335 gmx_mtop_atomloop_all_t aloop;
339 aloop = gmx_mtop_atomloop_all_init(mtop);
340 while (gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
341 if (atom->ptype == eptShell ||
342 atom->ptype == eptBond ||
343 atom->ptype == eptVSite) {
349 static gmx_bool nint_ftype(gmx_mtop_t *mtop,t_molinfo *mi,int ftype)
354 for(mb=0; mb<mtop->nmolblock; mb++) {
355 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
361 /* This routine reorders the molecule type array
362 * in the order of use in the molblocks,
363 * unused molecule types are deleted.
365 static void renumber_moltypes(gmx_mtop_t *sys,
366 int *nmolinfo,t_molinfo **molinfo)
372 snew(order,*nmolinfo);
374 for(mb=0; mb<sys->nmolblock; mb++) {
375 for(i=0; i<norder; i++) {
376 if (order[i] == sys->molblock[mb].type) {
381 /* This type did not occur yet, add it */
382 order[norder] = sys->molblock[mb].type;
383 /* Renumber the moltype in the topology */
386 sys->molblock[mb].type = i;
389 /* We still need to reorder the molinfo structs */
391 for(mi=0; mi<*nmolinfo; mi++) {
392 for(i=0; i<norder; i++) {
393 if (order[i] == mi) {
398 done_mi(&(*molinfo)[mi]);
400 minew[i] = (*molinfo)[mi];
409 static void molinfo2mtop(int nmi,t_molinfo *mi,gmx_mtop_t *mtop)
414 mtop->nmoltype = nmi;
415 snew(mtop->moltype,nmi);
416 for(m=0; m<nmi; m++) {
417 molt = &mtop->moltype[m];
418 molt->name = mi[m].name;
419 molt->atoms = mi[m].atoms;
420 /* ilists are copied later */
421 molt->cgs = mi[m].cgs;
422 molt->excls = mi[m].excls;
427 new_status(const char *topfile,const char *topppfile,const char *confin,
428 t_gromppopts *opts,t_inputrec *ir,gmx_bool bZero,
429 gmx_bool bGenVel,gmx_bool bVerbose,t_state *state,
430 gpp_atomtype_t atype,gmx_mtop_t *sys,
431 int *nmi,t_molinfo **mi,t_params plist[],
432 int *comb,double *reppow,real *fudgeQQ,
436 t_molinfo *molinfo=NULL;
438 gmx_molblock_t *molblock,*molbs;
440 int mb,i,nrmols,nmismatch;
443 char warn_buf[STRLEN];
447 /* Set gmx_boolean for GB */
448 if(ir->implicit_solvent)
451 /* TOPOLOGY processing */
452 sys->name = do_top(bVerbose,topfile,topppfile,opts,bZero,&(sys->symtab),
453 plist,comb,reppow,fudgeQQ,
454 atype,&nrmols,&molinfo,ir,
455 &nmolblock,&molblock,bGB,
459 snew(sys->molblock,nmolblock);
462 for(mb=0; mb<nmolblock; mb++) {
463 if (sys->nmolblock > 0 &&
464 molblock[mb].type == sys->molblock[sys->nmolblock-1].type) {
465 /* Merge consecutive blocks with the same molecule type */
466 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
467 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
468 } else if (molblock[mb].nmol > 0) {
469 /* Add a new molblock to the topology */
470 molbs = &sys->molblock[sys->nmolblock];
471 *molbs = molblock[mb];
472 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
473 molbs->nposres_xA = 0;
474 molbs->nposres_xB = 0;
475 sys->natoms += molbs->nmol*molbs->natoms_mol;
479 if (sys->nmolblock == 0) {
480 gmx_fatal(FARGS,"No molecules were defined in the system");
483 renumber_moltypes(sys,&nrmols,&molinfo);
486 convert_harmonics(nrmols,molinfo,atype);
488 if (ir->eDisre == edrNone) {
489 i = rm_interactions(F_DISRES,nrmols,molinfo);
491 set_warning_line(wi,"unknown",-1);
492 sprintf(warn_buf,"disre = no, removed %d distance restraints",i);
493 warning_note(wi,warn_buf);
496 if (opts->bOrire == FALSE) {
497 i = rm_interactions(F_ORIRES,nrmols,molinfo);
499 set_warning_line(wi,"unknown",-1);
500 sprintf(warn_buf,"orire = no, removed %d orientation restraints",i);
501 warning_note(wi,warn_buf);
505 /* Copy structures from msys to sys */
506 molinfo2mtop(nrmols,molinfo,sys);
508 gmx_mtop_finalize(sys);
510 /* COORDINATE file processing */
512 fprintf(stderr,"processing coordinates...\n");
514 get_stx_coordnum(confin,&state->natoms);
515 if (state->natoms != sys->natoms)
516 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
517 " does not match topology (%s, %d)",
518 confin,state->natoms,topfile,sys->natoms);
520 /* make space for coordinates and velocities */
523 init_t_atoms(confat,state->natoms,FALSE);
524 init_state(state,state->natoms,0,0,0,0);
525 read_stx_conf(confin,title,confat,state->x,state->v,NULL,state->box);
526 /* This call fixes the box shape for runs with pressure scaling */
527 set_box_rel(ir,state);
529 nmismatch = check_atom_names(topfile, confin, sys, confat);
530 free_t_atoms(confat,TRUE);
534 sprintf(buf,"%d non-matching atom name%s\n"
535 "atom names from %s will be used\n"
536 "atom names from %s will be ignored\n",
537 nmismatch,(nmismatch == 1) ? "" : "s",topfile,confin);
541 fprintf(stderr,"double-checking input for internal consistency...\n");
542 double_check(ir,state->box,nint_ftype(sys,molinfo,F_CONSTR),wi);
547 gmx_mtop_atomloop_all_t aloop;
550 snew(mass,state->natoms);
551 aloop = gmx_mtop_atomloop_all_init(sys);
552 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
556 if (opts->seed == -1) {
557 opts->seed = make_seed();
558 fprintf(stderr,"Setting gen_seed to %d\n",opts->seed);
560 maxwell_speed(opts->tempi,opts->seed,sys,state->v);
562 stop_cm(stdout,state->natoms,mass,state->x,state->v);
570 static void copy_state(const char *slog,t_trxframe *fr,
571 gmx_bool bReadVel,t_state *state,
576 if (fr->not_ok & FRAME_NOT_OK)
578 gmx_fatal(FARGS,"Can not start from an incomplete frame");
582 gmx_fatal(FARGS,"Did not find a frame with coordinates in file %s",
586 for(i=0; i<state->natoms; i++)
588 copy_rvec(fr->x[i],state->x[i]);
594 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
596 for(i=0; i<state->natoms; i++)
598 copy_rvec(fr->v[i],state->v[i]);
603 copy_mat(fr->box,state->box);
606 *use_time = fr->time;
609 static void cont_status(const char *slog,const char *ener,
610 gmx_bool bNeedVel,gmx_bool bGenVel, real fr_time,
611 t_inputrec *ir,t_state *state,
613 const output_env_t oenv)
614 /* If fr_time == -1 read the last frame available which is complete */
622 bReadVel = (bNeedVel && !bGenVel);
625 "Reading Coordinates%s and Box size from old trajectory\n",
626 bReadVel ? ", Velocities" : "");
629 fprintf(stderr,"Will read whole trajectory\n");
633 fprintf(stderr,"Will read till time %g\n",fr_time);
639 fprintf(stderr,"Velocities generated: "
640 "ignoring velocities in input trajectory\n");
642 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
646 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X | TRX_NEED_V);
652 "WARNING: Did not find a frame with velocities in file %s,\n"
653 " all velocities will be set to zero!\n\n",slog);
654 for(i=0; i<sys->natoms; i++)
656 clear_rvec(state->v[i]);
659 /* Search for a frame without velocities */
661 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
665 state->natoms = fr.natoms;
667 if (sys->natoms != state->natoms)
669 gmx_fatal(FARGS,"Number of atoms in Topology "
670 "is not the same as in Trajectory");
672 copy_state(slog,&fr,bReadVel,state,&use_time);
674 /* Find the appropriate frame */
675 while ((fr_time == -1 || fr.time < fr_time) &&
676 read_next_frame(oenv,fp,&fr))
678 copy_state(slog,&fr,bReadVel,state,&use_time);
683 /* Set the relative box lengths for preserving the box shape.
684 * Note that this call can lead to differences in the last bit
685 * with respect to using tpbconv to create a [TT].tpx[tt] file.
687 set_box_rel(ir,state);
689 fprintf(stderr,"Using frame at t = %g ps\n",use_time);
690 fprintf(stderr,"Starting time for run is %g ps\n",ir->init_t);
692 if ((ir->epc != epcNO || ir->etc ==etcNOSEHOOVER) && ener)
694 get_enx_state(ener,use_time,&sys->groups,ir,state);
695 preserve_box_shape(ir,state->box_rel,state->boxv);
699 static void read_posres(gmx_mtop_t *mtop,t_molinfo *molinfo,gmx_bool bTopB,
701 int rc_scaling, int ePBC,
705 gmx_bool bFirst = TRUE, *hadAtom;
711 int natoms,npbcdim=0;
712 char warn_buf[STRLEN],title[STRLEN];
713 int a,i,ai,j,k,mb,nat_molb;
714 gmx_molblock_t *molb;
718 get_stx_coordnum(fn,&natoms);
719 if (natoms != mtop->natoms) {
720 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);
721 warning(wi,warn_buf);
725 init_t_atoms(&dumat,natoms,FALSE);
726 read_stx_conf(fn,title,&dumat,x,v,NULL,box);
728 npbcdim = ePBC2npbcdim(ePBC);
730 if (rc_scaling != erscNO) {
731 copy_mat(box,invbox);
732 for(j=npbcdim; j<DIM; j++) {
733 clear_rvec(invbox[j]);
736 m_inv_ur0(invbox,invbox);
739 /* Copy the reference coordinates to mtop */
743 snew(hadAtom,natoms);
744 for(mb=0; mb<mtop->nmolblock; mb++) {
745 molb = &mtop->molblock[mb];
746 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
747 pr = &(molinfo[molb->type].plist[F_POSRES]);
748 prfb = &(molinfo[molb->type].plist[F_FBPOSRES]);
749 if (pr->nr > 0 || prfb->nr > 0) {
750 atom = mtop->moltype[molb->type].atoms.atom;
751 for(i=0; (i<pr->nr); i++) {
754 gmx_fatal(FARGS,"Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
755 ai+1,*molinfo[molb->type].name,fn,natoms);
758 if (rc_scaling == erscCOM) {
759 /* Determine the center of mass of the posres reference coordinates */
760 for(j=0; j<npbcdim; j++) {
761 sum[j] += atom[ai].m*x[a+ai][j];
763 totmass += atom[ai].m;
766 /* Same for flat-bottomed posres, but do not count an atom twice for COM */
767 for(i=0; (i<prfb->nr); i++) {
768 ai=prfb->param[i].AI;
770 gmx_fatal(FARGS,"Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
771 ai+1,*molinfo[molb->type].name,fn,natoms);
773 if (rc_scaling == erscCOM && hadAtom[ai] == FALSE) {
774 /* Determine the center of mass of the posres reference coordinates */
775 for(j=0; j<npbcdim; j++) {
776 sum[j] += atom[ai].m*x[a+ai][j];
778 totmass += atom[ai].m;
782 molb->nposres_xA = nat_molb;
783 snew(molb->posres_xA,molb->nposres_xA);
784 for(i=0; i<nat_molb; i++) {
785 copy_rvec(x[a+i],molb->posres_xA[i]);
788 molb->nposres_xB = nat_molb;
789 snew(molb->posres_xB,molb->nposres_xB);
790 for(i=0; i<nat_molb; i++) {
791 copy_rvec(x[a+i],molb->posres_xB[i]);
797 if (rc_scaling == erscCOM) {
799 gmx_fatal(FARGS,"The total mass of the position restraint atoms is 0");
800 for(j=0; j<npbcdim; j++)
801 com[j] = sum[j]/totmass;
802 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]);
805 if (rc_scaling != erscNO) {
806 for(mb=0; mb<mtop->nmolblock; mb++) {
807 molb = &mtop->molblock[mb];
808 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
809 if (molb->nposres_xA > 0 || molb->nposres_xB > 0) {
810 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
811 for(i=0; i<nat_molb; i++) {
812 for(j=0; j<npbcdim; j++) {
813 if (rc_scaling == erscALL) {
814 /* Convert from Cartesian to crystal coordinates */
815 xp[i][j] *= invbox[j][j];
816 for(k=j+1; k<npbcdim; k++) {
817 xp[i][j] += invbox[k][j]*xp[i][k];
819 } else if (rc_scaling == erscCOM) {
820 /* Subtract the center of mass */
828 if (rc_scaling == erscCOM) {
829 /* Convert the COM from Cartesian to crystal coordinates */
830 for(j=0; j<npbcdim; j++) {
831 com[j] *= invbox[j][j];
832 for(k=j+1; k<npbcdim; k++) {
833 com[j] += invbox[k][j]*com[k];
839 free_t_atoms(&dumat,TRUE);
845 static void gen_posres(gmx_mtop_t *mtop,t_molinfo *mi,
846 char *fnA, char *fnB,
847 int rc_scaling, int ePBC,
853 read_posres (mtop,mi,FALSE,fnA,rc_scaling,ePBC,com,wi);
854 if (strcmp(fnA,fnB) != 0) {
855 read_posres(mtop,mi,TRUE ,fnB,rc_scaling,ePBC,comB,wi);
859 static void set_wall_atomtype(gpp_atomtype_t at,t_gromppopts *opts,
860 t_inputrec *ir,warninp_t wi)
863 char warn_buf[STRLEN];
867 fprintf(stderr,"Searching the wall atom type(s)\n");
869 for(i=0; i<ir->nwall; i++)
871 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i],at);
872 if (ir->wall_atomtype[i] == NOTSET)
874 sprintf(warn_buf,"Specified wall atom type %s is not defined",opts->wall_atomtype[i]);
875 warning_error(wi,warn_buf);
880 static int nrdf_internal(t_atoms *atoms)
885 for(i=0; i<atoms->nr; i++) {
886 /* Vsite ptype might not be set here yet, so also check the mass */
887 if ((atoms->atom[i].ptype == eptAtom ||
888 atoms->atom[i].ptype == eptNucleus)
889 && atoms->atom[i].m > 0) {
894 case 0: nrdf = 0; break;
895 case 1: nrdf = 0; break;
896 case 2: nrdf = 1; break;
897 default: nrdf = nmass*3 - 6; break;
920 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
921 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
928 y2[i] = y2[i]*y2[i+1]+u[i];
934 interpolate1d( double xmin,
947 a = (xmin+(ix+1)*dx-x)/dx;
948 b = (x-xmin-ix*dx)/dx;
950 *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;
951 *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];
956 setup_cmap (int grid_spacing,
959 gmx_cmap_t * cmap_grid)
961 double *tmp_u,*tmp_u2,*tmp_yy,*tmp_y1,*tmp_t2,*tmp_grid;
963 int i,j,k,ii,jj,kk,idx;
965 double dx,xmin,v,v1,v2,v12;
968 snew(tmp_u,2*grid_spacing);
969 snew(tmp_u2,2*grid_spacing);
970 snew(tmp_yy,2*grid_spacing);
971 snew(tmp_y1,2*grid_spacing);
972 snew(tmp_t2,2*grid_spacing*2*grid_spacing);
973 snew(tmp_grid,2*grid_spacing*2*grid_spacing);
975 dx = 360.0/grid_spacing;
976 xmin = -180.0-dx*grid_spacing/2;
980 /* Compute an offset depending on which cmap we are using
981 * Offset will be the map number multiplied with the
982 * grid_spacing * grid_spacing * 2
984 offset = kk * grid_spacing * grid_spacing * 2;
986 for(i=0;i<2*grid_spacing;i++)
988 ii=(i+grid_spacing-grid_spacing/2)%grid_spacing;
990 for(j=0;j<2*grid_spacing;j++)
992 jj=(j+grid_spacing-grid_spacing/2)%grid_spacing;
993 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
997 for(i=0;i<2*grid_spacing;i++)
999 spline1d(dx,&(tmp_grid[2*grid_spacing*i]),2*grid_spacing,tmp_u,&(tmp_t2[2*grid_spacing*i]));
1002 for(i=grid_spacing/2;i<grid_spacing+grid_spacing/2;i++)
1004 ii = i-grid_spacing/2;
1007 for(j=grid_spacing/2;j<grid_spacing+grid_spacing/2;j++)
1009 jj = j-grid_spacing/2;
1012 for(k=0;k<2*grid_spacing;k++)
1014 interpolate1d(xmin,dx,&(tmp_grid[2*grid_spacing*k]),
1015 &(tmp_t2[2*grid_spacing*k]),psi,&tmp_yy[k],&tmp_y1[k]);
1018 spline1d(dx,tmp_yy,2*grid_spacing,tmp_u,tmp_u2);
1019 interpolate1d(xmin,dx,tmp_yy,tmp_u2,phi,&v,&v1);
1020 spline1d(dx,tmp_y1,2*grid_spacing,tmp_u,tmp_u2);
1021 interpolate1d(xmin,dx,tmp_y1,tmp_u2,phi,&v2,&v12);
1023 idx = ii*grid_spacing+jj;
1024 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1025 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1026 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1027 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1033 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1037 cmap_grid->ngrid = ngrid;
1038 cmap_grid->grid_spacing = grid_spacing;
1039 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1041 snew(cmap_grid->cmapdata,ngrid);
1043 for(i=0;i<cmap_grid->ngrid;i++)
1045 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1050 static int count_constraints(gmx_mtop_t *mtop,t_molinfo *mi,warninp_t wi)
1052 int count,count_mol,i,mb;
1053 gmx_molblock_t *molb;
1058 for(mb=0; mb<mtop->nmolblock; mb++) {
1060 molb = &mtop->molblock[mb];
1061 plist = mi[molb->type].plist;
1063 for(i=0; i<F_NRE; i++) {
1065 count_mol += 3*plist[i].nr;
1066 else if (interaction_function[i].flags & IF_CONSTRAINT)
1067 count_mol += plist[i].nr;
1070 if (count_mol > nrdf_internal(&mi[molb->type].atoms)) {
1072 "Molecule type '%s' has %d constraints.\n"
1073 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1074 *mi[molb->type].name,count_mol,
1075 nrdf_internal(&mi[molb->type].atoms));
1078 count += molb->nmol*count_mol;
1084 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1086 int i,nmiss,natoms,mt;
1088 const t_atoms *atoms;
1091 for(mt=0;mt<sys->nmoltype;mt++)
1093 atoms = &sys->moltype[mt].atoms;
1096 for(i=0;i<natoms;i++)
1098 q = atoms->atom[i].q;
1099 if ((get_atomtype_radius(atoms->atom[i].type,atype) == 0 ||
1100 get_atomtype_vol(atoms->atom[i].type,atype) == 0 ||
1101 get_atomtype_surftens(atoms->atom[i].type,atype) == 0 ||
1102 get_atomtype_gb_radius(atoms->atom[i].type,atype) == 0 ||
1103 get_atomtype_S_hct(atoms->atom[i].type,atype) == 0) &&
1106 fprintf(stderr,"\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1107 get_atomtype_name(atoms->atom[i].type,atype),q);
1115 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);
1120 static void check_gbsa_params(t_inputrec *ir,gpp_atomtype_t atype)
1124 /* If we are doing GBSA, check that we got the parameters we need
1125 * This checking is to see if there are GBSA paratmeters for all
1126 * atoms in the force field. To go around this for testing purposes
1127 * comment out the nerror++ counter temporarily
1130 for(i=0;i<get_atomtype_ntypes(atype);i++)
1132 if (get_atomtype_radius(i,atype) < 0 ||
1133 get_atomtype_vol(i,atype) < 0 ||
1134 get_atomtype_surftens(i,atype) < 0 ||
1135 get_atomtype_gb_radius(i,atype) < 0 ||
1136 get_atomtype_S_hct(i,atype) < 0)
1138 fprintf(stderr,"\nGB parameter(s) missing or negative for atom type '%s'\n",
1139 get_atomtype_name(i,atype));
1146 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);
1151 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1154 real verletbuf_drift,
1159 verletbuf_list_setup_t ls;
1162 char warn_buf[STRLEN];
1165 for(i=0; i<ir->opts.ngtc; i++)
1167 if (ir->opts.ref_t[i] < 0)
1169 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.");
1173 ref_T = max(ref_T,ir->opts.ref_t[i]);
1177 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n",verletbuf_drift,ref_T);
1179 for(i=0; i<ir->opts.ngtc; i++)
1181 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1183 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.",
1184 ir->opts.nrdf[i],ir->opts.ref_t[i],ref_T);
1185 warning_note(wi,warn_buf);
1189 /* Calculate the buffer size for simple atom vs atoms list */
1190 ls.cluster_size_i = 1;
1191 ls.cluster_size_j = 1;
1192 calc_verlet_buffer_size(mtop,det(box),ir,verletbuf_drift,
1193 &ls,&n_nonlin_vsite,&rlist_1x1);
1195 /* Set the pair-list buffer size in ir */
1196 verletbuf_get_list_setup(FALSE,&ls);
1197 calc_verlet_buffer_size(mtop,det(box),ir,verletbuf_drift,
1198 &ls,&n_nonlin_vsite,&ir->rlist);
1200 if (n_nonlin_vsite > 0)
1202 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);
1203 warning_note(wi,warn_buf);
1206 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1207 1,1,rlist_1x1,rlist_1x1-max(ir->rvdw,ir->rcoulomb));
1209 ir->rlistlong = ir->rlist;
1210 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1211 ls.cluster_size_i,ls.cluster_size_j,
1212 ir->rlist,ir->rlist-max(ir->rvdw,ir->rcoulomb));
1214 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC,box))
1216 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)));
1220 int cmain (int argc, char *argv[])
1222 static const char *desc[] = {
1223 "The gromacs preprocessor",
1224 "reads a molecular topology file, checks the validity of the",
1225 "file, expands the topology from a molecular description to an atomic",
1226 "description. The topology file contains information about",
1227 "molecule types and the number of molecules, the preprocessor",
1228 "copies each molecule as needed. ",
1229 "There is no limitation on the number of molecule types. ",
1230 "Bonds and bond-angles can be converted into constraints, separately",
1231 "for hydrogens and heavy atoms.",
1232 "Then a coordinate file is read and velocities can be generated",
1233 "from a Maxwellian distribution if requested.",
1234 "[TT]grompp[tt] also reads parameters for the [TT]mdrun[tt] ",
1235 "(eg. number of MD steps, time step, cut-off), and others such as",
1236 "NEMD parameters, which are corrected so that the net acceleration",
1238 "Eventually a binary file is produced that can serve as the sole input",
1239 "file for the MD program.[PAR]",
1241 "[TT]grompp[tt] uses the atom names from the topology file. The atom names",
1242 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1243 "warnings when they do not match the atom names in the topology.",
1244 "Note that the atom names are irrelevant for the simulation as",
1245 "only the atom types are used for generating interaction parameters.[PAR]",
1247 "[TT]grompp[tt] uses a built-in preprocessor to resolve includes, macros, ",
1248 "etc. The preprocessor supports the following keywords:[PAR]",
1249 "#ifdef VARIABLE[BR]",
1250 "#ifndef VARIABLE[BR]",
1253 "#define VARIABLE[BR]",
1254 "#undef VARIABLE[BR]"
1255 "#include \"filename\"[BR]",
1256 "#include <filename>[PAR]",
1257 "The functioning of these statements in your topology may be modulated by",
1258 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1259 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1260 "include = -I/home/john/doe[tt][BR]",
1261 "For further information a C-programming textbook may help you out.",
1262 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1263 "topology file written out so that you can verify its contents.[PAR]",
1265 /* cpp has been unnecessary for some time, hasn't it?
1266 "If your system does not have a C-preprocessor, you can still",
1267 "use [TT]grompp[tt], but you do not have access to the features ",
1268 "from the cpp. Command line options to the C-preprocessor can be given",
1269 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1272 "When using position restraints a file with restraint coordinates",
1273 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1274 "with respect to the conformation from the [TT]-c[tt] option.",
1275 "For free energy calculation the the coordinates for the B topology",
1276 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1277 "those of the A topology.[PAR]",
1279 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1280 "The last frame with coordinates and velocities will be read,",
1281 "unless the [TT]-time[tt] option is used. Only if this information",
1282 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1283 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1284 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1285 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1288 "[TT]grompp[tt] can be used to restart simulations (preserving",
1289 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1290 "However, for simply changing the number of run steps to extend",
1291 "a run, using [TT]tpbconv[tt] is more convenient than [TT]grompp[tt].",
1292 "You then supply the old checkpoint file directly to [TT]mdrun[tt]",
1293 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1294 "like output frequency, then supplying the checkpoint file to",
1295 "[TT]grompp[tt] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1296 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1298 "By default, all bonded interactions which have constant energy due to",
1299 "virtual site constructions will be removed. If this constant energy is",
1300 "not zero, this will result in a shift in the total energy. All bonded",
1301 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1302 "all constraints for distances which will be constant anyway because",
1303 "of virtual site constructions will be removed. If any constraints remain",
1304 "which involve virtual sites, a fatal error will result.[PAR]"
1306 "To verify your run input file, please take note of all warnings",
1307 "on the screen, and correct where necessary. Do also look at the contents",
1308 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1309 "the input that [TT]grompp[tt] has read. If in doubt, you can start [TT]grompp[tt]",
1310 "with the [TT]-debug[tt] option which will give you more information",
1311 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1312 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1313 "program. [TT]gmxcheck[tt] can be used to compare the contents of two",
1314 "run input files.[PAR]"
1316 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1317 "by [TT]grompp[tt] that otherwise halt output. In some cases, warnings are",
1318 "harmless, but usually they are not. The user is advised to carefully",
1319 "interpret the output messages before attempting to bypass them with",
1326 gpp_atomtype_t atype;
1328 int natoms,nvsite,comb,mt;
1332 real max_spacing,fudgeQQ;
1334 char fn[STRLEN],fnB[STRLEN];
1335 const char *mdparin;
1337 gmx_bool bNeedVel,bGenVel;
1338 gmx_bool have_atomnumber;
1340 t_params *gb_plist = NULL;
1341 gmx_genborn_t *born = NULL;
1343 gmx_bool bVerbose = FALSE;
1345 char warn_buf[STRLEN];
1348 { efMDP, NULL, NULL, ffREAD },
1349 { efMDP, "-po", "mdout", ffWRITE },
1350 { efSTX, "-c", NULL, ffREAD },
1351 { efSTX, "-r", NULL, ffOPTRD },
1352 { efSTX, "-rb", NULL, ffOPTRD },
1353 { efNDX, NULL, NULL, ffOPTRD },
1354 { efTOP, NULL, NULL, ffREAD },
1355 { efTOP, "-pp", "processed", ffOPTWR },
1356 { efTPX, "-o", NULL, ffWRITE },
1357 { efTRN, "-t", NULL, ffOPTRD },
1358 { efEDR, "-e", NULL, ffOPTRD },
1359 { efTRN, "-ref","rotref", ffOPTRW }
1361 #define NFILE asize(fnm)
1363 /* Command line options */
1364 static gmx_bool bRenum=TRUE;
1365 static gmx_bool bRmVSBds=TRUE,bZero=FALSE;
1366 static int i,maxwarn=0;
1367 static real fr_time=-1;
1369 { "-v", FALSE, etBOOL,{&bVerbose},
1370 "Be loud and noisy" },
1371 { "-time", FALSE, etREAL, {&fr_time},
1372 "Take frame at or first after this time." },
1373 { "-rmvsbds",FALSE, etBOOL, {&bRmVSBds},
1374 "Remove constant bonded interactions with virtual sites" },
1375 { "-maxwarn", FALSE, etINT, {&maxwarn},
1376 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1377 { "-zero", FALSE, etBOOL, {&bZero},
1378 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1379 { "-renum", FALSE, etBOOL, {&bRenum},
1380 "Renumber atomtypes and minimize number of atomtypes" }
1383 CopyRight(stderr,argv[0]);
1385 /* Initiate some variables */
1390 /* Parse the command line */
1391 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
1392 asize(desc),desc,0,NULL,&oenv);
1394 wi = init_warning(TRUE,maxwarn);
1396 /* PARAMETER file processing */
1397 mdparin = opt2fn("-f",NFILE,fnm);
1398 set_warning_line(wi,mdparin,-1);
1399 get_ir(mdparin,opt2fn("-po",NFILE,fnm),ir,opts,wi);
1402 fprintf(stderr,"checking input for internal consistency...\n");
1403 check_ir(mdparin,ir,opts,wi);
1405 if (ir->ld_seed == -1) {
1406 ir->ld_seed = make_seed();
1407 fprintf(stderr,"Setting the LD random seed to %d\n",ir->ld_seed);
1410 if (ir->expandedvals->lmc_seed == -1) {
1411 ir->expandedvals->lmc_seed = make_seed();
1412 fprintf(stderr,"Setting the lambda MC random seed to %d\n",ir->expandedvals->lmc_seed);
1415 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1416 bGenVel = (bNeedVel && opts->bGenVel);
1417 if (bGenVel && ir->bContinuation)
1420 "Generating velocities is inconsistent with attempting "
1421 "to continue a previous run. Choose only one of "
1422 "gen-vel = yes and continuation = yes.");
1423 warning_error(wi, warn_buf);
1429 atype = init_atomtype();
1431 pr_symtab(debug,0,"Just opened",&sys->symtab);
1433 strcpy(fn,ftp2fn(efTOP,NFILE,fnm));
1434 if (!gmx_fexist(fn))
1435 gmx_fatal(FARGS,"%s does not exist",fn);
1436 new_status(fn,opt2fn_null("-pp",NFILE,fnm),opt2fn("-c",NFILE,fnm),
1437 opts,ir,bZero,bGenVel,bVerbose,&state,
1438 atype,sys,&nmi,&mi,plist,&comb,&reppow,&fudgeQQ,
1443 pr_symtab(debug,0,"After new_status",&sys->symtab);
1445 if (ir->cutoff_scheme == ecutsVERLET)
1447 fprintf(stderr,"Removing all charge groups because cutoff-scheme=%s\n",
1448 ecutscheme_names[ir->cutoff_scheme]);
1450 /* Remove all charge groups */
1451 gmx_mtop_remove_chargegroups(sys);
1454 if (count_constraints(sys,mi,wi) && (ir->eConstrAlg == econtSHAKE)) {
1455 if (ir->eI == eiCG || ir->eI == eiLBFGS) {
1456 sprintf(warn_buf,"Can not do %s with %s, use %s",
1457 EI(ir->eI),econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1458 warning_error(wi,warn_buf);
1460 if (ir->bPeriodicMols) {
1461 sprintf(warn_buf,"Can not do periodic molecules with %s, use %s",
1462 econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1463 warning_error(wi,warn_buf);
1467 if ( EI_SD (ir->eI) && ir->etc != etcNO ) {
1468 warning_note(wi,"Temperature coupling is ignored with SD integrators.");
1471 /* If we are doing QM/MM, check that we got the atom numbers */
1472 have_atomnumber = TRUE;
1473 for (i=0; i<get_atomtype_ntypes(atype); i++) {
1474 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i,atype) >= 0);
1476 if (!have_atomnumber && ir->bQMMM)
1480 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1481 "field you are using does not contain atom numbers fields. This is an\n"
1482 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1483 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1484 "an integer just before the mass column in ffXXXnb.itp.\n"
1485 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1489 if ((ir->adress->const_wf>1) || (ir->adress->const_wf<0)) {
1490 warning_error(wi,"AdResS contant weighting function should be between 0 and 1\n\n");
1492 /** \TODO check size of ex+hy width against box size */
1495 /* Check for errors in the input now, since they might cause problems
1496 * during processing further down.
1498 check_warning_error(wi,FARGS);
1500 if (opt2bSet("-r",NFILE,fnm))
1501 sprintf(fn,"%s",opt2fn("-r",NFILE,fnm));
1503 sprintf(fn,"%s",opt2fn("-c",NFILE,fnm));
1504 if (opt2bSet("-rb",NFILE,fnm))
1505 sprintf(fnB,"%s",opt2fn("-rb",NFILE,fnm));
1509 if (nint_ftype(sys,mi,F_POSRES) > 0 || nint_ftype(sys,mi,F_FBPOSRES) > 0)
1513 fprintf(stderr,"Reading position restraint coords from %s",fn);
1514 if (strcmp(fn,fnB) == 0)
1516 fprintf(stderr,"\n");
1520 fprintf(stderr," and %s\n",fnB);
1523 gen_posres(sys,mi,fn,fnB,
1524 ir->refcoord_scaling,ir->ePBC,
1525 ir->posres_com,ir->posres_comB,
1530 /* set parameters for virtual site construction (not for vsiten) */
1531 for(mt=0; mt<sys->nmoltype; mt++) {
1533 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1535 /* now throw away all obsolete bonds, angles and dihedrals: */
1536 /* note: constraints are ALWAYS removed */
1538 for(mt=0; mt<sys->nmoltype; mt++) {
1539 clean_vsite_bondeds(mi[mt].plist,sys->moltype[mt].atoms.nr,bRmVSBds);
1543 /* If we are using CMAP, setup the pre-interpolation grid */
1546 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1547 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap,&sys->ffparams.cmap_grid);
1550 set_wall_atomtype(atype,opts,ir,wi);
1552 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1553 ntype = get_atomtype_ntypes(atype);
1556 if (ir->implicit_solvent != eisNO)
1558 /* Now we have renumbered the atom types, we can check the GBSA params */
1559 check_gbsa_params(ir,atype);
1561 /* Check that all atoms that have charge and/or LJ-parameters also have
1562 * sensible GB-parameters
1564 check_gbsa_params_charged(sys,atype);
1567 /* PELA: Copy the atomtype data to the topology atomtype list */
1568 copy_atomtype_atomtypes(atype,&(sys->atomtypes));
1571 pr_symtab(debug,0,"After renum_atype",&sys->symtab);
1574 fprintf(stderr,"converting bonded parameters...\n");
1576 ntype = get_atomtype_ntypes(atype);
1577 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1580 pr_symtab(debug,0,"After convert_params",&sys->symtab);
1582 /* set ptype to VSite for virtual sites */
1583 for(mt=0; mt<sys->nmoltype; mt++) {
1584 set_vsites_ptype(FALSE,&sys->moltype[mt]);
1587 pr_symtab(debug,0,"After virtual sites",&sys->symtab);
1589 /* Check velocity for virtual sites and shells */
1591 check_vel(sys,state.v);
1597 for(i=0; i<sys->nmoltype; i++) {
1598 check_cg_sizes(ftp2fn(efTOP,NFILE,fnm),&sys->moltype[i].cgs,wi);
1601 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1603 check_bonds_timestep(sys,ir->delta_t,wi);
1606 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1608 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.");
1611 check_warning_error(wi,FARGS);
1614 fprintf(stderr,"initialising group options...\n");
1615 do_index(mdparin,ftp2fn_null(efNDX,NFILE,fnm),
1617 bGenVel ? state.v : NULL,
1620 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1623 if (EI_DYNAMICS(ir->eI) &&
1624 !(EI_MD(ir->eI) && ir->etc==etcNO) &&
1625 inputrec2nboundeddim(ir) == 3)
1627 set_verlet_buffer(sys,ir,state.box,ir->verletbuf_drift,wi);
1631 /* Init the temperature coupling state */
1632 init_gtc_state(&state,ir->opts.ngtc,0,ir->opts.nhchainlength); /* need to add nnhpres here? */
1635 fprintf(stderr,"Checking consistency between energy and charge groups...\n");
1636 check_eg_vs_cg(sys);
1639 pr_symtab(debug,0,"After index",&sys->symtab);
1640 triple_check(mdparin,ir,sys,wi);
1641 close_symtab(&sys->symtab);
1643 pr_symtab(debug,0,"After close",&sys->symtab);
1645 /* make exclusions between QM atoms */
1647 if (ir->QMMMscheme==eQMMMschemenormal && ir->ns_type == ensSIMPLE ){
1648 gmx_fatal(FARGS,"electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1651 generate_qmexcl(sys,ir,wi);
1655 if (ftp2bSet(efTRN,NFILE,fnm)) {
1657 fprintf(stderr,"getting data from old trajectory ...\n");
1658 cont_status(ftp2fn(efTRN,NFILE,fnm),ftp2fn_null(efEDR,NFILE,fnm),
1659 bNeedVel,bGenVel,fr_time,ir,&state,sys,oenv);
1662 if (ir->ePBC==epbcXY && ir->nwall!=2)
1664 clear_rvec(state.box[ZZ]);
1667 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1669 set_warning_line(wi,mdparin,-1);
1670 check_chargegroup_radii(sys,ir,state.x,wi);
1673 if (EEL_FULL(ir->coulombtype)) {
1674 /* Calculate the optimal grid dimensions */
1675 copy_mat(state.box,box);
1676 if (ir->ePBC==epbcXY && ir->nwall==2)
1677 svmul(ir->wall_ewald_zfac,box[ZZ],box[ZZ]);
1678 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1680 /* Mark fourier_spacing as not used */
1681 ir->fourier_spacing = 0;
1683 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1685 set_warning_line(wi,mdparin,-1);
1686 warning_error(wi,"Some of the Fourier grid sizes are set, but all of them need to be set.");
1688 max_spacing = calc_grid(stdout,box,ir->fourier_spacing,
1689 &(ir->nkx),&(ir->nky),&(ir->nkz));
1692 if (ir->ePull != epullNO)
1693 set_pull_init(ir,sys,state.x,state.box,oenv,opts->pull_start);
1697 set_reference_positions(ir->rot,sys,state.x,state.box,
1698 opt2fn("-ref",NFILE,fnm),opt2bSet("-ref",NFILE,fnm),
1702 /* reset_multinr(sys); */
1704 if (EEL_PME(ir->coulombtype)) {
1705 float ratio = pme_load_estimate(sys,ir,state.box);
1706 fprintf(stderr,"Estimate for the relative computational load of the PME mesh part: %.2f\n",ratio);
1707 /* With free energy we might need to do PME both for the A and B state
1708 * charges. This will double the cost, but the optimal performance will
1709 * then probably be at a slightly larger cut-off and grid spacing.
1711 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1712 (ir->efep != efepNO && ratio > 2.0/3.0)) {
1714 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1715 "and for highly parallel simulations between 0.25 and 0.33,\n"
1716 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1717 if (ir->efep != efepNO) {
1719 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1725 char warn_buf[STRLEN];
1726 double cio = compute_io(ir,sys->natoms,&sys->groups,F_NRE,1);
1727 sprintf(warn_buf,"This run will generate roughly %.0f Mb of data",cio);
1729 set_warning_line(wi,mdparin,-1);
1730 warning_note(wi,warn_buf);
1732 printf("%s\n",warn_buf);
1736 /* MRS: eventually figure out better logic for initializing the fep
1737 values that makes declaring the lambda and declaring the state not
1738 potentially conflict if not handled correctly. */
1739 if (ir->efep != efepNO)
1741 state.fep_state = ir->fepvals->init_fep_state;
1742 for (i=0;i<efptNR;i++)
1744 /* init_lambda trumps state definitions*/
1745 if (ir->fepvals->init_lambda >= 0)
1747 state.lambda[i] = ir->fepvals->init_lambda;
1751 if (ir->fepvals->all_lambda[i] == NULL)
1753 gmx_fatal(FARGS,"Values of lambda not set for a free energy calculation!");
1757 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1764 fprintf(stderr,"writing run input file...\n");
1766 done_warning(wi,FARGS);
1768 write_tpx_state(ftp2fn(efTPX,NFILE,fnm),ir,&state,sys);