2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
42 #include <sys/types.h>
66 #include "sortwater.h"
68 #include "gmx_fatal.h"
74 #include "vsite_parm.h"
80 #include "compute_io.h"
81 #include "gpp_atomtype.h"
82 #include "gpp_tomorse.h"
83 #include "mtop_util.h"
85 #include "calc_verletbuf.h"
87 static int rm_interactions(int ifunc,int nrmols,t_molinfo mols[])
92 /* For all the molecule types */
93 for(i=0; i<nrmols; i++) {
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)
108 gmx_incons("comparing atom names");
112 for(mb=0; mb<mtop->nmolblock; mb++) {
113 tat = &mtop->moltype[mtop->molblock[mb].type].atoms;
114 for(m=0; m<mtop->molblock[mb].nmol; m++) {
115 for(j=0; j < tat->nr; j++) {
116 if (strcmp( *(tat->atomname[j]) , *(at->atomname[i]) ) != 0) {
117 if (nmismatch < MAXMISMATCH) {
119 "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
120 i+1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
121 } else if (nmismatch == MAXMISMATCH) {
122 fprintf(stderr,"(more than %d non-matching atom names)\n",MAXMISMATCH);
134 static void check_eg_vs_cg(gmx_mtop_t *mtop)
136 int astart,mb,m,cg,j,firstj;
137 unsigned char firsteg,eg;
140 /* Go through all the charge groups and make sure all their
141 * atoms are in the same energy group.
145 for(mb=0; mb<mtop->nmolblock; mb++) {
146 molt = &mtop->moltype[mtop->molblock[mb].type];
147 for(m=0; m<mtop->molblock[mb].nmol; m++) {
148 for(cg=0; cg<molt->cgs.nr;cg++) {
149 /* Get the energy group of the first atom in this charge group */
150 firstj = astart + molt->cgs.index[cg];
151 firsteg = ggrpnr(&mtop->groups,egcENER,firstj);
152 for(j=molt->cgs.index[cg]+1;j<molt->cgs.index[cg+1];j++) {
153 eg = ggrpnr(&mtop->groups,egcENER,astart+j);
155 gmx_fatal(FARGS,"atoms %d and %d in charge group %d of molecule type '%s' are in different energy groups",
156 firstj+1,astart+j+1,cg+1,*molt->name);
160 astart += molt->atoms.nr;
165 static void check_cg_sizes(const char *topfn,t_block *cgs,warninp_t wi)
168 char warn_buf[STRLEN];
171 for(cg=0; cg<cgs->nr; cg++)
173 maxsize = max(maxsize,cgs->index[cg+1]-cgs->index[cg]);
176 if (maxsize > MAX_CHARGEGROUP_SIZE)
178 gmx_fatal(FARGS,"The largest charge group contains %d atoms. The maximum is %d.",maxsize,MAX_CHARGEGROUP_SIZE);
180 else if (maxsize > 10)
182 set_warning_line(wi,topfn,-1);
184 "The largest charge group contains %d atoms.\n"
185 "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"
186 "For efficiency and accuracy, charge group should consist of a few atoms.\n"
187 "For all-atom force fields use: CH3, CH2, CH, NH2, NH, OH, CO2, CO, etc.",
189 warning_note(wi,warn_buf);
193 static void check_bonds_timestep(gmx_mtop_t *mtop,double dt,warninp_t wi)
195 /* This check is not intended to ensure accurate integration,
196 * rather it is to signal mistakes in the mdp settings.
197 * A common mistake is to forget to turn on constraints
198 * for MD after energy minimization with flexible bonds.
199 * This check can also detect too large time steps for flexible water
200 * models, but such errors will often be masked by the constraints
201 * mdp options, which turns flexible water into water with bond constraints,
202 * but without an angle constraint. Unfortunately such incorrect use
203 * of water models can not easily be detected without checking
204 * for specific model names.
206 * The stability limit of leap-frog or velocity verlet is 4.44 steps
207 * per oscillational period.
208 * But accurate bonds distributions are lost far before that limit.
209 * To allow relatively common schemes (although not common with Gromacs)
210 * of dt=1 fs without constraints and dt=2 fs with only H-bond constraints
211 * we set the note limit to 10.
213 int min_steps_warn=5;
214 int min_steps_note=10;
217 gmx_moltype_t *moltype,*w_moltype;
219 t_ilist *ilist,*ilb,*ilc,*ils;
221 int i,a1,a2,w_a1,w_a2,j;
222 real twopi2,limit2,fc,re,m1,m2,period2,w_period2;
223 gmx_bool bFound,bWater,bWarn;
224 char warn_buf[STRLEN];
226 ip = mtop->ffparams.iparams;
228 twopi2 = sqr(2*M_PI);
230 limit2 = sqr(min_steps_note*dt);
236 for(molt=0; molt<mtop->nmoltype; molt++)
238 moltype = &mtop->moltype[molt];
239 atom = moltype->atoms.atom;
240 ilist = moltype->ilist;
241 ilc = &ilist[F_CONSTR];
242 ils = &ilist[F_SETTLE];
243 for(ftype=0; ftype<F_NRE; ftype++)
245 if (!(ftype == F_BONDS || ftype == F_G96BONDS || ftype == F_HARMONIC))
251 for(i=0; i<ilb->nr; i+=3)
253 fc = ip[ilb->iatoms[i]].harmonic.krA;
254 re = ip[ilb->iatoms[i]].harmonic.rA;
255 if (ftype == F_G96BONDS)
257 /* Convert squared sqaure fc to harmonic fc */
260 a1 = ilb->iatoms[i+1];
261 a2 = ilb->iatoms[i+2];
264 if (fc > 0 && m1 > 0 && m2 > 0)
266 period2 = twopi2*m1*m2/((m1 + m2)*fc);
270 period2 = GMX_FLOAT_MAX;
274 fprintf(debug,"fc %g m1 %g m2 %g period %g\n",
275 fc,m1,m2,sqrt(period2));
277 if (period2 < limit2)
280 for(j=0; j<ilc->nr; j+=3)
282 if ((ilc->iatoms[j+1] == a1 && ilc->iatoms[j+2] == a2) ||
283 (ilc->iatoms[j+1] == a2 && ilc->iatoms[j+2] == a1))
288 for(j=0; j<ils->nr; j+=4)
290 if ((a1 == ils->iatoms[j+1] || a1 == ils->iatoms[j+2] || a1 == ils->iatoms[j+3]) &&
291 (a2 == ils->iatoms[j+1] || a2 == ils->iatoms[j+2] || a2 == ils->iatoms[j+3]))
297 (w_moltype == NULL || period2 < w_period2))
309 if (w_moltype != NULL)
311 bWarn = (w_period2 < sqr(min_steps_warn*dt));
312 /* A check that would recognize most water models */
313 bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' &&
314 w_moltype->atoms.nr <= 5);
315 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"
318 w_a1+1,*w_moltype->atoms.atomname[w_a1],
319 w_a2+1,*w_moltype->atoms.atomname[w_a2],
320 sqrt(w_period2),bWarn ? min_steps_warn : min_steps_note,dt,
322 "Maybe you asked for fexible water." :
323 "Maybe you forgot to change the constraints mdp option.");
326 warning(wi,warn_buf);
330 warning_note(wi,warn_buf);
335 static void check_vel(gmx_mtop_t *mtop,rvec v[])
337 gmx_mtop_atomloop_all_t aloop;
341 aloop = gmx_mtop_atomloop_all_init(mtop);
342 while (gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
343 if (atom->ptype == eptShell ||
344 atom->ptype == eptBond ||
345 atom->ptype == eptVSite) {
351 static gmx_bool nint_ftype(gmx_mtop_t *mtop,t_molinfo *mi,int ftype)
356 for(mb=0; mb<mtop->nmolblock; mb++) {
357 nint += mtop->molblock[mb].nmol*mi[mtop->molblock[mb].type].plist[ftype].nr;
363 /* This routine reorders the molecule type array
364 * in the order of use in the molblocks,
365 * unused molecule types are deleted.
367 static void renumber_moltypes(gmx_mtop_t *sys,
368 int *nmolinfo,t_molinfo **molinfo)
374 snew(order,*nmolinfo);
376 for(mb=0; mb<sys->nmolblock; mb++) {
377 for(i=0; i<norder; i++) {
378 if (order[i] == sys->molblock[mb].type) {
383 /* This type did not occur yet, add it */
384 order[norder] = sys->molblock[mb].type;
385 /* Renumber the moltype in the topology */
388 sys->molblock[mb].type = i;
391 /* We still need to reorder the molinfo structs */
393 for(mi=0; mi<*nmolinfo; mi++) {
394 for(i=0; i<norder; i++) {
395 if (order[i] == mi) {
400 done_mi(&(*molinfo)[mi]);
402 minew[i] = (*molinfo)[mi];
411 static void molinfo2mtop(int nmi,t_molinfo *mi,gmx_mtop_t *mtop)
416 mtop->nmoltype = nmi;
417 snew(mtop->moltype,nmi);
418 for(m=0; m<nmi; m++) {
419 molt = &mtop->moltype[m];
420 molt->name = mi[m].name;
421 molt->atoms = mi[m].atoms;
422 /* ilists are copied later */
423 molt->cgs = mi[m].cgs;
424 molt->excls = mi[m].excls;
429 new_status(const char *topfile,const char *topppfile,const char *confin,
430 t_gromppopts *opts,t_inputrec *ir,gmx_bool bZero,
431 gmx_bool bGenVel,gmx_bool bVerbose,t_state *state,
432 gpp_atomtype_t atype,gmx_mtop_t *sys,
433 int *nmi,t_molinfo **mi,t_params plist[],
434 int *comb,double *reppow,real *fudgeQQ,
438 t_molinfo *molinfo=NULL;
440 gmx_molblock_t *molblock,*molbs;
442 int mb,i,nrmols,nmismatch;
445 char warn_buf[STRLEN];
449 /* Set gmx_boolean for GB */
450 if(ir->implicit_solvent)
453 /* TOPOLOGY processing */
454 sys->name = do_top(bVerbose,topfile,topppfile,opts,bZero,&(sys->symtab),
455 plist,comb,reppow,fudgeQQ,
456 atype,&nrmols,&molinfo,ir,
457 &nmolblock,&molblock,bGB,
461 snew(sys->molblock,nmolblock);
464 for(mb=0; mb<nmolblock; mb++) {
465 if (sys->nmolblock > 0 &&
466 molblock[mb].type == sys->molblock[sys->nmolblock-1].type) {
467 /* Merge consecutive blocks with the same molecule type */
468 sys->molblock[sys->nmolblock-1].nmol += molblock[mb].nmol;
469 sys->natoms += molblock[mb].nmol*sys->molblock[sys->nmolblock-1].natoms_mol;
470 } else if (molblock[mb].nmol > 0) {
471 /* Add a new molblock to the topology */
472 molbs = &sys->molblock[sys->nmolblock];
473 *molbs = molblock[mb];
474 molbs->natoms_mol = molinfo[molbs->type].atoms.nr;
475 molbs->nposres_xA = 0;
476 molbs->nposres_xB = 0;
477 sys->natoms += molbs->nmol*molbs->natoms_mol;
481 if (sys->nmolblock == 0) {
482 gmx_fatal(FARGS,"No molecules were defined in the system");
485 renumber_moltypes(sys,&nrmols,&molinfo);
488 convert_harmonics(nrmols,molinfo,atype);
490 if (ir->eDisre == edrNone) {
491 i = rm_interactions(F_DISRES,nrmols,molinfo);
493 set_warning_line(wi,"unknown",-1);
494 sprintf(warn_buf,"disre = no, removed %d distance restraints",i);
495 warning_note(wi,warn_buf);
498 if (opts->bOrire == FALSE) {
499 i = rm_interactions(F_ORIRES,nrmols,molinfo);
501 set_warning_line(wi,"unknown",-1);
502 sprintf(warn_buf,"orire = no, removed %d orientation restraints",i);
503 warning_note(wi,warn_buf);
507 /* Copy structures from msys to sys */
508 molinfo2mtop(nrmols,molinfo,sys);
510 gmx_mtop_finalize(sys);
512 /* COORDINATE file processing */
514 fprintf(stderr,"processing coordinates...\n");
516 get_stx_coordnum(confin,&state->natoms);
517 if (state->natoms != sys->natoms)
518 gmx_fatal(FARGS,"number of coordinates in coordinate file (%s, %d)\n"
519 " does not match topology (%s, %d)",
520 confin,state->natoms,topfile,sys->natoms);
522 /* make space for coordinates and velocities */
525 init_t_atoms(confat,state->natoms,FALSE);
526 init_state(state,state->natoms,0,0,0,0);
527 read_stx_conf(confin,title,confat,state->x,state->v,NULL,state->box);
528 /* This call fixes the box shape for runs with pressure scaling */
529 set_box_rel(ir,state);
531 nmismatch = check_atom_names(topfile, confin, sys, confat);
532 free_t_atoms(confat,TRUE);
536 sprintf(buf,"%d non-matching atom name%s\n"
537 "atom names from %s will be used\n"
538 "atom names from %s will be ignored\n",
539 nmismatch,(nmismatch == 1) ? "" : "s",topfile,confin);
543 fprintf(stderr,"double-checking input for internal consistency...\n");
544 double_check(ir,state->box,nint_ftype(sys,molinfo,F_CONSTR),wi);
549 gmx_mtop_atomloop_all_t aloop;
552 snew(mass,state->natoms);
553 aloop = gmx_mtop_atomloop_all_init(sys);
554 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
558 if (opts->seed == -1) {
559 opts->seed = make_seed();
560 fprintf(stderr,"Setting gen_seed to %d\n",opts->seed);
562 maxwell_speed(opts->tempi,opts->seed,sys,state->v);
564 stop_cm(stdout,state->natoms,mass,state->x,state->v);
572 static void copy_state(const char *slog,t_trxframe *fr,
573 gmx_bool bReadVel,t_state *state,
578 if (fr->not_ok & FRAME_NOT_OK)
580 gmx_fatal(FARGS,"Can not start from an incomplete frame");
584 gmx_fatal(FARGS,"Did not find a frame with coordinates in file %s",
588 for(i=0; i<state->natoms; i++)
590 copy_rvec(fr->x[i],state->x[i]);
596 gmx_incons("Trajecory frame unexpectedly does not contain velocities");
598 for(i=0; i<state->natoms; i++)
600 copy_rvec(fr->v[i],state->v[i]);
605 copy_mat(fr->box,state->box);
608 *use_time = fr->time;
611 static void cont_status(const char *slog,const char *ener,
612 gmx_bool bNeedVel,gmx_bool bGenVel, real fr_time,
613 t_inputrec *ir,t_state *state,
615 const output_env_t oenv)
616 /* If fr_time == -1 read the last frame available which is complete */
624 bReadVel = (bNeedVel && !bGenVel);
627 "Reading Coordinates%s and Box size from old trajectory\n",
628 bReadVel ? ", Velocities" : "");
631 fprintf(stderr,"Will read whole trajectory\n");
635 fprintf(stderr,"Will read till time %g\n",fr_time);
641 fprintf(stderr,"Velocities generated: "
642 "ignoring velocities in input trajectory\n");
644 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
648 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X | TRX_NEED_V);
654 "WARNING: Did not find a frame with velocities in file %s,\n"
655 " all velocities will be set to zero!\n\n",slog);
656 for(i=0; i<sys->natoms; i++)
658 clear_rvec(state->v[i]);
661 /* Search for a frame without velocities */
663 read_first_frame(oenv,&fp,slog,&fr,TRX_NEED_X);
667 state->natoms = fr.natoms;
669 if (sys->natoms != state->natoms)
671 gmx_fatal(FARGS,"Number of atoms in Topology "
672 "is not the same as in Trajectory");
674 copy_state(slog,&fr,bReadVel,state,&use_time);
676 /* Find the appropriate frame */
677 while ((fr_time == -1 || fr.time < fr_time) &&
678 read_next_frame(oenv,fp,&fr))
680 copy_state(slog,&fr,bReadVel,state,&use_time);
685 /* Set the relative box lengths for preserving the box shape.
686 * Note that this call can lead to differences in the last bit
687 * with respect to using tpbconv to create a [TT].tpx[tt] file.
689 set_box_rel(ir,state);
691 fprintf(stderr,"Using frame at t = %g ps\n",use_time);
692 fprintf(stderr,"Starting time for run is %g ps\n",ir->init_t);
694 if ((ir->epc != epcNO || ir->etc ==etcNOSEHOOVER) && ener)
696 get_enx_state(ener,use_time,&sys->groups,ir,state);
697 preserve_box_shape(ir,state->box_rel,state->boxv);
701 static void read_posres(gmx_mtop_t *mtop,t_molinfo *molinfo,gmx_bool bTopB,
703 int rc_scaling, int ePBC,
707 gmx_bool bFirst = TRUE;
713 int natoms,npbcdim=0;
714 char warn_buf[STRLEN],title[STRLEN];
715 int a,i,ai,j,k,mb,nat_molb;
716 gmx_molblock_t *molb;
720 get_stx_coordnum(fn,&natoms);
721 if (natoms != mtop->natoms) {
722 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);
723 warning(wi,warn_buf);
727 init_t_atoms(&dumat,natoms,FALSE);
728 read_stx_conf(fn,title,&dumat,x,v,NULL,box);
730 npbcdim = ePBC2npbcdim(ePBC);
732 if (rc_scaling != erscNO) {
733 copy_mat(box,invbox);
734 for(j=npbcdim; j<DIM; j++) {
735 clear_rvec(invbox[j]);
738 m_inv_ur0(invbox,invbox);
741 /* Copy the reference coordinates to mtop */
745 for(mb=0; mb<mtop->nmolblock; mb++) {
746 molb = &mtop->molblock[mb];
747 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
748 pr = &(molinfo[molb->type].plist[F_POSRES]);
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);
757 if (rc_scaling == erscCOM) {
758 /* Determine the center of mass of the posres reference coordinates */
759 for(j=0; j<npbcdim; j++) {
760 sum[j] += atom[ai].m*x[a+ai][j];
762 totmass += atom[ai].m;
766 molb->nposres_xA = nat_molb;
767 snew(molb->posres_xA,molb->nposres_xA);
768 for(i=0; i<nat_molb; i++) {
769 copy_rvec(x[a+i],molb->posres_xA[i]);
772 molb->nposres_xB = nat_molb;
773 snew(molb->posres_xB,molb->nposres_xB);
774 for(i=0; i<nat_molb; i++) {
775 copy_rvec(x[a+i],molb->posres_xB[i]);
781 if (rc_scaling == erscCOM) {
783 gmx_fatal(FARGS,"The total mass of the position restraint atoms is 0");
784 for(j=0; j<npbcdim; j++)
785 com[j] = sum[j]/totmass;
786 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]);
789 if (rc_scaling != erscNO) {
790 for(mb=0; mb<mtop->nmolblock; mb++) {
791 molb = &mtop->molblock[mb];
792 nat_molb = molb->nmol*mtop->moltype[molb->type].atoms.nr;
793 if (molb->nposres_xA > 0 || molb->nposres_xB > 0) {
794 xp = (!bTopB ? molb->posres_xA : molb->posres_xB);
795 for(i=0; i<nat_molb; i++) {
796 for(j=0; j<npbcdim; j++) {
797 if (rc_scaling == erscALL) {
798 /* Convert from Cartesian to crystal coordinates */
799 xp[i][j] *= invbox[j][j];
800 for(k=j+1; k<npbcdim; k++) {
801 xp[i][j] += invbox[k][j]*xp[i][k];
803 } else if (rc_scaling == erscCOM) {
804 /* Subtract the center of mass */
812 if (rc_scaling == erscCOM) {
813 /* Convert the COM from Cartesian to crystal coordinates */
814 for(j=0; j<npbcdim; j++) {
815 com[j] *= invbox[j][j];
816 for(k=j+1; k<npbcdim; k++) {
817 com[j] += invbox[k][j]*com[k];
823 free_t_atoms(&dumat,TRUE);
828 static void gen_posres(gmx_mtop_t *mtop,t_molinfo *mi,
829 char *fnA, char *fnB,
830 int rc_scaling, int ePBC,
836 read_posres (mtop,mi,FALSE,fnA,rc_scaling,ePBC,com,wi);
837 if (strcmp(fnA,fnB) != 0) {
838 read_posres(mtop,mi,TRUE ,fnB,rc_scaling,ePBC,comB,wi);
842 static void set_wall_atomtype(gpp_atomtype_t at,t_gromppopts *opts,
843 t_inputrec *ir,warninp_t wi)
846 char warn_buf[STRLEN];
850 fprintf(stderr,"Searching the wall atom type(s)\n");
852 for(i=0; i<ir->nwall; i++)
854 ir->wall_atomtype[i] = get_atomtype_type(opts->wall_atomtype[i],at);
855 if (ir->wall_atomtype[i] == NOTSET)
857 sprintf(warn_buf,"Specified wall atom type %s is not defined",opts->wall_atomtype[i]);
858 warning_error(wi,warn_buf);
863 static int nrdf_internal(t_atoms *atoms)
868 for(i=0; i<atoms->nr; i++) {
869 /* Vsite ptype might not be set here yet, so also check the mass */
870 if ((atoms->atom[i].ptype == eptAtom ||
871 atoms->atom[i].ptype == eptNucleus)
872 && atoms->atom[i].m > 0) {
877 case 0: nrdf = 0; break;
878 case 1: nrdf = 0; break;
879 case 2: nrdf = 1; break;
880 default: nrdf = nmass*3 - 6; break;
903 q = (y[i+1]-2.0*y[i]+y[i-1])/dx;
904 u[i] = (3.0*q/dx-0.5*u[i-1])/p;
911 y2[i] = y2[i]*y2[i+1]+u[i];
917 interpolate1d( double xmin,
930 a = (xmin+(ix+1)*dx-x)/dx;
931 b = (x-xmin-ix*dx)/dx;
933 *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;
934 *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];
939 setup_cmap (int grid_spacing,
942 gmx_cmap_t * cmap_grid)
944 double *tmp_u,*tmp_u2,*tmp_yy,*tmp_y1,*tmp_t2,*tmp_grid;
946 int i,j,k,ii,jj,kk,idx;
948 double dx,xmin,v,v1,v2,v12;
951 snew(tmp_u,2*grid_spacing);
952 snew(tmp_u2,2*grid_spacing);
953 snew(tmp_yy,2*grid_spacing);
954 snew(tmp_y1,2*grid_spacing);
955 snew(tmp_t2,2*grid_spacing*2*grid_spacing);
956 snew(tmp_grid,2*grid_spacing*2*grid_spacing);
958 dx = 360.0/grid_spacing;
959 xmin = -180.0-dx*grid_spacing/2;
963 /* Compute an offset depending on which cmap we are using
964 * Offset will be the map number multiplied with the
965 * grid_spacing * grid_spacing * 2
967 offset = kk * grid_spacing * grid_spacing * 2;
969 for(i=0;i<2*grid_spacing;i++)
971 ii=(i+grid_spacing-grid_spacing/2)%grid_spacing;
973 for(j=0;j<2*grid_spacing;j++)
975 jj=(j+grid_spacing-grid_spacing/2)%grid_spacing;
976 tmp_grid[i*grid_spacing*2+j] = grid[offset+ii*grid_spacing+jj];
980 for(i=0;i<2*grid_spacing;i++)
982 spline1d(dx,&(tmp_grid[2*grid_spacing*i]),2*grid_spacing,tmp_u,&(tmp_t2[2*grid_spacing*i]));
985 for(i=grid_spacing/2;i<grid_spacing+grid_spacing/2;i++)
987 ii = i-grid_spacing/2;
990 for(j=grid_spacing/2;j<grid_spacing+grid_spacing/2;j++)
992 jj = j-grid_spacing/2;
995 for(k=0;k<2*grid_spacing;k++)
997 interpolate1d(xmin,dx,&(tmp_grid[2*grid_spacing*k]),
998 &(tmp_t2[2*grid_spacing*k]),psi,&tmp_yy[k],&tmp_y1[k]);
1001 spline1d(dx,tmp_yy,2*grid_spacing,tmp_u,tmp_u2);
1002 interpolate1d(xmin,dx,tmp_yy,tmp_u2,phi,&v,&v1);
1003 spline1d(dx,tmp_y1,2*grid_spacing,tmp_u,tmp_u2);
1004 interpolate1d(xmin,dx,tmp_y1,tmp_u2,phi,&v2,&v12);
1006 idx = ii*grid_spacing+jj;
1007 cmap_grid->cmapdata[kk].cmap[idx*4] = grid[offset+ii*grid_spacing+jj];
1008 cmap_grid->cmapdata[kk].cmap[idx*4+1] = v1;
1009 cmap_grid->cmapdata[kk].cmap[idx*4+2] = v2;
1010 cmap_grid->cmapdata[kk].cmap[idx*4+3] = v12;
1016 void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
1020 cmap_grid->ngrid = ngrid;
1021 cmap_grid->grid_spacing = grid_spacing;
1022 nelem = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
1024 snew(cmap_grid->cmapdata,ngrid);
1026 for(i=0;i<cmap_grid->ngrid;i++)
1028 snew(cmap_grid->cmapdata[i].cmap,4*nelem);
1033 static int count_constraints(gmx_mtop_t *mtop,t_molinfo *mi,warninp_t wi)
1035 int count,count_mol,i,mb;
1036 gmx_molblock_t *molb;
1041 for(mb=0; mb<mtop->nmolblock; mb++) {
1043 molb = &mtop->molblock[mb];
1044 plist = mi[molb->type].plist;
1046 for(i=0; i<F_NRE; i++) {
1048 count_mol += 3*plist[i].nr;
1049 else if (interaction_function[i].flags & IF_CONSTRAINT)
1050 count_mol += plist[i].nr;
1053 if (count_mol > nrdf_internal(&mi[molb->type].atoms)) {
1055 "Molecule type '%s' has %d constraints.\n"
1056 "For stability and efficiency there should not be more constraints than internal number of degrees of freedom: %d.\n",
1057 *mi[molb->type].name,count_mol,
1058 nrdf_internal(&mi[molb->type].atoms));
1061 count += molb->nmol*count_mol;
1067 static void check_gbsa_params_charged(gmx_mtop_t *sys, gpp_atomtype_t atype)
1069 int i,nmiss,natoms,mt;
1071 const t_atoms *atoms;
1074 for(mt=0;mt<sys->nmoltype;mt++)
1076 atoms = &sys->moltype[mt].atoms;
1079 for(i=0;i<natoms;i++)
1081 q = atoms->atom[i].q;
1082 if ((get_atomtype_radius(atoms->atom[i].type,atype) == 0 ||
1083 get_atomtype_vol(atoms->atom[i].type,atype) == 0 ||
1084 get_atomtype_surftens(atoms->atom[i].type,atype) == 0 ||
1085 get_atomtype_gb_radius(atoms->atom[i].type,atype) == 0 ||
1086 get_atomtype_S_hct(atoms->atom[i].type,atype) == 0) &&
1089 fprintf(stderr,"\nGB parameter(s) zero for atom type '%s' while charge is %g\n",
1090 get_atomtype_name(atoms->atom[i].type,atype),q);
1098 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);
1103 static void check_gbsa_params(t_inputrec *ir,gpp_atomtype_t atype)
1107 /* If we are doing GBSA, check that we got the parameters we need
1108 * This checking is to see if there are GBSA paratmeters for all
1109 * atoms in the force field. To go around this for testing purposes
1110 * comment out the nerror++ counter temporarily
1113 for(i=0;i<get_atomtype_ntypes(atype);i++)
1115 if (get_atomtype_radius(i,atype) < 0 ||
1116 get_atomtype_vol(i,atype) < 0 ||
1117 get_atomtype_surftens(i,atype) < 0 ||
1118 get_atomtype_gb_radius(i,atype) < 0 ||
1119 get_atomtype_S_hct(i,atype) < 0)
1121 fprintf(stderr,"\nGB parameter(s) missing or negative for atom type '%s'\n",
1122 get_atomtype_name(i,atype));
1129 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);
1134 static void set_verlet_buffer(const gmx_mtop_t *mtop,
1137 real verletbuf_drift,
1142 verletbuf_list_setup_t ls;
1145 char warn_buf[STRLEN];
1148 for(i=0; i<ir->opts.ngtc; i++)
1150 if (ir->opts.ref_t[i] < 0)
1152 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.");
1156 ref_T = max(ref_T,ir->opts.ref_t[i]);
1160 printf("Determining Verlet buffer for an energy drift of %g kJ/mol/ps at %g K\n",verletbuf_drift,ref_T);
1162 for(i=0; i<ir->opts.ngtc; i++)
1164 if (ir->opts.ref_t[i] >= 0 && ir->opts.ref_t[i] != ref_T)
1166 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.",
1167 ir->opts.nrdf[i],ir->opts.ref_t[i],ref_T);
1168 warning_note(wi,warn_buf);
1172 /* Calculate the buffer size for simple atom vs atoms list */
1173 ls.cluster_size_i = 1;
1174 ls.cluster_size_j = 1;
1175 calc_verlet_buffer_size(mtop,det(box),ir,verletbuf_drift,
1176 &ls,&n_nonlin_vsite,&rlist_1x1);
1178 /* Set the pair-list buffer size in ir */
1179 verletbuf_get_list_setup(FALSE,&ls);
1180 calc_verlet_buffer_size(mtop,det(box),ir,verletbuf_drift,
1181 &ls,&n_nonlin_vsite,&ir->rlist);
1183 if (n_nonlin_vsite > 0)
1185 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);
1186 warning_note(wi,warn_buf);
1189 printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n",
1190 1,1,rlist_1x1,rlist_1x1-max(ir->rvdw,ir->rcoulomb));
1192 ir->rlistlong = ir->rlist;
1193 printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
1194 ls.cluster_size_i,ls.cluster_size_j,
1195 ir->rlist,ir->rlist-max(ir->rvdw,ir->rcoulomb));
1197 if (sqr(ir->rlistlong) >= max_cutoff2(ir->ePBC,box))
1199 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)));
1203 int cmain (int argc, char *argv[])
1205 static const char *desc[] = {
1206 "The gromacs preprocessor",
1207 "reads a molecular topology file, checks the validity of the",
1208 "file, expands the topology from a molecular description to an atomic",
1209 "description. The topology file contains information about",
1210 "molecule types and the number of molecules, the preprocessor",
1211 "copies each molecule as needed. ",
1212 "There is no limitation on the number of molecule types. ",
1213 "Bonds and bond-angles can be converted into constraints, separately",
1214 "for hydrogens and heavy atoms.",
1215 "Then a coordinate file is read and velocities can be generated",
1216 "from a Maxwellian distribution if requested.",
1217 "[TT]grompp[tt] also reads parameters for the [TT]mdrun[tt] ",
1218 "(eg. number of MD steps, time step, cut-off), and others such as",
1219 "NEMD parameters, which are corrected so that the net acceleration",
1221 "Eventually a binary file is produced that can serve as the sole input",
1222 "file for the MD program.[PAR]",
1224 "[TT]grompp[tt] uses the atom names from the topology file. The atom names",
1225 "in the coordinate file (option [TT]-c[tt]) are only read to generate",
1226 "warnings when they do not match the atom names in the topology.",
1227 "Note that the atom names are irrelevant for the simulation as",
1228 "only the atom types are used for generating interaction parameters.[PAR]",
1230 "[TT]grompp[tt] uses a built-in preprocessor to resolve includes, macros, ",
1231 "etc. The preprocessor supports the following keywords:[PAR]",
1232 "#ifdef VARIABLE[BR]",
1233 "#ifndef VARIABLE[BR]",
1236 "#define VARIABLE[BR]",
1237 "#undef VARIABLE[BR]"
1238 "#include \"filename\"[BR]",
1239 "#include <filename>[PAR]",
1240 "The functioning of these statements in your topology may be modulated by",
1241 "using the following two flags in your [TT].mdp[tt] file:[PAR]",
1242 "[TT]define = -DVARIABLE1 -DVARIABLE2[BR]",
1243 "include = -I/home/john/doe[tt][BR]",
1244 "For further information a C-programming textbook may help you out.",
1245 "Specifying the [TT]-pp[tt] flag will get the pre-processed",
1246 "topology file written out so that you can verify its contents.[PAR]",
1248 /* cpp has been unnecessary for some time, hasn't it?
1249 "If your system does not have a C-preprocessor, you can still",
1250 "use [TT]grompp[tt], but you do not have access to the features ",
1251 "from the cpp. Command line options to the C-preprocessor can be given",
1252 "in the [TT].mdp[tt] file. See your local manual (man cpp).[PAR]",
1255 "When using position restraints a file with restraint coordinates",
1256 "can be supplied with [TT]-r[tt], otherwise restraining will be done",
1257 "with respect to the conformation from the [TT]-c[tt] option.",
1258 "For free energy calculation the the coordinates for the B topology",
1259 "can be supplied with [TT]-rb[tt], otherwise they will be equal to",
1260 "those of the A topology.[PAR]",
1262 "Starting coordinates can be read from trajectory with [TT]-t[tt].",
1263 "The last frame with coordinates and velocities will be read,",
1264 "unless the [TT]-time[tt] option is used. Only if this information",
1265 "is absent will the coordinates in the [TT]-c[tt] file be used.",
1266 "Note that these velocities will not be used when [TT]gen_vel = yes[tt]",
1267 "in your [TT].mdp[tt] file. An energy file can be supplied with",
1268 "[TT]-e[tt] to read Nose-Hoover and/or Parrinello-Rahman coupling",
1271 "[TT]grompp[tt] can be used to restart simulations (preserving",
1272 "continuity) by supplying just a checkpoint file with [TT]-t[tt].",
1273 "However, for simply changing the number of run steps to extend",
1274 "a run, using [TT]tpbconv[tt] is more convenient than [TT]grompp[tt].",
1275 "You then supply the old checkpoint file directly to [TT]mdrun[tt]",
1276 "with [TT]-cpi[tt]. If you wish to change the ensemble or things",
1277 "like output frequency, then supplying the checkpoint file to",
1278 "[TT]grompp[tt] with [TT]-t[tt] along with a new [TT].mdp[tt] file",
1279 "with [TT]-f[tt] is the recommended procedure.[PAR]",
1281 "By default, all bonded interactions which have constant energy due to",
1282 "virtual site constructions will be removed. If this constant energy is",
1283 "not zero, this will result in a shift in the total energy. All bonded",
1284 "interactions can be kept by turning off [TT]-rmvsbds[tt]. Additionally,",
1285 "all constraints for distances which will be constant anyway because",
1286 "of virtual site constructions will be removed. If any constraints remain",
1287 "which involve virtual sites, a fatal error will result.[PAR]"
1289 "To verify your run input file, please take note of all warnings",
1290 "on the screen, and correct where necessary. Do also look at the contents",
1291 "of the [TT]mdout.mdp[tt] file; this contains comment lines, as well as",
1292 "the input that [TT]grompp[tt] has read. If in doubt, you can start [TT]grompp[tt]",
1293 "with the [TT]-debug[tt] option which will give you more information",
1294 "in a file called [TT]grompp.log[tt] (along with real debug info). You",
1295 "can see the contents of the run input file with the [TT]gmxdump[tt]",
1296 "program. [TT]gmxcheck[tt] can be used to compare the contents of two",
1297 "run input files.[PAR]"
1299 "The [TT]-maxwarn[tt] option can be used to override warnings printed",
1300 "by [TT]grompp[tt] that otherwise halt output. In some cases, warnings are",
1301 "harmless, but usually they are not. The user is advised to carefully",
1302 "interpret the output messages before attempting to bypass them with",
1309 gpp_atomtype_t atype;
1311 int natoms,nvsite,comb,mt;
1315 real max_spacing,fudgeQQ;
1317 char fn[STRLEN],fnB[STRLEN];
1318 const char *mdparin;
1320 gmx_bool bNeedVel,bGenVel;
1321 gmx_bool have_atomnumber;
1323 t_params *gb_plist = NULL;
1324 gmx_genborn_t *born = NULL;
1326 gmx_bool bVerbose = FALSE;
1328 char warn_buf[STRLEN];
1331 { efMDP, NULL, NULL, ffREAD },
1332 { efMDP, "-po", "mdout", ffWRITE },
1333 { efSTX, "-c", NULL, ffREAD },
1334 { efSTX, "-r", NULL, ffOPTRD },
1335 { efSTX, "-rb", NULL, ffOPTRD },
1336 { efNDX, NULL, NULL, ffOPTRD },
1337 { efTOP, NULL, NULL, ffREAD },
1338 { efTOP, "-pp", "processed", ffOPTWR },
1339 { efTPX, "-o", NULL, ffWRITE },
1340 { efTRN, "-t", NULL, ffOPTRD },
1341 { efEDR, "-e", NULL, ffOPTRD },
1342 { efTRN, "-ref","rotref", ffOPTRW }
1344 #define NFILE asize(fnm)
1346 /* Command line options */
1347 static gmx_bool bRenum=TRUE;
1348 static gmx_bool bRmVSBds=TRUE,bZero=FALSE;
1349 static int i,maxwarn=0;
1350 static real fr_time=-1;
1352 { "-v", FALSE, etBOOL,{&bVerbose},
1353 "Be loud and noisy" },
1354 { "-time", FALSE, etREAL, {&fr_time},
1355 "Take frame at or first after this time." },
1356 { "-rmvsbds",FALSE, etBOOL, {&bRmVSBds},
1357 "Remove constant bonded interactions with virtual sites" },
1358 { "-maxwarn", FALSE, etINT, {&maxwarn},
1359 "Number of allowed warnings during input processing. Not for normal use and may generate unstable systems" },
1360 { "-zero", FALSE, etBOOL, {&bZero},
1361 "Set parameters for bonded interactions without defaults to zero instead of generating an error" },
1362 { "-renum", FALSE, etBOOL, {&bRenum},
1363 "Renumber atomtypes and minimize number of atomtypes" }
1366 CopyRight(stderr,argv[0]);
1368 /* Initiate some variables */
1373 /* Parse the command line */
1374 parse_common_args(&argc,argv,0,NFILE,fnm,asize(pa),pa,
1375 asize(desc),desc,0,NULL,&oenv);
1377 wi = init_warning(TRUE,maxwarn);
1379 /* PARAMETER file processing */
1380 mdparin = opt2fn("-f",NFILE,fnm);
1381 set_warning_line(wi,mdparin,-1);
1382 get_ir(mdparin,opt2fn("-po",NFILE,fnm),ir,opts,wi);
1385 fprintf(stderr,"checking input for internal consistency...\n");
1386 check_ir(mdparin,ir,opts,wi);
1388 if (ir->ld_seed == -1) {
1389 ir->ld_seed = make_seed();
1390 fprintf(stderr,"Setting the LD random seed to %d\n",ir->ld_seed);
1393 if (ir->expandedvals->lmc_seed == -1) {
1394 ir->expandedvals->lmc_seed = make_seed();
1395 fprintf(stderr,"Setting the lambda MC random seed to %d\n",ir->expandedvals->lmc_seed);
1398 bNeedVel = EI_STATE_VELOCITY(ir->eI);
1399 bGenVel = (bNeedVel && opts->bGenVel);
1400 if (bGenVel && ir->bContinuation)
1403 "Generating velocities is inconsistent with attempting "
1404 "to continue a previous run. Choose only one of "
1405 "gen-vel = yes and continuation = yes.");
1406 warning_error(wi, warn_buf);
1412 atype = init_atomtype();
1414 pr_symtab(debug,0,"Just opened",&sys->symtab);
1416 strcpy(fn,ftp2fn(efTOP,NFILE,fnm));
1417 if (!gmx_fexist(fn))
1418 gmx_fatal(FARGS,"%s does not exist",fn);
1419 new_status(fn,opt2fn_null("-pp",NFILE,fnm),opt2fn("-c",NFILE,fnm),
1420 opts,ir,bZero,bGenVel,bVerbose,&state,
1421 atype,sys,&nmi,&mi,plist,&comb,&reppow,&fudgeQQ,
1426 pr_symtab(debug,0,"After new_status",&sys->symtab);
1428 if (ir->cutoff_scheme == ecutsVERLET)
1430 fprintf(stderr,"Removing all charge groups because cutoff-scheme=%s\n",
1431 ecutscheme_names[ir->cutoff_scheme]);
1433 /* Remove all charge groups */
1434 gmx_mtop_remove_chargegroups(sys);
1437 if (count_constraints(sys,mi,wi) && (ir->eConstrAlg == econtSHAKE)) {
1438 if (ir->eI == eiCG || ir->eI == eiLBFGS) {
1439 sprintf(warn_buf,"Can not do %s with %s, use %s",
1440 EI(ir->eI),econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1441 warning_error(wi,warn_buf);
1443 if (ir->bPeriodicMols) {
1444 sprintf(warn_buf,"Can not do periodic molecules with %s, use %s",
1445 econstr_names[econtSHAKE],econstr_names[econtLINCS]);
1446 warning_error(wi,warn_buf);
1450 if ( EI_SD (ir->eI) && ir->etc != etcNO ) {
1451 warning_note(wi,"Temperature coupling is ignored with SD integrators.");
1454 /* If we are doing QM/MM, check that we got the atom numbers */
1455 have_atomnumber = TRUE;
1456 for (i=0; i<get_atomtype_ntypes(atype); i++) {
1457 have_atomnumber = have_atomnumber && (get_atomtype_atomnumber(i,atype) >= 0);
1459 if (!have_atomnumber && ir->bQMMM)
1463 "It appears as if you are trying to run a QM/MM calculation, but the force\n"
1464 "field you are using does not contain atom numbers fields. This is an\n"
1465 "optional field (introduced in Gromacs 3.3) for general runs, but mandatory\n"
1466 "for QM/MM. The good news is that it is easy to add - put the atom number as\n"
1467 "an integer just before the mass column in ffXXXnb.itp.\n"
1468 "NB: United atoms have the same atom numbers as normal ones.\n\n");
1472 if ((ir->adress->const_wf>1) || (ir->adress->const_wf<0)) {
1473 warning_error(wi,"AdResS contant weighting function should be between 0 and 1\n\n");
1475 /** \TODO check size of ex+hy width against box size */
1478 /* Check for errors in the input now, since they might cause problems
1479 * during processing further down.
1481 check_warning_error(wi,FARGS);
1483 if (opt2bSet("-r",NFILE,fnm))
1484 sprintf(fn,"%s",opt2fn("-r",NFILE,fnm));
1486 sprintf(fn,"%s",opt2fn("-c",NFILE,fnm));
1487 if (opt2bSet("-rb",NFILE,fnm))
1488 sprintf(fnB,"%s",opt2fn("-rb",NFILE,fnm));
1492 if (nint_ftype(sys,mi,F_POSRES) > 0)
1496 fprintf(stderr,"Reading position restraint coords from %s",fn);
1497 if (strcmp(fn,fnB) == 0)
1499 fprintf(stderr,"\n");
1503 fprintf(stderr," and %s\n",fnB);
1506 gen_posres(sys,mi,fn,fnB,
1507 ir->refcoord_scaling,ir->ePBC,
1508 ir->posres_com,ir->posres_comB,
1513 /* set parameters for virtual site construction (not for vsiten) */
1514 for(mt=0; mt<sys->nmoltype; mt++) {
1516 set_vsites(bVerbose, &sys->moltype[mt].atoms, atype, mi[mt].plist);
1518 /* now throw away all obsolete bonds, angles and dihedrals: */
1519 /* note: constraints are ALWAYS removed */
1521 for(mt=0; mt<sys->nmoltype; mt++) {
1522 clean_vsite_bondeds(mi[mt].plist,sys->moltype[mt].atoms.nr,bRmVSBds);
1526 /* If we are using CMAP, setup the pre-interpolation grid */
1529 init_cmap_grid(&sys->ffparams.cmap_grid, plist->nc, plist->grid_spacing);
1530 setup_cmap(plist->grid_spacing, plist->nc, plist->cmap,&sys->ffparams.cmap_grid);
1533 set_wall_atomtype(atype,opts,ir,wi);
1535 renum_atype(plist, sys, ir->wall_atomtype, atype, bVerbose);
1536 ntype = get_atomtype_ntypes(atype);
1539 if (ir->implicit_solvent != eisNO)
1541 /* Now we have renumbered the atom types, we can check the GBSA params */
1542 check_gbsa_params(ir,atype);
1544 /* Check that all atoms that have charge and/or LJ-parameters also have
1545 * sensible GB-parameters
1547 check_gbsa_params_charged(sys,atype);
1550 /* PELA: Copy the atomtype data to the topology atomtype list */
1551 copy_atomtype_atomtypes(atype,&(sys->atomtypes));
1554 pr_symtab(debug,0,"After renum_atype",&sys->symtab);
1557 fprintf(stderr,"converting bonded parameters...\n");
1559 ntype = get_atomtype_ntypes(atype);
1560 convert_params(ntype, plist, mi, comb, reppow, fudgeQQ, sys);
1563 pr_symtab(debug,0,"After convert_params",&sys->symtab);
1565 /* set ptype to VSite for virtual sites */
1566 for(mt=0; mt<sys->nmoltype; mt++) {
1567 set_vsites_ptype(FALSE,&sys->moltype[mt]);
1570 pr_symtab(debug,0,"After virtual sites",&sys->symtab);
1572 /* Check velocity for virtual sites and shells */
1574 check_vel(sys,state.v);
1580 for(i=0; i<sys->nmoltype; i++) {
1581 check_cg_sizes(ftp2fn(efTOP,NFILE,fnm),&sys->moltype[i].cgs,wi);
1584 if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
1586 check_bonds_timestep(sys,ir->delta_t,wi);
1589 if (EI_ENERGY_MINIMIZATION(ir->eI) && 0 == ir->nsteps)
1591 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.");
1594 check_warning_error(wi,FARGS);
1597 fprintf(stderr,"initialising group options...\n");
1598 do_index(mdparin,ftp2fn_null(efNDX,NFILE,fnm),
1600 bGenVel ? state.v : NULL,
1603 if (ir->cutoff_scheme == ecutsVERLET && ir->verletbuf_drift > 0 &&
1606 if (EI_DYNAMICS(ir->eI) &&
1607 !(EI_MD(ir->eI) && ir->etc==etcNO) &&
1608 inputrec2nboundeddim(ir) == 3)
1610 set_verlet_buffer(sys,ir,state.box,ir->verletbuf_drift,wi);
1614 /* Init the temperature coupling state */
1615 init_gtc_state(&state,ir->opts.ngtc,0,ir->opts.nhchainlength); /* need to add nnhpres here? */
1618 fprintf(stderr,"Checking consistency between energy and charge groups...\n");
1619 check_eg_vs_cg(sys);
1622 pr_symtab(debug,0,"After index",&sys->symtab);
1623 triple_check(mdparin,ir,sys,wi);
1624 close_symtab(&sys->symtab);
1626 pr_symtab(debug,0,"After close",&sys->symtab);
1628 /* make exclusions between QM atoms */
1630 if (ir->QMMMscheme==eQMMMschemenormal && ir->ns_type == ensSIMPLE ){
1631 gmx_fatal(FARGS,"electrostatic embedding only works with grid neighboursearching, use ns-type=grid instead\n");
1634 generate_qmexcl(sys,ir,wi);
1638 if (ftp2bSet(efTRN,NFILE,fnm)) {
1640 fprintf(stderr,"getting data from old trajectory ...\n");
1641 cont_status(ftp2fn(efTRN,NFILE,fnm),ftp2fn_null(efEDR,NFILE,fnm),
1642 bNeedVel,bGenVel,fr_time,ir,&state,sys,oenv);
1645 if (ir->ePBC==epbcXY && ir->nwall!=2)
1647 clear_rvec(state.box[ZZ]);
1650 if (ir->cutoff_scheme != ecutsVERLET && ir->rlist > 0)
1652 set_warning_line(wi,mdparin,-1);
1653 check_chargegroup_radii(sys,ir,state.x,wi);
1656 if (EEL_FULL(ir->coulombtype)) {
1657 /* Calculate the optimal grid dimensions */
1658 copy_mat(state.box,box);
1659 if (ir->ePBC==epbcXY && ir->nwall==2)
1660 svmul(ir->wall_ewald_zfac,box[ZZ],box[ZZ]);
1661 if (ir->nkx > 0 && ir->nky > 0 && ir->nkz > 0)
1663 /* Mark fourier_spacing as not used */
1664 ir->fourier_spacing = 0;
1666 else if (ir->nkx != 0 && ir->nky != 0 && ir->nkz != 0)
1668 set_warning_line(wi,mdparin,-1);
1669 warning_error(wi,"Some of the Fourier grid sizes are set, but all of them need to be set.");
1671 max_spacing = calc_grid(stdout,box,ir->fourier_spacing,
1672 &(ir->nkx),&(ir->nky),&(ir->nkz));
1675 if (ir->ePull != epullNO)
1676 set_pull_init(ir,sys,state.x,state.box,oenv,opts->pull_start);
1680 set_reference_positions(ir->rot,sys,state.x,state.box,
1681 opt2fn("-ref",NFILE,fnm),opt2bSet("-ref",NFILE,fnm),
1685 /* reset_multinr(sys); */
1687 if (EEL_PME(ir->coulombtype)) {
1688 float ratio = pme_load_estimate(sys,ir,state.box);
1689 fprintf(stderr,"Estimate for the relative computational load of the PME mesh part: %.2f\n",ratio);
1690 /* With free energy we might need to do PME both for the A and B state
1691 * charges. This will double the cost, but the optimal performance will
1692 * then probably be at a slightly larger cut-off and grid spacing.
1694 if ((ir->efep == efepNO && ratio > 1.0/2.0) ||
1695 (ir->efep != efepNO && ratio > 2.0/3.0)) {
1697 "The optimal PME mesh load for parallel simulations is below 0.5\n"
1698 "and for highly parallel simulations between 0.25 and 0.33,\n"
1699 "for higher performance, increase the cut-off and the PME grid spacing.\n");
1700 if (ir->efep != efepNO) {
1702 "For free energy simulations, the optimal load limit increases from 0.5 to 0.667\n");
1708 char warn_buf[STRLEN];
1709 double cio = compute_io(ir,sys->natoms,&sys->groups,F_NRE,1);
1710 sprintf(warn_buf,"This run will generate roughly %.0f Mb of data",cio);
1712 set_warning_line(wi,mdparin,-1);
1713 warning_note(wi,warn_buf);
1715 printf("%s\n",warn_buf);
1719 /* MRS: eventually figure out better logic for initializing the fep
1720 values that makes declaring the lambda and declaring the state not
1721 potentially conflict if not handled correctly. */
1722 if (ir->efep != efepNO)
1724 state.fep_state = ir->fepvals->init_fep_state;
1725 for (i=0;i<efptNR;i++)
1727 /* init_lambda trumps state definitions*/
1728 if (ir->fepvals->init_lambda >= 0)
1730 state.lambda[i] = ir->fepvals->init_lambda;
1734 if (ir->fepvals->all_lambda[i] == NULL)
1736 gmx_fatal(FARGS,"Values of lambda not set for a free energy calculation!");
1740 state.lambda[i] = ir->fepvals->all_lambda[i][state.fep_state];
1747 fprintf(stderr,"writing run input file...\n");
1749 done_warning(wi,FARGS);
1751 write_tpx_state(ftp2fn(efTPX,NFILE,fnm),ir,&state,sys);