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"
84 static int rm_interactions(int ifunc,int nrmols,t_molinfo mols[])
89 /* For all the molecule types */
90 for(i=0; i<nrmols; i++) {
91 n += mols[i].plist[ifunc].nr;
92 mols[i].plist[ifunc].nr=0;
97 static int check_atom_names(const char *fn1, const char *fn2,
98 gmx_mtop_t *mtop, t_atoms *at)
100 int mb,m,i,j,nmismatch;
102 #define MAXMISMATCH 20
104 if (mtop->natoms != at->nr)
105 gmx_incons("comparing atom names");
109 for(mb=0; mb<mtop->nmolblock; mb++) {
110 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
111 for(m=0; m<mtop->molblock[mb].nmol; m++) {
112 for(j=0; j < tat->nr; j++) {
113 if (strcmp( *(tat->atomname[j]) , *(at->atomname[i]) ) != 0) {
114 if (nmismatch < MAXMISMATCH) {
116 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
117 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
118 } else if (nmismatch == MAXMISMATCH) {
119 fprintf(stderr,"(more than %d non-matching atom names)\n",MAXMISMATCH);
131 static void check_eg_vs_cg(gmx_mtop_t *mtop)
133 int astart,mb,m,cg,j,firstj;
134 unsigned char firsteg,eg;
137 /* Go through all the charge groups and make sure all their
138 * atoms are in the same energy group.
142 for(mb=0; mb<mtop->nmolblock; mb++) {
143 molt = &mtop->moltype[mtop->molblock[mb].type];
144 for(m=0; m<mtop->molblock[mb].nmol; m++) {
145 for(cg=0; cg<molt->cgs.nr;cg++) {
146 /* Get the energy group of the first atom in this charge group */
147 firstj = astart + molt->cgs.index[cg];
148 firsteg = ggrpnr(&mtop->groups,egcENER,firstj);
149 for(j=molt->cgs.index[cg]+1;j<molt->cgs.index[cg+1];j++) {
150 eg = ggrpnr(&mtop->groups,egcENER,astart+j);
152 gmx_fatal(FARGS,"atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
153 firstj+1,astart+j+1,cg+1,*molt->name);
157 astart += molt->atoms.nr;
162 static void check_cg_sizes(const char *topfn,t_block *cgs,warninp_t wi)
165 char warn_buf[STRLEN];
168 for(cg=0; cg<cgs->nr; cg++)
170 maxsize = max(maxsize,cgs->index[cg+1]-cgs->index[cg]);
173 if (maxsize > MAX_CHARGEGROUP_SIZE)
175 gmx_fatal(FARGS,"The largest charge group contains %d atoms. The maximum is %d.",maxsize,MAX_CHARGEGROUP_SIZE);
177 else if (maxsize > 10)
179 set_warning_line(wi,topfn,-1);
181 "The largest charge group contains %d atoms.\n"
182 "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"
183 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
184 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
186 warning_note(wi,warn_buf);
190 static void check_bonds_timestep(gmx_mtop_t *mtop,double dt,warninp_t wi)
192 /* This check is not intended to ensure accurate integration,
193 * rather it is to signal mistakes in the mdp settings.
194 * A common mistake is to forget to turn on constraints
195 * for MD after energy minimization with flexible bonds.
196 * This check can also detect too large time steps for flexible water
197 * models, but such errors will often be masked by the constraints
198 * mdp options, which turns flexible water into water with bond constraints,
199 * but without an angle constraint. Unfortunately such incorrect use
200 * of water models can not easily be detected without checking
201 * for specific model names.
203 * The stability limit of leap-frog or velocity verlet is 4.44 steps
204 * per oscillational period.
205 * But accurate bonds distributions are lost far before that limit.
206 * To allow relatively common schemes (although not common with Gromacs)
207 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
208 * we set the note limit to 10.
210 int min_steps_warn=5;
211 int min_steps_note=10;
214 gmx_moltype_t *moltype,*w_moltype;
216 t_ilist *ilist,*ilb,*ilc,*ils;
218 int i,a1,a2,w_a1,w_a2,j;
219 real twopi2,limit2,fc,re,m1,m2,period2,w_period2;
220 gmx_bool bFound,bWater,bWarn;
221 char warn_buf[STRLEN];
223 ip = mtop->ffparams.iparams;
225 twopi2 = sqr(2*M_PI);
227 limit2 = sqr(min_steps_note*dt);
233 for(molt=0; molt<mtop->nmoltype; molt++)
235 moltype = &mtop->moltype[molt];
236 atom = moltype->atoms.atom;
237 ilist = moltype->ilist;
238 ilc = &ilist[F_CONSTR];
239 ils = &ilist[F_SETTLE];
240 for(ftype=0; ftype<F_NRE; ftype++)
242 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
248 for(i=0; i<ilb->nr; i+=3)
250 fc = ip[ilb->iatoms[i]].harmonic.krA;
251 re = ip[ilb->iatoms[i]].harmonic.rA;
252 if (ftype == F_G96BONDS)
254 /* Convert squared sqaure fc to harmonic fc */
257 a1 = ilb->iatoms[i+1];
258 a2 = ilb->iatoms[i+2];
261 if (fc > 0 && m1 > 0 && m2 > 0)
263 period2 = twopi2*m1*m2/((m1 + m2)*fc);
267 period2 = GMX_FLOAT_MAX;
271 fprintf(debug,"fc %g m1 %g m2 %g period %g\n",
272 fc,m1,m2,sqrt(period2));
274 if (period2 < limit2)
277 for(j=0; j<ilc->nr; j+=3)
279 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
280 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
285 for(j=0; j<ils->nr; j+=4)
287 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
288 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
294 (w_moltype == NULL || period2 < w_period2))
306 if (w_moltype != NULL)
308 bWarn = (w_period2 < sqr(min_steps_warn*dt));
309 /* A check that would recognize most water models */
310 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
311 w_moltype->atoms.nr <= 5);
312 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"
315 w_a1+1,*w_moltype->atoms.atomname[w_a1],
316 w_a2+1,*w_moltype->atoms.atomname[w_a2],
317 sqrt(w_period2),bWarn ? min_steps_warn : min_steps_note,dt,
319 "Maybe you asked for fexible water." :
320 "Maybe you forgot to change the constraints mdp option.");
323 warning(wi,warn_buf);
327 warning_note(wi,warn_buf);
332 static void check_vel(gmx_mtop_t *mtop,rvec v[])
334 gmx_mtop_atomloop_all_t aloop;
338 aloop = gmx_mtop_atomloop_all_init(mtop);
339 while (gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
340 if (atom->ptype == eptShell ||
341 atom->ptype == eptBond ||
342 atom->ptype == eptVSite) {
348 static gmx_bool nint_ftype(gmx_mtop_t *mtop,t_molinfo *mi,int ftype)
353 for(mb=0; mb<mtop->nmolblock; mb++) {
354 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
360 /* This routine reorders the molecule type array
361 * in the order of use in the molblocks,
362 * unused molecule types are deleted.
364 static void renumber_moltypes(gmx_mtop_t *sys,
365 int *nmolinfo,t_molinfo **molinfo)
371 snew(order,*nmolinfo);
373 for(mb=0; mb<sys->nmolblock; mb++) {
374 for(i=0; i<norder; i++) {
375 if (order[i] == sys->molblock[mb].type) {
380 /* This type did not occur yet, add it */
381 order[norder] = sys->molblock[mb].type;
382 /* Renumber the moltype in the topology */
385 sys->molblock[mb].type = i;
388 /* We still need to reorder the molinfo structs */
390 for(mi=0; mi<*nmolinfo; mi++) {
391 for(i=0; i<norder; i++) {
392 if (order[i] == mi) {
397 done_mi(&(*molinfo)[mi]);
399 minew[i] = (*molinfo)[mi];
408 static void molinfo2mtop(int nmi,t_molinfo *mi,gmx_mtop_t *mtop)
413 mtop->nmoltype = nmi;
414 snew(mtop->moltype,nmi);
415 for(m=0; m<nmi; m++) {
416 molt = &mtop->moltype[m];
417 molt->name = mi[m].name;
418 molt->atoms = mi[m].atoms;
419 /* ilists are copied later */
420 molt->cgs = mi[m].cgs;
421 molt->excls = mi[m].excls;
426 new_status(const char *topfile,const char *topppfile,const char *confin,
427 t_gromppopts *opts,t_inputrec *ir,gmx_bool bZero,
428 gmx_bool bGenVel,gmx_bool bVerbose,t_state *state,
429 gpp_atomtype_t atype,gmx_mtop_t *sys,
430 int *nmi,t_molinfo **mi,t_params plist[],
431 int *comb,double *reppow,real *fudgeQQ,
435 t_molinfo *molinfo=NULL;
437 gmx_molblock_t *molblock,*molbs;
439 int mb,i,nrmols,nmismatch;
442 char warn_buf[STRLEN];
446 /* Set gmx_boolean for GB */
447 if(ir->implicit_solvent)
450 /* TOPOLOGY processing */
451 sys->name = do_top(bVerbose,topfile,topppfile,opts,bZero,&(sys->symtab),
452 plist,comb,reppow,fudgeQQ,
453 atype,&nrmols,&molinfo,ir,
454 &nmolblock,&molblock,bGB,
458 snew(sys->molblock,nmolblock);
461 for(mb=0; mb<nmolblock; mb++) {
462 if (sys->nmolblock > 0 &&
463 molblock[mb].type == sys->molblock[sys->nmolblock-1].type) {
464 /* Merge consecutive blocks with the same molecule type */
465 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
466 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
467 } else if (molblock[mb].nmol > 0) {
468 /* Add a new molblock to the topology */
469 molbs = &sys->molblock[sys->nmolblock];
470 *molbs = molblock[mb];
471 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
472 molbs->nposres_xA = 0;
473 molbs->nposres_xB = 0;
474 sys->natoms += molbs->nmol*molbs->natoms_mol;
478 if (sys->nmolblock == 0) {
479 gmx_fatal(FARGS,"No molecules were defined in the system");
482 renumber_moltypes(sys,&nrmols,&molinfo);
485 convert_harmonics(nrmols,molinfo,atype);
487 if (ir->eDisre == edrNone) {
488 i = rm_interactions(F_DISRES,nrmols,molinfo);
490 set_warning_line(wi,"unknown",-1);
491 sprintf(warn_buf,"disre = no, removed %d distance restraints",i);
492 warning_note(wi,warn_buf);
495 if (opts->bOrire == FALSE) {
496 i = rm_interactions(F_ORIRES,nrmols,molinfo);
498 set_warning_line(wi,"unknown",-1);
499 sprintf(warn_buf,"orire = no, removed %d orientation restraints",i);
500 warning_note(wi,warn_buf);
504 /* Copy structures from msys to sys */
505 molinfo2mtop(nrmols,molinfo,sys);
507 gmx_mtop_finalize(sys);
509 /* COORDINATE file processing */
511 fprintf(stderr,"processing coordinates...\n");
513 get_stx_coordnum(confin,&state->natoms);
514 if (state->natoms != sys->natoms)
515 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
516 " does not match topology (%s, %d)",
517 confin,state->natoms,topfile,sys->natoms);
519 /* make space for coordinates and velocities */
522 init_t_atoms(confat,state->natoms,FALSE);
523 init_state(state,state->natoms,0,0,0,0);
524 read_stx_conf(confin,title,confat,state->x,state->v,NULL,state->box);
525 /* This call fixes the box shape for runs with pressure scaling */
526 set_box_rel(ir,state);
528 nmismatch = check_atom_names(topfile, confin, sys, confat);
529 free_t_atoms(confat,TRUE);
533 sprintf(buf,"%d non-matching atom name%s\n"
534 "atom names from %s will be used\n"
535 "atom names from %s will be ignored\n",
536 nmismatch,(nmismatch == 1) ? "" : "s",topfile,confin);
540 fprintf(stderr,"double-checking input for internal consistency...\n");
541 double_check(ir,state->box,nint_ftype(sys,molinfo,F_CONSTR),wi);
546 gmx_mtop_atomloop_all_t aloop;
549 snew(mass,state->natoms);
550 aloop = gmx_mtop_atomloop_all_init(sys);
551 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
555 if (opts->seed == -1) {
556 opts->seed = make_seed();
557 fprintf(stderr,"Setting gen_seed to %d\n",opts->seed);
559 maxwell_speed(opts->tempi,opts->seed,sys,state->v);
561 stop_cm(stdout,state->natoms,mass,state->x,state->v);
569 static void copy_state(const char *slog,t_trxframe *fr,
570 gmx_bool bReadVel,t_state *state,
575 if (fr->not_ok & FRAME_NOT_OK)
577 gmx_fatal(FARGS,"Can not start from an incomplete frame");
581 gmx_fatal(FARGS,"Did not find a frame with coordinates in file %s",
585 for(i=0; i<state->natoms; i++)
587 copy_rvec(fr->x[i],state->x[i]);
593 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
595 for(i=0; i<state->natoms; i++)
597 copy_rvec(fr->v[i],state->v[i]);
602 copy_mat(fr->box,state->box);
605 *use_time = fr->time;
608 static void cont_status(const char *slog,const char *ener,
609 gmx_bool bNeedVel,gmx_bool bGenVel, real fr_time,
610 t_inputrec *ir,t_state *state,
612 const output_env_t oenv)
613 /* If fr_time == -1 read the last frame available which is complete */
621 bReadVel = (bNeedVel && !bGenVel);
624 "Reading Coordinates%s and Box size from old trajectory\n",
625 bReadVel ? ", Velocities" : "");
628 fprintf(stderr,"Will read whole trajectory\n");
632 fprintf(stderr,"Will read till time %g\n",fr_time);
638 fprintf(stderr,"Velocities generated: "
639 "ignoring velocities in input trajectory\n");
641 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
645 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X | TRX_NEED_V);
651 "WARNING: Did not find a frame with velocities in file %s,\n"
652 " all velocities will be set to zero!\n\n",slog);
653 for(i=0; i<sys->natoms; i++)
655 clear_rvec(state->v[i]);
658 /* Search for a frame without velocities */
660 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
664 state->natoms = fr.natoms;
666 if (sys->natoms != state->natoms)
668 gmx_fatal(FARGS,"Number of atoms in Topology "
669 "is not the same as in Trajectory");
671 copy_state(slog,&fr,bReadVel,state,&use_time);
673 /* Find the appropriate frame */
674 while ((fr_time == -1 || fr.time < fr_time) &&
675 read_next_frame(oenv,fp,&fr))
677 copy_state(slog,&fr,bReadVel,state,&use_time);
682 /* Set the relative box lengths for preserving the box shape.
683 * Note that this call can lead to differences in the last bit
684 * with respect to using tpbconv to create a [TT].tpx[tt] file.
686 set_box_rel(ir,state);
688 fprintf(stderr,"Using frame at t = %g ps\n",use_time);
689 fprintf(stderr,"Starting time for run is %g ps\n",ir->init_t);
691 if ((ir->epc != epcNO || ir->etc ==etcNOSEHOOVER) && ener)
693 get_enx_state(ener,use_time,&sys->groups,ir,state);
694 preserve_box_shape(ir,state->box_rel,state->boxv);
698 static void read_posres(gmx_mtop_t *mtop,t_molinfo *molinfo,gmx_bool bTopB,
700 int rc_scaling, int ePBC,
704 gmx_bool bFirst = TRUE;
710 int natoms,npbcdim=0;
711 char warn_buf[STRLEN],title[STRLEN];
712 int a,i,ai,j,k,mb,nat_molb;
713 gmx_molblock_t *molb;
717 get_stx_coordnum(fn,&natoms);
718 if (natoms != mtop->natoms) {
719 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);
720 warning(wi,warn_buf);
724 init_t_atoms(&dumat,natoms,FALSE);
725 read_stx_conf(fn,title,&dumat,x,v,NULL,box);
727 npbcdim = ePBC2npbcdim(ePBC);
729 if (rc_scaling != erscNO) {
730 copy_mat(box,invbox);
731 for(j=npbcdim; j<DIM; j++) {
732 clear_rvec(invbox[j]);
735 m_inv_ur0(invbox,invbox);
738 /* Copy the reference coordinates to mtop */
742 for(mb=0; mb<mtop->nmolblock; mb++) {
743 molb = &mtop->molblock[mb];
744 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
745 pr = &(molinfo[molb->type].plist[F_POSRES]);
747 atom = mtop->moltype[molb->type].atoms.atom;
748 for(i=0; (i<pr->nr); i++) {
751 gmx_fatal(FARGS,"Position restraint atom index (%d) in moltype '%s' is larger than number of atoms in %s (%d).\n",
752 ai+1,*molinfo[molb->type].name,fn,natoms);
754 if (rc_scaling == erscCOM) {
755 /* Determine the center of mass of the posres reference coordinates */
756 for(j=0; j<npbcdim; j++) {
757 sum[j] += atom[ai].m*x[a+ai][j];
759 totmass += atom[ai].m;
763 molb->nposres_xA = nat_molb;
764 snew(molb->posres_xA,molb->nposres_xA);
765 for(i=0; i<nat_molb; i++) {
766 copy_rvec(x[a+i],molb->posres_xA[i]);
769 molb->nposres_xB = nat_molb;
770 snew(molb->posres_xB,molb->nposres_xB);
771 for(i=0; i<nat_molb; i++) {
772 copy_rvec(x[a+i],molb->posres_xB[i]);
778 if (rc_scaling == erscCOM) {
780 gmx_fatal(FARGS,"The total mass of the position restraint atoms is 0");
781 for(j=0; j<npbcdim; j++)
782 com[j] = sum[j]/totmass;
783 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]);
786 if (rc_scaling != erscNO) {
787 for(mb=0; mb<mtop->nmolblock; mb++) {
788 molb = &mtop->molblock[mb];
789 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
790 if (molb->nposres_xA > 0 || molb->nposres_xB > 0) {
791 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
792 for(i=0; i<nat_molb; i++) {
793 for(j=0; j<npbcdim; j++) {
794 if (rc_scaling == erscALL) {
795 /* Convert from Cartesian to crystal coordinates */
796 xp[i][j] *= invbox[j][j];
797 for(k=j+1; k<npbcdim; k++) {
798 xp[i][j] += invbox[k][j]*xp[i][k];
800 } else if (rc_scaling == erscCOM) {
801 /* Subtract the center of mass */
809 if (rc_scaling == erscCOM) {
810 /* Convert the COM from Cartesian to crystal coordinates */
811 for(j=0; j<npbcdim; j++) {
812 com[j] *= invbox[j][j];
813 for(k=j+1; k<npbcdim; k++) {
814 com[j] += invbox[k][j]*com[k];
820 free_t_atoms(&dumat,TRUE);
825 static void gen_posres(gmx_mtop_t *mtop,t_molinfo *mi,
826 char *fnA, char *fnB,
827 int rc_scaling, int ePBC,
833 read_posres (mtop,mi,FALSE,fnA,rc_scaling,ePBC,com,wi);
834 if (strcmp(fnA,fnB) != 0) {
835 read_posres(mtop,mi,TRUE ,fnB,rc_scaling,ePBC,comB,wi);
839 static void set_wall_atomtype(gpp_atomtype_t at,t_gromppopts *opts,
840 t_inputrec *ir,warninp_t wi)
843 char warn_buf[STRLEN];
847 fprintf(stderr,"Searching the wall atom type(s)\n");
849 for(i=0; i<ir->nwall; i++)
851 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i],at);
852 if (ir->wall_atomtype[i] == NOTSET)
854 sprintf(warn_buf,"Specified wall atom type %s is not defined",opts->wall_atomtype[i]);
855 warning_error(wi,warn_buf);
860 static int nrdf_internal(t_atoms *atoms)
865 for(i=0; i<atoms->nr; i++) {
866 /* Vsite ptype might not be set here yet, so also check the mass */
867 if ((atoms->atom[i].ptype == eptAtom ||
868 atoms->atom[i].ptype == eptNucleus)
869 && atoms->atom[i].m > 0) {
874 case 0: nrdf = 0; break;
875 case 1: nrdf = 0; break;
876 case 2: nrdf = 1; break;
877 default: nrdf = nmass*3 - 6; break;
900 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
901 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
908 y2[i] = y2[i]*y2[i+1]+u[i];
914 interpolate1d( double xmin,
927 a = (xmin+(ix+1)*dx-x)/dx;
928 b = (x-xmin-ix*dx)/dx;
930 *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;
931 *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];
936 setup_cmap (int grid_spacing,
939 gmx_cmap_t * cmap_grid)
941 double *tmp_u,*tmp_u2,*tmp_yy,*tmp_y1,*tmp_t2,*tmp_grid;
943 int i,j,k,ii,jj,kk,idx;
945 double dx,xmin,v,v1,v2,v12;
948 snew(tmp_u,2*grid_spacing);
949 snew(tmp_u2,2*grid_spacing);
950 snew(tmp_yy,2*grid_spacing);
951 snew(tmp_y1,2*grid_spacing);
952 snew(tmp_t2,2*grid_spacing*2*grid_spacing);
953 snew(tmp_grid,2*grid_spacing*2*grid_spacing);
955 dx = 360.0/grid_spacing;
956 xmin = -180.0-dx*grid_spacing/2;
960 /* Compute an offset depending on which cmap we are using
961 * Offset will be the map number multiplied with the
962 * grid_spacing * grid_spacing * 2
964 offset = kk * grid_spacing * grid_spacing * 2;
966 for(i=0;i<2*grid_spacing;i++)
968 ii=(i+grid_spacing-grid_spacing/2)%grid_spacing;
970 for(j=0;j<2*grid_spacing;j++)
972 jj=(j+grid_spacing-grid_spacing/2)%grid_spacing;
973 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
977 for(i=0;i<2*grid_spacing;i++)
979 spline1d(dx,&(tmp_grid[2*grid_spacing*i]),2*grid_spacing,tmp_u,&(tmp_t2[2*grid_spacing*i]));
982 for(i=grid_spacing/2;i<grid_spacing+grid_spacing/2;i++)
984 ii = i-grid_spacing/2;
987 for(j=grid_spacing/2;j<grid_spacing+grid_spacing/2;j++)
989 jj = j-grid_spacing/2;
992 for(k=0;k<2*grid_spacing;k++)
994 interpolate1d(xmin,dx,&(tmp_grid[2*grid_spacing*k]),
995 &(tmp_t2[2*grid_spacing*k]),psi,&tmp_yy[k],&tmp_y1[k]);
998 spline1d(dx,tmp_yy,2*grid_spacing,tmp_u,tmp_u2);
999 interpolate1d(xmin,dx,tmp_yy,tmp_u2,phi,&v,&v1);
1000 spline1d(dx,tmp_y1,2*grid_spacing,tmp_u,tmp_u2);
1001 interpolate1d(xmin,dx,tmp_y1,tmp_u2,phi,&v2,&v12);
1003 idx = ii*grid_spacing+jj;
1004 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1005 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1006 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1007 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1013 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1017 cmap_grid->ngrid = ngrid;
1018 cmap_grid->grid_spacing = grid_spacing;
1019 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1021 snew(cmap_grid->cmapdata,ngrid);
1023 for(i=0;i<cmap_grid->ngrid;i++)
1025 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1030 static int count_constraints(gmx_mtop_t *mtop,t_molinfo *mi,warninp_t wi)
1032 int count,count_mol,i,mb;
1033 gmx_molblock_t *molb;
1038 for(mb=0; mb<mtop->nmolblock; mb++) {
1040 molb = &mtop->molblock[mb];
1041 plist = mi[molb->type].plist;
1043 for(i=0; i<F_NRE; i++) {
1045 count_mol += 3*plist[i].nr;
1046 else if (interaction_function[i].flags & IF_CONSTRAINT)
1047 count_mol += plist[i].nr;
1050 if (count_mol > nrdf_internal(&mi[molb->type].atoms)) {
1052 "Molecule type '%s' has %d constraints.\n"
1053 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1054 *mi[molb->type].name,count_mol,
1055 nrdf_internal(&mi[molb->type].atoms));
1058 count += molb->nmol*count_mol;
1064 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1066 int i,nmiss,natoms,mt;
1068 const t_atoms *atoms;
1071 for(mt=0;mt<sys->nmoltype;mt++)
1073 atoms = &sys->moltype[mt].atoms;
1076 for(i=0;i<natoms;i++)
1078 q = atoms->atom[i].q;
1079 if ((get_atomtype_radius(atoms->atom[i].type,atype) == 0 ||
1080 get_atomtype_vol(atoms->atom[i].type,atype) == 0 ||
1081 get_atomtype_surftens(atoms->atom[i].type,atype) == 0 ||
1082 get_atomtype_gb_radius(atoms->atom[i].type,atype) == 0 ||
1083 get_atomtype_S_hct(atoms->atom[i].type,atype) == 0) &&
1086 fprintf(stderr,"\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1087 get_atomtype_name(atoms->atom[i].type,atype),q);
1095 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);
1100 static void check_gbsa_params(t_inputrec *ir,gpp_atomtype_t atype)
1104 /* If we are doing GBSA, check that we got the parameters we need
1105 * This checking is to see if there are GBSA paratmeters for all
1106 * atoms in the force field. To go around this for testing purposes
1107 * comment out the nerror++ counter temporarily
1110 for(i=0;i<get_atomtype_ntypes(atype);i++)
1112 if (get_atomtype_radius(i,atype) < 0 ||
1113 get_atomtype_vol(i,atype) < 0 ||
1114 get_atomtype_surftens(i,atype) < 0 ||
1115 get_atomtype_gb_radius(i,atype) < 0 ||
1116 get_atomtype_S_hct(i,atype) < 0)
1118 fprintf(stderr,"\nGB parameter(s) missing or negative for atom type '%s'\n",
1119 get_atomtype_name(i,atype));
1126 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);
1131 static void check_settle(gmx_mtop_t *sys)
1135 nra = interaction_function[F_SETTLE].nratoms;
1136 for(i=0; (i<sys->nmoltype); i++)
1138 for(j=0; (j<sys->moltype[i].ilist[F_SETTLE].nr); j+=nra+1)
1140 cgj1 = sys->moltype[i].cgs.index[j+1];
1142 gmx_fatal(FARGS,"For SETTLE you need to have all atoms involved in one charge group. Please fix your topology.");
1147 int main (int argc, char *argv[])
1149 static const char *desc[] = {
1150 "The gromacs preprocessor",
1151 "reads a molecular topology file, checks the validity of the",
1152 "file, expands the topology from a molecular description to an atomic",
1153 "description. The topology file contains information about",
1154 "molecule types and the number of molecules, the preprocessor",
1155 "copies each molecule as needed. ",
1156 "There is no limitation on the number of molecule types. ",
1157 "Bonds and bond-angles can be converted into constraints, separately",
1158 "for hydrogens and heavy atoms.",
1159 "Then a coordinate file is read and velocities can be generated",
1160 "from a Maxwellian distribution if requested.",
1161 "[TT]grompp[tt] also reads parameters for the [TT]mdrun[tt] ",
1162 "(eg. number of MD steps, time step, cut-off), and others such as",
1163 "NEMD parameters, which are corrected so that the net acceleration",
1165 "Eventually a binary file is produced that can serve as the sole input",
1166 "file for the MD program.[PAR]",
1168 "[TT]grompp[tt] uses the atom names from the topology file. The atom names",
1169 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1170 "warnings when they do not match the atom names in the topology.",
1171 "Note that the atom names are irrelevant for the simulation as",
1172 "only the atom types are used for generating interaction parameters.[PAR]",
1174 "[TT]grompp[tt] uses a built-in preprocessor to resolve includes, macros, ",
1175 "etc. The preprocessor supports the following keywords:[PAR]",
1176 "#ifdef VARIABLE[BR]",
1177 "#ifndef VARIABLE[BR]",
1180 "#define VARIABLE[BR]",
1181 "#undef VARIABLE[BR]"
1182 "#include \"filename\"[BR]",
1183 "#include <filename>[PAR]",
1184 "The functioning of these statements in your topology may be modulated by",
1185 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1186 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1187 "include = -I/home/john/doe[tt][BR]",
1188 "For further information a C-programming textbook may help you out.",
1189 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1190 "topology file written out so that you can verify its contents.[PAR]",
1192 /* cpp has been unnecessary for some time, hasn't it?
1193 "If your system does not have a C-preprocessor, you can still",
1194 "use [TT]grompp[tt], but you do not have access to the features ",
1195 "from the cpp. Command line options to the C-preprocessor can be given",
1196 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1199 "When using position restraints a file with restraint coordinates",
1200 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1201 "with respect to the conformation from the [TT]-c[tt] option.",
1202 "For free energy calculation the the coordinates for the B topology",
1203 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1204 "those of the A topology.[PAR]",
1206 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1207 "The last frame with coordinates and velocities will be read,",
1208 "unless the [TT]-time[tt] option is used. Only if this information",
1209 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1210 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1211 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1212 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1215 "[TT]grompp[tt] can be used to restart simulations (preserving",
1216 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1217 "However, for simply changing the number of run steps to extend",
1218 "a run, using [TT]tpbconv[tt] is more convenient than [TT]grompp[tt].",
1219 "You then supply the old checkpoint file directly to [TT]mdrun[tt]",
1220 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1221 "like output frequency, then supplying the checkpoint file to",
1222 "[TT]grompp[tt] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1223 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1225 "By default, all bonded interactions which have constant energy due to",
1226 "virtual site constructions will be removed. If this constant energy is",
1227 "not zero, this will result in a shift in the total energy. All bonded",
1228 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1229 "all constraints for distances which will be constant anyway because",
1230 "of virtual site constructions will be removed. If any constraints remain",
1231 "which involve virtual sites, a fatal error will result.[PAR]"
1233 "To verify your run input file, please take note of all warnings",
1234 "on the screen, and correct where necessary. Do also look at the contents",
1235 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1236 "the input that [TT]grompp[tt] has read. If in doubt, you can start [TT]grompp[tt]",
1237 "with the [TT]-debug[tt] option which will give you more information",
1238 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1239 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1240 "program. [TT]gmxcheck[tt] can be used to compare the contents of two",
1241 "run input files.[PAR]"
1243 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1244 "by [TT]grompp[tt] that otherwise halt output. In some cases, warnings are",
1245 "harmless, but usually they are not. The user is advised to carefully",
1246 "interpret the output messages before attempting to bypass them with",
1253 gpp_atomtype_t atype;
1255 int natoms,nvsite,comb,mt;
1259 real max_spacing,fudgeQQ;
1261 char fn[STRLEN],fnB[STRLEN];
1262 const char *mdparin;
1264 gmx_bool bNeedVel,bGenVel;
1265 gmx_bool have_atomnumber;
1267 t_params *gb_plist = NULL;
1268 gmx_genborn_t *born = NULL;
1270 gmx_bool bVerbose = FALSE;
1272 char warn_buf[STRLEN];
1275 { efMDP, NULL, NULL, ffREAD },
1276 { efMDP, "-po", "mdout", ffWRITE },
1277 { efSTX, "-c", NULL, ffREAD },
1278 { efSTX, "-r", NULL, ffOPTRD },
1279 { efSTX, "-rb", NULL, ffOPTRD },
1280 { efNDX, NULL, NULL, ffOPTRD },
1281 { efTOP, NULL, NULL, ffREAD },
1282 { efTOP, "-pp", "processed", ffOPTWR },
1283 { efTPX, "-o", NULL, ffWRITE },
1284 { efTRN, "-t", NULL, ffOPTRD },
1285 { efEDR, "-e", NULL, ffOPTRD },
1286 { efTRN, "-ref","rotref", ffOPTRW }
1288 #define NFILE asize(fnm)
1290 /* Command line options */
1291 static gmx_bool bRenum=TRUE;
1292 static gmx_bool bRmVSBds=TRUE,bZero=FALSE;
1293 static int i,maxwarn=0;
1294 static real fr_time=-1;
1296 { "-v", FALSE, etBOOL,{&bVerbose},
1297 "Be loud and noisy" },
1298 { "-time", FALSE, etREAL, {&fr_time},
1299 "Take frame at or first after this time." },
1300 { "-rmvsbds",FALSE, etBOOL, {&bRmVSBds},
1301 "Remove constant bonded interactions with virtual sites" },
1302 { "-maxwarn", FALSE, etINT, {&maxwarn},
1303 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1304 { "-zero", FALSE, etBOOL, {&bZero},
1305 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1306 { "-renum", FALSE, etBOOL, {&bRenum},
1307 "Renumber atomtypes and minimize number of atomtypes" }
1310 CopyRight(stderr,argv[0]);
1312 /* Initiate some variables */
1317 /* Parse the command line */
1318 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
1319 asize(desc),desc,0,NULL,&oenv);
1321 wi = init_warning(TRUE,maxwarn);
1323 /* PARAMETER file processing */
1324 mdparin = opt2fn("-f",NFILE,fnm);
1325 set_warning_line(wi,mdparin,-1);
1326 get_ir(mdparin,opt2fn("-po",NFILE,fnm),ir,opts,wi);
1329 fprintf(stderr,"checking input for internal consistency...\n");
1330 check_ir(mdparin,ir,opts,wi);
1332 if (ir->ld_seed == -1) {
1333 ir->ld_seed = make_seed();
1334 fprintf(stderr,"Setting the LD random seed to %d\n",ir->ld_seed);
1337 if (ir->expandedvals->lmc_seed == -1) {
1338 ir->expandedvals->lmc_seed = make_seed();
1339 fprintf(stderr,"Setting the lambda MC random seed to %d\n",ir->expandedvals->lmc_seed);
1342 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1343 bGenVel = (bNeedVel && opts->bGenVel);
1348 atype = init_atomtype();
1350 pr_symtab(debug,0,"Just opened",&sys->symtab);
1352 strcpy(fn,ftp2fn(efTOP,NFILE,fnm));
1353 if (!gmx_fexist(fn))
1354 gmx_fatal(FARGS,"%s does not exist",fn);
1355 new_status(fn,opt2fn_null("-pp",NFILE,fnm),opt2fn("-c",NFILE,fnm),
1356 opts,ir,bZero,bGenVel,bVerbose,&state,
1357 atype,sys,&nmi,&mi,plist,&comb,&reppow,&fudgeQQ,
1362 pr_symtab(debug,0,"After new_status",&sys->symtab);
1364 if (count_constraints(sys,mi,wi) && (ir->eConstrAlg == econtSHAKE)) {
1365 if (ir->eI == eiCG || ir->eI == eiLBFGS) {
1366 sprintf(warn_buf,"Can not do %s with %s, use %s",
1367 EI(ir->eI),econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1368 warning_error(wi,warn_buf);
1370 if (ir->bPeriodicMols) {
1371 sprintf(warn_buf,"Can not do periodic molecules with %s, use %s",
1372 econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1373 warning_error(wi,warn_buf);
1377 if ( EI_SD (ir->eI) && ir->etc != etcNO ) {
1378 warning_note(wi,"Temperature coupling is ignored with SD integrators.");
1381 /* If we are doing QM/MM, check that we got the atom numbers */
1382 have_atomnumber = TRUE;
1383 for (i=0; i<get_atomtype_ntypes(atype); i++) {
1384 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i,atype) >= 0);
1386 if (!have_atomnumber && ir->bQMMM)
1390 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1391 "field you are using does not contain atom numbers fields. This is an\n"
1392 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1393 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1394 "an integer just before the mass column in ffXXXnb.itp.\n"
1395 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1399 if ((ir->adress->const_wf>1) || (ir->adress->const_wf<0)) {
1400 warning_error(wi,"AdResS contant weighting function should be between 0 and 1\n\n");
1402 /** \TODO check size of ex+hy width against box size */
1405 /* Check for errors in the input now, since they might cause problems
1406 * during processing further down.
1408 check_warning_error(wi,FARGS);
1410 if (opt2bSet("-r",NFILE,fnm))
1411 sprintf(fn,"%s",opt2fn("-r",NFILE,fnm));
1413 sprintf(fn,"%s",opt2fn("-c",NFILE,fnm));
1414 if (opt2bSet("-rb",NFILE,fnm))
1415 sprintf(fnB,"%s",opt2fn("-rb",NFILE,fnm));
1419 if (nint_ftype(sys,mi,F_POSRES) > 0)
1423 fprintf(stderr,"Reading position restraint coords from %s",fn);
1424 if (strcmp(fn,fnB) == 0)
1426 fprintf(stderr,"\n");
1430 fprintf(stderr," and %s\n",fnB);
1433 gen_posres(sys,mi,fn,fnB,
1434 ir->refcoord_scaling,ir->ePBC,
1435 ir->posres_com,ir->posres_comB,
1440 /* set parameters for virtual site construction (not for vsiten) */
1441 for(mt=0; mt<sys->nmoltype; mt++) {
1443 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1445 /* now throw away all obsolete bonds, angles and dihedrals: */
1446 /* note: constraints are ALWAYS removed */
1448 for(mt=0; mt<sys->nmoltype; mt++) {
1449 clean_vsite_bondeds(mi[mt].plist,sys->moltype[mt].atoms.nr,bRmVSBds);
1453 /* If we are using CMAP, setup the pre-interpolation grid */
1456 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1457 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap,&sys->ffparams.cmap_grid);
1460 set_wall_atomtype(atype,opts,ir,wi);
1462 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1463 ntype = get_atomtype_ntypes(atype);
1466 if (ir->implicit_solvent != eisNO)
1468 /* Now we have renumbered the atom types, we can check the GBSA params */
1469 check_gbsa_params(ir,atype);
1471 /* Check that all atoms that have charge and/or LJ-parameters also have
1472 * sensible GB-parameters
1474 check_gbsa_params_charged(sys,atype);
1477 /* PELA: Copy the atomtype data to the topology atomtype list */
1478 copy_atomtype_atomtypes(atype,&(sys->atomtypes));
1481 pr_symtab(debug,0,"After renum_atype",&sys->symtab);
1484 fprintf(stderr,"converting bonded parameters...\n");
1486 ntype = get_atomtype_ntypes(atype);
1487 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1490 pr_symtab(debug,0,"After convert_params",&sys->symtab);
1492 /* set ptype to VSite for virtual sites */
1493 for(mt=0; mt<sys->nmoltype; mt++) {
1494 set_vsites_ptype(FALSE,&sys->moltype[mt]);
1497 pr_symtab(debug,0,"After virtual sites",&sys->symtab);
1499 /* Check velocity for virtual sites and shells */
1501 check_vel(sys,state.v);
1504 /* check for charge groups in settles */
1510 for(i=0; i<sys->nmoltype; i++) {
1511 check_cg_sizes(ftp2fn(efTOP,NFILE,fnm),&sys->moltype[i].cgs,wi);
1514 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1516 check_bonds_timestep(sys,ir->delta_t,wi);
1519 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1521 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.");
1524 check_warning_error(wi,FARGS);
1527 fprintf(stderr,"initialising group options...\n");
1528 do_index(mdparin,ftp2fn_null(efNDX,NFILE,fnm),
1530 bGenVel ? state.v : NULL,
1533 /* Init the temperature coupling state */
1534 init_gtc_state(&state,ir->opts.ngtc,0,ir->opts.nhchainlength); /* need to add nnhpres here? */
1537 fprintf(stderr,"Checking consistency between energy and charge groups...\n");
1538 check_eg_vs_cg(sys);
1541 pr_symtab(debug,0,"After index",&sys->symtab);
1542 triple_check(mdparin,ir,sys,wi);
1543 close_symtab(&sys->symtab);
1545 pr_symtab(debug,0,"After close",&sys->symtab);
1547 /* make exclusions between QM atoms */
1549 if (ir->QMMMscheme==eQMMMschemenormal && ir->ns_type == ensSIMPLE ){
1550 gmx_fatal(FARGS,"electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1553 generate_qmexcl(sys,ir,wi);
1557 if (ftp2bSet(efTRN,NFILE,fnm)) {
1559 fprintf(stderr,"getting data from old trajectory ...\n");
1560 cont_status(ftp2fn(efTRN,NFILE,fnm),ftp2fn_null(efEDR,NFILE,fnm),
1561 bNeedVel,bGenVel,fr_time,ir,&state,sys,oenv);
1564 if (ir->ePBC==epbcXY && ir->nwall!=2)
1566 clear_rvec(state.box[ZZ]);
1571 set_warning_line(wi,mdparin,-1);
1572 check_chargegroup_radii(sys,ir,state.x,wi);
1575 if (EEL_FULL(ir->coulombtype)) {
1576 /* Calculate the optimal grid dimensions */
1577 copy_mat(state.box,box);
1578 if (ir->ePBC==epbcXY && ir->nwall==2)
1579 svmul(ir->wall_ewald_zfac,box[ZZ],box[ZZ]);
1580 max_spacing = calc_grid(stdout,box,opts->fourierspacing,
1581 &(ir->nkx),&(ir->nky),&(ir->nkz));
1584 if (ir->ePull != epullNO)
1585 set_pull_init(ir,sys,state.x,state.box,oenv,opts->pull_start);
1589 set_reference_positions(ir->rot,sys,state.x,state.box,
1590 opt2fn("-ref",NFILE,fnm),opt2bSet("-ref",NFILE,fnm),
1594 /* reset_multinr(sys); */
1596 if (EEL_PME(ir->coulombtype)) {
1597 float ratio = pme_load_estimate(sys,ir,state.box);
1598 fprintf(stderr,"Estimate for the relative computational load of the PME mesh part: %.2f\n",ratio);
1599 /* With free energy we might need to do PME both for the A and B state
1600 * charges. This will double the cost, but the optimal performance will
1601 * then probably be at a slightly larger cut-off and grid spacing.
1603 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1604 (ir->efep != efepNO && ratio > 2.0/3.0)) {
1606 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1607 "and for highly parallel simulations between 0.25 and 0.33,\n"
1608 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1609 if (ir->efep != efepNO) {
1611 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1617 char warn_buf[STRLEN];
1618 double cio = compute_io(ir,sys->natoms,&sys->groups,F_NRE,1);
1619 sprintf(warn_buf,"This run will generate roughly %.0f Mb of data",cio);
1621 set_warning_line(wi,mdparin,-1);
1622 warning_note(wi,warn_buf);
1624 printf("%s\n",warn_buf);
1628 /* MRS: eventually figure out better logic for initializing the fep
1629 values that makes declaring the lambda and declaring the state not
1630 potentially conflict if not handled correctly. */
1631 if (ir->efep != efepNO)
1633 state.fep_state = ir->fepvals->init_fep_state;
1634 for (i=0;i<efptNR;i++)
1636 /* init_lambda trumps state definitions*/
1637 if (ir->fepvals->init_lambda >= 0)
1639 state.lambda[i] = ir->fepvals->init_lambda;
1643 if (ir->fepvals->all_lambda[i] == NULL)
1645 gmx_fatal(FARGS,"Values of lambda not set for a free energy calculation!");
1649 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1656 fprintf(stderr,"writing run input file...\n");
1658 done_warning(wi,FARGS);
1660 write_tpx_state(ftp2fn(efTPX,NFILE,fnm),ir,&state,sys);