Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / mdlib / tpi.c
blob2a5ae3d18e828df7208bc2ab40c26de121538a0e
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
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
33 * And Hey:
34 * GROwing Monsters And Cloning Shrimps
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <string.h>
41 #include <time.h>
42 #include <math.h>
43 #include "sysstuff.h"
44 #include "string2.h"
45 #include "network.h"
46 #include "confio.h"
47 #include "copyrite.h"
48 #include "smalloc.h"
49 #include "nrnb.h"
50 #include "main.h"
51 #include "chargegroup.h"
52 #include "force.h"
53 #include "macros.h"
54 #include "random.h"
55 #include "names.h"
56 #include "gmx_fatal.h"
57 #include "txtdump.h"
58 #include "typedefs.h"
59 #include "update.h"
60 #include "random.h"
61 #include "constr.h"
62 #include "vec.h"
63 #include "statutil.h"
64 #include "tgroup.h"
65 #include "mdebin.h"
66 #include "vsite.h"
67 #include "force.h"
68 #include "mdrun.h"
69 #include "domdec.h"
70 #include "partdec.h"
71 #include "gmx_random.h"
72 #include "physics.h"
73 #include "xvgr.h"
74 #include "mdatoms.h"
75 #include "ns.h"
76 #include "gmx_wallcycle.h"
77 #include "mtop_util.h"
78 #include "gmxfio.h"
79 #include "pme.h"
80 #include "gbutil.h"
82 #if ( defined(GMX_IA32_SSE) || defined(GMX_X86_64_SSE) || defined(GMX_X86_64_SSE2) )
83 #if defined(GMX_DOUBLE)
84 #include "gmx_sse2_double.h"
85 #else
86 #include "gmx_sse2_single.h"
87 #endif
88 #endif
91 static void global_max(t_commrec *cr,int *n)
93 int *sum,i;
95 snew(sum,cr->nnodes);
96 sum[cr->nodeid] = *n;
97 gmx_sumi(cr->nnodes,sum,cr);
98 for(i=0; i<cr->nnodes; i++)
99 *n = max(*n,sum[i]);
101 sfree(sum);
104 static void realloc_bins(double **bin,int *nbin,int nbin_new)
106 int i;
108 if (nbin_new != *nbin) {
109 srenew(*bin,nbin_new);
110 for(i=*nbin; i<nbin_new; i++)
111 (*bin)[i] = 0;
112 *nbin = nbin_new;
116 double do_tpi(FILE *fplog,t_commrec *cr,
117 int nfile, const t_filenm fnm[],
118 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
119 int nstglobalcomm,
120 gmx_vsite_t *vsite,gmx_constr_t constr,
121 int stepout,
122 t_inputrec *inputrec,
123 gmx_mtop_t *top_global,t_fcdata *fcd,
124 t_state *state,
125 t_mdatoms *mdatoms,
126 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
127 gmx_edsam_t ed,
128 t_forcerec *fr,
129 int repl_ex_nst,int repl_ex_seed,
130 gmx_membed_t *membed,
131 real cpt_period,real max_hours,
132 const char *deviceOptions,
133 unsigned long Flags,
134 gmx_runtime_t *runtime)
136 const char *TPI="Test Particle Insertion";
137 gmx_localtop_t *top;
138 gmx_groups_t *groups;
139 gmx_enerdata_t *enerd;
140 rvec *f;
141 real lambda,t,temp,beta,drmax,epot;
142 double embU,sum_embU,*sum_UgembU,V,V_all,VembU_all;
143 t_trxstatus *status;
144 t_trxframe rerun_fr;
145 gmx_bool bDispCorr,bCharge,bRFExcl,bNotLastFrame,bStateChanged,bNS,bOurStep;
146 tensor force_vir,shake_vir,vir,pres;
147 int cg_tp,a_tp0,a_tp1,ngid,gid_tp,nener,e;
148 rvec *x_mol;
149 rvec mu_tot,x_init,dx,x_tp;
150 int nnodes,frame,nsteps,step;
151 int i,start,end;
152 gmx_rng_t tpi_rand;
153 FILE *fp_tpi=NULL;
154 char *ptr,*dump_pdb,**leg,str[STRLEN],str2[STRLEN];
155 double dbl,dump_ener;
156 gmx_bool bCavity;
157 int nat_cavity=0,d;
158 real *mass_cavity=NULL,mass_tot;
159 int nbin;
160 double invbinw,*bin,refvolshift,logV,bUlogV;
161 real dvdl,prescorr,enercorr,dvdlcorr;
162 gmx_bool bEnergyOutOfBounds;
163 const char *tpid_leg[2]={"direct","reweighted"};
165 /* Since there is no upper limit to the insertion energies,
166 * we need to set an upper limit for the distribution output.
168 real bU_bin_limit = 50;
169 real bU_logV_bin_limit = bU_bin_limit + 10;
171 nnodes = cr->nnodes;
173 top = gmx_mtop_generate_local_top(top_global,inputrec);
175 groups = &top_global->groups;
177 bCavity = (inputrec->eI == eiTPIC);
178 if (bCavity) {
179 ptr = getenv("GMX_TPIC_MASSES");
180 if (ptr == NULL) {
181 nat_cavity = 1;
182 } else {
183 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
184 * The center of mass of the last atoms is then used for TPIC.
186 nat_cavity = 0;
187 while (sscanf(ptr,"%lf%n",&dbl,&i) > 0) {
188 srenew(mass_cavity,nat_cavity+1);
189 mass_cavity[nat_cavity] = dbl;
190 fprintf(fplog,"mass[%d] = %f\n",
191 nat_cavity+1,mass_cavity[nat_cavity]);
192 nat_cavity++;
193 ptr += i;
195 if (nat_cavity == 0)
196 gmx_fatal(FARGS,"Found %d masses in GMX_TPIC_MASSES",nat_cavity);
201 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
202 state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
203 /* We never need full pbc for TPI */
204 fr->ePBC = epbcXYZ;
205 /* Determine the temperature for the Boltzmann weighting */
206 temp = inputrec->opts.ref_t[0];
207 if (fplog) {
208 for(i=1; (i<inputrec->opts.ngtc); i++) {
209 if (inputrec->opts.ref_t[i] != temp) {
210 fprintf(fplog,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
211 fprintf(stderr,"\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
214 fprintf(fplog,
215 "\n The temperature for test particle insertion is %.3f K\n\n",
216 temp);
218 beta = 1.0/(BOLTZ*temp);
220 /* Number of insertions per frame */
221 nsteps = inputrec->nsteps;
223 /* Use the same neighborlist with more insertions points
224 * in a sphere of radius drmax around the initial point
226 /* This should be a proper mdp parameter */
227 drmax = inputrec->rtpi;
229 /* An environment variable can be set to dump all configurations
230 * to pdb with an insertion energy <= this value.
232 dump_pdb = getenv("GMX_TPI_DUMP");
233 dump_ener = 0;
234 if (dump_pdb)
235 sscanf(dump_pdb,"%lf",&dump_ener);
237 atoms2md(top_global,inputrec,0,NULL,0,top_global->natoms,mdatoms);
238 update_mdatoms(mdatoms,inputrec->init_lambda);
240 snew(enerd,1);
241 init_enerdata(groups->grps[egcENER].nr,inputrec->n_flambda,enerd);
242 snew(f,top_global->natoms);
244 /* Print to log file */
245 runtime_start(runtime);
246 print_date_and_time(fplog,cr->nodeid,
247 "Started Test Particle Insertion",runtime);
248 wallcycle_start(wcycle,ewcRUN);
250 /* The last charge group is the group to be inserted */
251 cg_tp = top->cgs.nr - 1;
252 a_tp0 = top->cgs.index[cg_tp];
253 a_tp1 = top->cgs.index[cg_tp+1];
254 if (debug)
255 fprintf(debug,"TPI cg %d, atoms %d-%d\n",cg_tp,a_tp0,a_tp1);
256 if (a_tp1 - a_tp0 > 1 &&
257 (inputrec->rlist < inputrec->rcoulomb ||
258 inputrec->rlist < inputrec->rvdw))
259 gmx_fatal(FARGS,"Can not do TPI for multi-atom molecule with a twin-range cut-off");
260 snew(x_mol,a_tp1-a_tp0);
262 bDispCorr = (inputrec->eDispCorr != edispcNO);
263 bCharge = FALSE;
264 for(i=a_tp0; i<a_tp1; i++) {
265 /* Copy the coordinates of the molecule to be insterted */
266 copy_rvec(state->x[i],x_mol[i-a_tp0]);
267 /* Check if we need to print electrostatic energies */
268 bCharge |= (mdatoms->chargeA[i] != 0 ||
269 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
271 bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype!=eelRF_NEC);
273 calc_cgcm(fplog,cg_tp,cg_tp+1,&(top->cgs),state->x,fr->cg_cm);
274 if (bCavity) {
275 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog) {
276 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
277 fprintf(stderr,"WARNING: Your TPI molecule is not centered at 0,0,0\n");
279 } else {
280 /* Center the molecule to be inserted at zero */
281 for(i=0; i<a_tp1-a_tp0; i++)
282 rvec_dec(x_mol[i],fr->cg_cm[cg_tp]);
285 if (fplog) {
286 fprintf(fplog,"\nWill insert %d atoms %s partial charges\n",
287 a_tp1-a_tp0,bCharge ? "with" : "without");
289 fprintf(fplog,"\nWill insert %d times in each frame of %s\n",
290 nsteps,opt2fn("-rerun",nfile,fnm));
293 if (!bCavity)
295 if (inputrec->nstlist > 1)
297 if (drmax==0 && a_tp1-a_tp0==1)
299 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 if (fplog)
303 fprintf(fplog,"Will use the same neighborlist for %d insertions in a sphere of radius %f\n",inputrec->nstlist,drmax);
307 else
309 if (fplog)
311 fprintf(fplog,"Will insert randomly in a sphere of radius %f around the center of the cavity\n",drmax);
315 ngid = groups->grps[egcENER].nr;
316 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
317 nener = 1 + ngid;
318 if (bDispCorr)
319 nener += 1;
320 if (bCharge) {
321 nener += ngid;
322 if (bRFExcl)
323 nener += 1;
324 if (EEL_FULL(fr->eeltype))
325 nener += 1;
327 snew(sum_UgembU,nener);
329 /* Initialize random generator */
330 tpi_rand = gmx_rng_init(inputrec->ld_seed);
332 if (MASTER(cr)) {
333 fp_tpi = xvgropen(opt2fn("-tpi",nfile,fnm),
334 "TPI energies","Time (ps)",
335 "(kJ mol\\S-1\\N) / (nm\\S3\\N)",oenv);
336 xvgr_subtitle(fp_tpi,"f. are averages over one frame",oenv);
337 snew(leg,4+nener);
338 e = 0;
339 sprintf(str,"-kT log(<Ve\\S-\\betaU\\N>/<V>)");
340 leg[e++] = strdup(str);
341 sprintf(str,"f. -kT log<e\\S-\\betaU\\N>");
342 leg[e++] = strdup(str);
343 sprintf(str,"f. <e\\S-\\betaU\\N>");
344 leg[e++] = strdup(str);
345 sprintf(str,"f. V");
346 leg[e++] = strdup(str);
347 sprintf(str,"f. <Ue\\S-\\betaU\\N>");
348 leg[e++] = strdup(str);
349 for(i=0; i<ngid; i++) {
350 sprintf(str,"f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
351 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
352 leg[e++] = strdup(str);
354 if (bDispCorr) {
355 sprintf(str,"f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
356 leg[e++] = strdup(str);
358 if (bCharge) {
359 for(i=0; i<ngid; i++) {
360 sprintf(str,"f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
361 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
362 leg[e++] = strdup(str);
364 if (bRFExcl) {
365 sprintf(str,"f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
366 leg[e++] = strdup(str);
368 if (EEL_FULL(fr->eeltype)) {
369 sprintf(str,"f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
370 leg[e++] = strdup(str);
373 xvgr_legend(fp_tpi,4+nener,(const char**)leg,oenv);
374 for(i=0; i<4+nener; i++)
375 sfree(leg[i]);
376 sfree(leg);
378 clear_rvec(x_init);
379 V_all = 0;
380 VembU_all = 0;
382 invbinw = 10;
383 nbin = 10;
384 snew(bin,nbin);
386 bNotLastFrame = read_first_frame(oenv,&status,opt2fn("-rerun",nfile,fnm),
387 &rerun_fr,TRX_NEED_X);
388 frame = 0;
390 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
391 mdatoms->nr - (a_tp1 - a_tp0))
392 gmx_fatal(FARGS,"Number of atoms in trajectory (%d)%s "
393 "is not equal the number in the run input file (%d) "
394 "minus the number of atoms to insert (%d)\n",
395 rerun_fr.natoms,bCavity ? " minus one" : "",
396 mdatoms->nr,a_tp1-a_tp0);
398 refvolshift = log(det(rerun_fr.box));
400 #if ( defined(GMX_IA32_SSE) || defined(GMX_X86_64_SSE) || defined(GMX_X86_64_SSE2) )
401 /* Make sure we don't detect SSE overflow generated before this point */
402 gmx_mm_check_and_reset_overflow();
403 #endif
405 while (bNotLastFrame)
407 lambda = rerun_fr.lambda;
408 t = rerun_fr.time;
410 sum_embU = 0;
411 for(e=0; e<nener; e++)
413 sum_UgembU[e] = 0;
416 /* Copy the coordinates from the input trajectory */
417 for(i=0; i<rerun_fr.natoms; i++)
419 copy_rvec(rerun_fr.x[i],state->x[i]);
422 V = det(rerun_fr.box);
423 logV = log(V);
425 bStateChanged = TRUE;
426 bNS = TRUE;
427 for(step=0; step<nsteps; step++)
429 /* In parallel all nodes generate all random configurations.
430 * In that way the result is identical to a single cpu tpi run.
432 if (!bCavity)
434 /* Random insertion in the whole volume */
435 bNS = (step % inputrec->nstlist == 0);
436 if (bNS)
438 /* Generate a random position in the box */
439 x_init[XX] = gmx_rng_uniform_real(tpi_rand)*state->box[XX][XX];
440 x_init[YY] = gmx_rng_uniform_real(tpi_rand)*state->box[YY][YY];
441 x_init[ZZ] = gmx_rng_uniform_real(tpi_rand)*state->box[ZZ][ZZ];
443 if (inputrec->nstlist == 1)
445 copy_rvec(x_init,x_tp);
447 else
449 /* Generate coordinates within |dx|=drmax of x_init */
452 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
453 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
454 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
456 while (norm2(dx) > drmax*drmax);
457 rvec_add(x_init,dx,x_tp);
460 else
462 /* Random insertion around a cavity location
463 * given by the last coordinate of the trajectory.
465 if (step == 0)
467 if (nat_cavity == 1)
469 /* Copy the location of the cavity */
470 copy_rvec(rerun_fr.x[rerun_fr.natoms-1],x_init);
472 else
474 /* Determine the center of mass of the last molecule */
475 clear_rvec(x_init);
476 mass_tot = 0;
477 for(i=0; i<nat_cavity; i++)
479 for(d=0; d<DIM; d++)
481 x_init[d] +=
482 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
484 mass_tot += mass_cavity[i];
486 for(d=0; d<DIM; d++)
488 x_init[d] /= mass_tot;
492 /* Generate coordinates within |dx|=drmax of x_init */
495 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
496 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
497 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
499 while (norm2(dx) > drmax*drmax);
500 rvec_add(x_init,dx,x_tp);
503 if (a_tp1 - a_tp0 == 1)
505 /* Insert a single atom, just copy the insertion location */
506 copy_rvec(x_tp,state->x[a_tp0]);
508 else
510 /* Copy the coordinates from the top file */
511 for(i=a_tp0; i<a_tp1; i++)
513 copy_rvec(x_mol[i-a_tp0],state->x[i]);
515 /* Rotate the molecule randomly */
516 rotate_conf(a_tp1-a_tp0,state->x+a_tp0,NULL,
517 2*M_PI*gmx_rng_uniform_real(tpi_rand),
518 2*M_PI*gmx_rng_uniform_real(tpi_rand),
519 2*M_PI*gmx_rng_uniform_real(tpi_rand));
520 /* Shift to the insertion location */
521 for(i=a_tp0; i<a_tp1; i++)
523 rvec_inc(state->x[i],x_tp);
527 /* Check if this insertion belongs to this node */
528 bOurStep = TRUE;
529 if (PAR(cr))
531 switch (inputrec->eI)
533 case eiTPI:
534 bOurStep = ((step / inputrec->nstlist) % nnodes == cr->nodeid);
535 break;
536 case eiTPIC:
537 bOurStep = (step % nnodes == cr->nodeid);
538 break;
539 default:
540 gmx_fatal(FARGS,"Unknown integrator %s",ei_names[inputrec->eI]);
543 if (bOurStep)
545 /* Clear some matrix variables */
546 clear_mat(force_vir);
547 clear_mat(shake_vir);
548 clear_mat(vir);
549 clear_mat(pres);
551 /* Set the charge group center of mass of the test particle */
552 copy_rvec(x_init,fr->cg_cm[top->cgs.nr-1]);
554 /* Calc energy (no forces) on new positions.
555 * Since we only need the intermolecular energy
556 * and the RF exclusion terms of the inserted molecule occur
557 * within a single charge group we can pass NULL for the graph.
558 * This also avoids shifts that would move charge groups
559 * out of the box.
561 * Some checks above ensure than we can not have
562 * twin-range interactions together with nstlist > 1,
563 * therefore we do not need to remember the LR energies.
565 /* Make do_force do a single node force calculation */
566 cr->nnodes = 1;
567 do_force(fplog,cr,inputrec,
568 step,nrnb,wcycle,top,top_global,&top_global->groups,
569 rerun_fr.box,state->x,&state->hist,
570 f,force_vir,mdatoms,enerd,fcd,
571 lambda,NULL,fr,NULL,mu_tot,t,NULL,NULL,FALSE,
572 GMX_FORCE_NONBONDED |
573 (bNS ? GMX_FORCE_NS | GMX_FORCE_DOLR : 0) |
574 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
575 cr->nnodes = nnodes;
576 bStateChanged = FALSE;
577 bNS = FALSE;
579 /* Calculate long range corrections to pressure and energy */
580 calc_dispcorr(fplog,inputrec,fr,step,top_global->natoms,rerun_fr.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] += dvdlcorr;
588 epot = enerd->term[F_EPOT];
589 bEnergyOutOfBounds = FALSE;
590 #if ( defined(GMX_IA32_SSE) || defined(GMX_X86_64_SSE) || defined(GMX_X86_64_SSE2) )
591 /* With SSE the energy can overflow, check for this */
592 if (gmx_mm_check_and_reset_overflow())
594 if (debug)
596 fprintf(debug,"Found an SSE overflow, assuming the energy is out of bounds\n");
598 bEnergyOutOfBounds = TRUE;
600 #endif
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)
621 if (debug)
623 fprintf(debug,"\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n",t,step,epot);
625 embU = 0;
627 else
629 embU = exp(-beta*epot);
630 sum_embU += embU;
631 /* Determine the weighted energy contributions of each energy group */
632 e = 0;
633 sum_UgembU[e++] += epot*embU;
634 if (fr->bBHAM)
636 for(i=0; i<ngid; i++)
638 sum_UgembU[e++] +=
639 (enerd->grpp.ener[egBHAMSR][GID(i,gid_tp,ngid)] +
640 enerd->grpp.ener[egBHAMLR][GID(i,gid_tp,ngid)])*embU;
643 else
645 for(i=0; i<ngid; i++)
647 sum_UgembU[e++] +=
648 (enerd->grpp.ener[egLJSR][GID(i,gid_tp,ngid)] +
649 enerd->grpp.ener[egLJLR][GID(i,gid_tp,ngid)])*embU;
652 if (bDispCorr)
654 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
656 if (bCharge)
658 for(i=0; i<ngid; i++)
660 sum_UgembU[e++] +=
661 (enerd->grpp.ener[egCOULSR][GID(i,gid_tp,ngid)] +
662 enerd->grpp.ener[egCOULLR][GID(i,gid_tp,ngid)])*embU;
664 if (bRFExcl)
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)
677 bin[0]++;
679 else
681 i = (int)((bU_logV_bin_limit
682 - (beta*epot - logV + refvolshift))*invbinw
683 + 0.5);
684 if (i < 0)
686 i = 0;
688 if (i >= nbin)
690 realloc_bins(&bin,&nbin,i+10);
692 bin[i]++;
695 if (debug)
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);
711 if (PAR(cr))
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);
718 frame++;
719 V_all += V;
720 VembU_all += V*sum_embU/nsteps;
722 if (fp_tpi)
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,
734 sum_embU/nsteps,V);
735 for(e=0; e<nener; e++)
737 fprintf(fp_tpi," %12.5e",sum_UgembU[e]/nsteps);
739 fprintf(fp_tpi,"\n");
740 fflush(fp_tpi);
743 bNotLastFrame = read_next_frame(oenv, status,&rerun_fr);
744 } /* End of the loop */
745 runtime_end(runtime);
747 close_trj(status);
749 if (fp_tpi != NULL)
751 gmx_fio_fclose(fp_tpi);
754 if (fplog != NULL)
756 fprintf(fplog,"\n");
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 */
762 if (PAR(cr))
764 /* When running in parallel sum the bins over the processes */
765 i = nbin;
766 global_max(cr,&i);
767 realloc_bins(&bin,&nbin,i);
768 gmx_sumd(nbin,bin,cr);
770 if (MASTER(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",
782 bUlogV,
783 (int)(bin[i]+0.5),
784 bin[i]*exp(-bUlogV)*V_all/VembU_all);
786 gmx_fio_fclose(fp_tpi);
788 sfree(bin);
790 sfree(sum_UgembU);
792 runtime->nsteps_done = frame*inputrec->nsteps;
794 return 0;