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.
53 #include "chargegroup.h"
58 #include "gmx_fatal.h"
73 #include "gmx_random.h"
78 #include "gmx_wallcycle.h"
79 #include "mtop_util.h"
85 #include "gmx_x86_sse2.h"
89 static void global_max(t_commrec *cr,int *n)
95 gmx_sumi(cr->nnodes,sum,cr);
96 for(i=0; i<cr->nnodes; i++)
102 static void realloc_bins(double **bin,int *nbin,int nbin_new)
106 if (nbin_new != *nbin) {
107 srenew(*bin,nbin_new);
108 for(i=*nbin; i<nbin_new; i++)
114 double do_tpi(FILE *fplog,t_commrec *cr,
115 int nfile, const t_filenm fnm[],
116 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
118 gmx_vsite_t *vsite,gmx_constr_t constr,
120 t_inputrec *inputrec,
121 gmx_mtop_t *top_global,t_fcdata *fcd,
124 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
127 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
129 real cpt_period,real max_hours,
130 const char *deviceOptions,
132 gmx_runtime_t *runtime)
134 const char *TPI="Test Particle Insertion";
136 gmx_groups_t *groups;
137 gmx_enerdata_t *enerd;
139 real lambda,t,temp,beta,drmax,epot;
140 double embU,sum_embU,*sum_UgembU,V,V_all,VembU_all;
143 gmx_bool bDispCorr,bCharge,bRFExcl,bNotLastFrame,bStateChanged,bNS,bOurStep;
144 tensor force_vir,shake_vir,vir,pres;
145 int cg_tp,a_tp0,a_tp1,ngid,gid_tp,nener,e;
147 rvec mu_tot,x_init,dx,x_tp;
148 int nnodes,frame,nsteps,step;
152 char *ptr,*dump_pdb,**leg,str[STRLEN],str2[STRLEN];
153 double dbl,dump_ener;
156 real *mass_cavity=NULL,mass_tot;
158 double invbinw,*bin,refvolshift,logV,bUlogV;
159 real dvdl,prescorr,enercorr,dvdlcorr;
160 gmx_bool bEnergyOutOfBounds;
161 const char *tpid_leg[2]={"direct","reweighted"};
163 /* Since there is no upper limit to the insertion energies,
164 * we need to set an upper limit for the distribution output.
166 real bU_bin_limit = 50;
167 real bU_logV_bin_limit = bU_bin_limit + 10;
171 top = gmx_mtop_generate_local_top(top_global,inputrec);
173 groups = &top_global->groups;
175 bCavity = (inputrec->eI == eiTPIC);
177 ptr = getenv("GMX_TPIC_MASSES");
181 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
182 * The center of mass of the last atoms is then used for TPIC.
185 while (sscanf(ptr,"%lf%n",&dbl,&i) > 0) {
186 srenew(mass_cavity,nat_cavity+1);
187 mass_cavity[nat_cavity] = dbl;
188 fprintf(fplog,"mass[%d] = %f\n",
189 nat_cavity+1,mass_cavity[nat_cavity]);
194 gmx_fatal(FARGS,"Found %d masses in GMX_TPIC_MASSES",nat_cavity);
199 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
200 state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
201 /* We never need full pbc for TPI */
203 /* Determine the temperature for the Boltzmann weighting */
204 temp = inputrec->opts.ref_t[0];
206 for(i=1; (i<inputrec->opts.ngtc); i++) {
207 if (inputrec->opts.ref_t[i] != temp) {
208 fprintf(fplog,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
209 fprintf(stderr,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
213 "\n The temperature for test particle insertion is %.3f K\n\n",
216 beta = 1.0/(BOLTZ*temp);
218 /* Number of insertions per frame */
219 nsteps = inputrec->nsteps;
221 /* Use the same neighborlist with more insertions points
222 * in a sphere of radius drmax around the initial point
224 /* This should be a proper mdp parameter */
225 drmax = inputrec->rtpi;
227 /* An environment variable can be set to dump all configurations
228 * to pdb with an insertion energy <= this value.
230 dump_pdb = getenv("GMX_TPI_DUMP");
233 sscanf(dump_pdb,"%lf",&dump_ener);
235 atoms2md(top_global,inputrec,0,NULL,0,top_global->natoms,mdatoms);
236 update_mdatoms(mdatoms,inputrec->fepvals->init_lambda);
239 init_enerdata(groups->grps[egcENER].nr,inputrec->fepvals->n_lambda,enerd);
240 snew(f,top_global->natoms);
242 /* Print to log file */
243 runtime_start(runtime);
244 print_date_and_time(fplog,cr->nodeid,
245 "Started Test Particle Insertion",runtime);
246 wallcycle_start(wcycle,ewcRUN);
248 /* The last charge group is the group to be inserted */
249 cg_tp = top->cgs.nr - 1;
250 a_tp0 = top->cgs.index[cg_tp];
251 a_tp1 = top->cgs.index[cg_tp+1];
253 fprintf(debug,"TPI cg %d, atoms %d-%d\n",cg_tp,a_tp0,a_tp1);
254 if (a_tp1 - a_tp0 > 1 &&
255 (inputrec->rlist < inputrec->rcoulomb ||
256 inputrec->rlist < inputrec->rvdw))
257 gmx_fatal(FARGS,"Can not do TPI for multi-atom molecule with a twin-range cut-off");
258 snew(x_mol,a_tp1-a_tp0);
260 bDispCorr = (inputrec->eDispCorr != edispcNO);
262 for(i=a_tp0; i<a_tp1; i++) {
263 /* Copy the coordinates of the molecule to be insterted */
264 copy_rvec(state->x[i],x_mol[i-a_tp0]);
265 /* Check if we need to print electrostatic energies */
266 bCharge |= (mdatoms->chargeA[i] != 0 ||
267 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
269 bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype!=eelRF_NEC);
271 calc_cgcm(fplog,cg_tp,cg_tp+1,&(top->cgs),state->x,fr->cg_cm);
273 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog) {
274 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
275 fprintf(stderr,"WARNING: Your TPI molecule is not centered at 0,0,0\n");
278 /* Center the molecule to be inserted at zero */
279 for(i=0; i<a_tp1-a_tp0; i++)
280 rvec_dec(x_mol[i],fr->cg_cm[cg_tp]);
284 fprintf(fplog,"\nWill insert %d atoms %s partial charges\n",
285 a_tp1-a_tp0,bCharge ? "with" : "without");
287 fprintf(fplog,"\nWill insert %d times in each frame of %s\n",
288 nsteps,opt2fn("-rerun",nfile,fnm));
293 if (inputrec->nstlist > 1)
295 if (drmax==0 && a_tp1-a_tp0==1)
297 gmx_fatal(FARGS,"Re-using the neighborlist %d times for insertions of a single atom in a sphere of radius %f does not make sense",inputrec->nstlist,drmax);
301 fprintf(fplog,"Will use the same neighborlist for %d insertions in a sphere of radius %f\n",inputrec->nstlist,drmax);
309 fprintf(fplog,"Will insert randomly in a sphere of radius %f around the center of the cavity\n",drmax);
313 ngid = groups->grps[egcENER].nr;
314 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
322 if (EEL_FULL(fr->eeltype))
325 snew(sum_UgembU,nener);
327 /* Initialize random generator */
328 tpi_rand = gmx_rng_init(inputrec->ld_seed);
331 fp_tpi = xvgropen(opt2fn("-tpi",nfile,fnm),
332 "TPI energies","Time (ps)",
333 "(kJ mol\\S-1\\N) / (nm\\S3\\N)",oenv);
334 xvgr_subtitle(fp_tpi,"f. are averages over one frame",oenv);
337 sprintf(str,"-kT log(<Ve\\S-\\betaU\\N>/<V>)");
338 leg[e++] = strdup(str);
339 sprintf(str,"f. -kT log<e\\S-\\betaU\\N>");
340 leg[e++] = strdup(str);
341 sprintf(str,"f. <e\\S-\\betaU\\N>");
342 leg[e++] = strdup(str);
344 leg[e++] = strdup(str);
345 sprintf(str,"f. <Ue\\S-\\betaU\\N>");
346 leg[e++] = strdup(str);
347 for(i=0; i<ngid; i++) {
348 sprintf(str,"f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
349 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
350 leg[e++] = strdup(str);
353 sprintf(str,"f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
354 leg[e++] = strdup(str);
357 for(i=0; i<ngid; i++) {
358 sprintf(str,"f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
359 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
360 leg[e++] = strdup(str);
363 sprintf(str,"f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
364 leg[e++] = strdup(str);
366 if (EEL_FULL(fr->eeltype)) {
367 sprintf(str,"f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
368 leg[e++] = strdup(str);
371 xvgr_legend(fp_tpi,4+nener,(const char**)leg,oenv);
372 for(i=0; i<4+nener; i++)
384 bNotLastFrame = read_first_frame(oenv,&status,opt2fn("-rerun",nfile,fnm),
385 &rerun_fr,TRX_NEED_X);
388 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
389 mdatoms->nr - (a_tp1 - a_tp0))
390 gmx_fatal(FARGS,"Number of atoms in trajectory (%d)%s "
391 "is not equal the number in the run input file (%d) "
392 "minus the number of atoms to insert (%d)\n",
393 rerun_fr.natoms,bCavity ? " minus one" : "",
394 mdatoms->nr,a_tp1-a_tp0);
396 refvolshift = log(det(rerun_fr.box));
399 /* Make sure we don't detect SSE overflow generated before this point */
400 gmx_mm_check_and_reset_overflow();
403 while (bNotLastFrame)
405 lambda = rerun_fr.lambda;
409 for(e=0; e<nener; e++)
414 /* Copy the coordinates from the input trajectory */
415 for(i=0; i<rerun_fr.natoms; i++)
417 copy_rvec(rerun_fr.x[i],state->x[i]);
419 copy_mat(rerun_fr.box,state->box);
424 bStateChanged = TRUE;
426 for(step=0; step<nsteps; step++)
428 /* In parallel all nodes generate all random configurations.
429 * In that way the result is identical to a single cpu tpi run.
433 /* Random insertion in the whole volume */
434 bNS = (step % inputrec->nstlist == 0);
437 /* Generate a random position in the box */
438 x_init[XX] = gmx_rng_uniform_real(tpi_rand)*state->box[XX][XX];
439 x_init[YY] = gmx_rng_uniform_real(tpi_rand)*state->box[YY][YY];
440 x_init[ZZ] = gmx_rng_uniform_real(tpi_rand)*state->box[ZZ][ZZ];
442 if (inputrec->nstlist == 1)
444 copy_rvec(x_init,x_tp);
448 /* Generate coordinates within |dx|=drmax of x_init */
451 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
452 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
453 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
455 while (norm2(dx) > drmax*drmax);
456 rvec_add(x_init,dx,x_tp);
461 /* Random insertion around a cavity location
462 * given by the last coordinate of the trajectory.
468 /* Copy the location of the cavity */
469 copy_rvec(rerun_fr.x[rerun_fr.natoms-1],x_init);
473 /* Determine the center of mass of the last molecule */
476 for(i=0; i<nat_cavity; i++)
481 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
483 mass_tot += mass_cavity[i];
487 x_init[d] /= mass_tot;
491 /* Generate coordinates within |dx|=drmax of x_init */
494 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
495 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
496 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
498 while (norm2(dx) > drmax*drmax);
499 rvec_add(x_init,dx,x_tp);
502 if (a_tp1 - a_tp0 == 1)
504 /* Insert a single atom, just copy the insertion location */
505 copy_rvec(x_tp,state->x[a_tp0]);
509 /* Copy the coordinates from the top file */
510 for(i=a_tp0; i<a_tp1; i++)
512 copy_rvec(x_mol[i-a_tp0],state->x[i]);
514 /* Rotate the molecule randomly */
515 rotate_conf(a_tp1-a_tp0,state->x+a_tp0,NULL,
516 2*M_PI*gmx_rng_uniform_real(tpi_rand),
517 2*M_PI*gmx_rng_uniform_real(tpi_rand),
518 2*M_PI*gmx_rng_uniform_real(tpi_rand));
519 /* Shift to the insertion location */
520 for(i=a_tp0; i<a_tp1; i++)
522 rvec_inc(state->x[i],x_tp);
526 /* Check if this insertion belongs to this node */
530 switch (inputrec->eI)
533 bOurStep = ((step / inputrec->nstlist) % nnodes == cr->nodeid);
536 bOurStep = (step % nnodes == cr->nodeid);
539 gmx_fatal(FARGS,"Unknown integrator %s",ei_names[inputrec->eI]);
544 /* Clear some matrix variables */
545 clear_mat(force_vir);
546 clear_mat(shake_vir);
550 /* Set the charge group center of mass of the test particle */
551 copy_rvec(x_init,fr->cg_cm[top->cgs.nr-1]);
553 /* Calc energy (no forces) on new positions.
554 * Since we only need the intermolecular energy
555 * and the RF exclusion terms of the inserted molecule occur
556 * within a single charge group we can pass NULL for the graph.
557 * This also avoids shifts that would move charge groups
560 * Some checks above ensure than we can not have
561 * twin-range interactions together with nstlist > 1,
562 * therefore we do not need to remember the LR energies.
564 /* Make do_force do a single node force calculation */
566 do_force(fplog,cr,inputrec,
567 step,nrnb,wcycle,top,top_global,&top_global->groups,
568 state->box,state->x,&state->hist,
569 f,force_vir,mdatoms,enerd,fcd,
571 NULL,fr,NULL,mu_tot,t,NULL,NULL,FALSE,
572 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
573 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS | GMX_FORCE_DO_LR : 0) |
574 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
576 bStateChanged = FALSE;
579 /* Calculate long range corrections to pressure and energy */
580 calc_dispcorr(fplog,inputrec,fr,step,top_global->natoms,state->box,
581 lambda,pres,vir,&prescorr,&enercorr,&dvdlcorr);
582 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
583 enerd->term[F_DISPCORR] = enercorr;
584 enerd->term[F_EPOT] += enercorr;
585 enerd->term[F_PRES] += prescorr;
586 enerd->term[F_DVDL_VDW] += dvdlcorr;
588 epot = enerd->term[F_EPOT];
589 bEnergyOutOfBounds = FALSE;
591 /* With SSE the energy can overflow, check for this */
592 if (gmx_mm_check_and_reset_overflow())
596 fprintf(debug,"Found an SSE overflow, assuming the energy is out of bounds\n");
598 bEnergyOutOfBounds = TRUE;
601 /* If the compiler doesn't optimize this check away
602 * we catch the NAN energies.
603 * The epot>GMX_REAL_MAX check catches inf values,
604 * which should nicely result in embU=0 through the exp below,
605 * but it does not hurt to check anyhow.
607 /* Non-bonded Interaction usually diverge at r=0.
608 * With tabulated interaction functions the first few entries
609 * should be capped in a consistent fashion between
610 * repulsion, dispersion and Coulomb to avoid accidental
611 * negative values in the total energy.
612 * The table generation code in tables.c does this.
613 * With user tbales the user should take care of this.
615 if (epot != epot || epot > GMX_REAL_MAX)
617 bEnergyOutOfBounds = TRUE;
619 if (bEnergyOutOfBounds)
623 fprintf(debug,"\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n",t,step,epot);
629 embU = exp(-beta*epot);
631 /* Determine the weighted energy contributions of each energy group */
633 sum_UgembU[e++] += epot*embU;
636 for(i=0; i<ngid; i++)
639 (enerd->grpp.ener[egBHAMSR][GID(i,gid_tp,ngid)] +
640 enerd->grpp.ener[egBHAMLR][GID(i,gid_tp,ngid)])*embU;
645 for(i=0; i<ngid; i++)
648 (enerd->grpp.ener[egLJSR][GID(i,gid_tp,ngid)] +
649 enerd->grpp.ener[egLJLR][GID(i,gid_tp,ngid)])*embU;
654 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
658 for(i=0; i<ngid; i++)
661 (enerd->grpp.ener[egCOULSR][GID(i,gid_tp,ngid)] +
662 enerd->grpp.ener[egCOULLR][GID(i,gid_tp,ngid)])*embU;
666 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
668 if (EEL_FULL(fr->eeltype))
670 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
675 if (embU == 0 || beta*epot > bU_bin_limit)
681 i = (int)((bU_logV_bin_limit
682 - (beta*epot - logV + refvolshift))*invbinw
690 realloc_bins(&bin,&nbin,i+10);
697 fprintf(debug,"TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
698 step,epot,x_tp[XX],x_tp[YY],x_tp[ZZ]);
701 if (dump_pdb && epot <= dump_ener)
703 sprintf(str,"t%g_step%d.pdb",t,step);
704 sprintf(str2,"t: %f step %d ener: %f",t,step,epot);
705 write_sto_conf_mtop(str,str2,top_global,state->x,state->v,
706 inputrec->ePBC,state->box);
713 /* When running in parallel sum the energies over the processes */
714 gmx_sumd(1, &sum_embU, cr);
715 gmx_sumd(nener,sum_UgembU,cr);
720 VembU_all += V*sum_embU/nsteps;
724 if (bVerbose || frame%10==0 || frame<10)
726 fprintf(stderr,"mu %10.3e <mu> %10.3e\n",
727 -log(sum_embU/nsteps)/beta,-log(VembU_all/V_all)/beta);
730 fprintf(fp_tpi,"%10.3f %12.5e %12.5e %12.5e %12.5e",
732 VembU_all==0 ? 20/beta : -log(VembU_all/V_all)/beta,
733 sum_embU==0 ? 20/beta : -log(sum_embU/nsteps)/beta,
735 for(e=0; e<nener; e++)
737 fprintf(fp_tpi," %12.5e",sum_UgembU[e]/nsteps);
739 fprintf(fp_tpi,"\n");
743 bNotLastFrame = read_next_frame(oenv, status,&rerun_fr);
744 } /* End of the loop */
745 runtime_end(runtime);
751 gmx_fio_fclose(fp_tpi);
757 fprintf(fplog," <V> = %12.5e nm^3\n",V_all/frame);
758 fprintf(fplog," <mu> = %12.5e kJ/mol\n",-log(VembU_all/V_all)/beta);
761 /* Write the Boltzmann factor histogram */
764 /* When running in parallel sum the bins over the processes */
767 realloc_bins(&bin,&nbin,i);
768 gmx_sumd(nbin,bin,cr);
772 fp_tpi = xvgropen(opt2fn("-tpid",nfile,fnm),
773 "TPI energy distribution",
774 "\\betaU - log(V/<V>)","count",oenv);
775 sprintf(str,"number \\betaU > %g: %9.3e",bU_bin_limit,bin[0]);
776 xvgr_subtitle(fp_tpi,str,oenv);
777 xvgr_legend(fp_tpi,2,(const char **)tpid_leg,oenv);
778 for(i=nbin-1; i>0; i--)
780 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
781 fprintf(fp_tpi,"%6.2f %10d %12.5e\n",
784 bin[i]*exp(-bUlogV)*V_all/VembU_all);
786 gmx_fio_fclose(fp_tpi);
792 runtime->nsteps_done = frame*inputrec->nsteps;