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 * GROwing Monsters And Cloning Shrimps
51 #include "chargegroup.h"
56 #include "gmx_fatal.h"
71 #include "gmx_random.h"
76 #include "gmx_wallcycle.h"
77 #include "mtop_util.h"
83 #include "gmx_x86_sse2.h"
87 static void global_max(t_commrec *cr, int *n)
91 snew(sum, cr->nnodes);
93 gmx_sumi(cr->nnodes, sum, cr);
94 for (i = 0; i < cr->nnodes; i++)
102 static void realloc_bins(double **bin, int *nbin, int nbin_new)
106 if (nbin_new != *nbin)
108 srenew(*bin, nbin_new);
109 for (i = *nbin; i < nbin_new; i++)
117 double do_tpi(FILE *fplog, t_commrec *cr,
118 int nfile, const t_filenm fnm[],
119 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
121 gmx_vsite_t *vsite, gmx_constr_t constr,
123 t_inputrec *inputrec,
124 gmx_mtop_t *top_global, t_fcdata *fcd,
127 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
130 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
132 real cpt_period, real max_hours,
133 const char *deviceOptions,
135 gmx_runtime_t *runtime)
137 const char *TPI = "Test Particle Insertion";
139 gmx_groups_t *groups;
140 gmx_enerdata_t *enerd;
142 real lambda, t, temp, beta, drmax, epot;
143 double embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
146 gmx_bool bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS, bOurStep;
147 tensor force_vir, shake_vir, vir, pres;
148 int cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
150 rvec mu_tot, x_init, dx, x_tp;
151 int nnodes, frame, nsteps, step;
155 char *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
156 double dbl, dump_ener;
158 int nat_cavity = 0, d;
159 real *mass_cavity = NULL, mass_tot;
161 double invbinw, *bin, refvolshift, logV, bUlogV;
162 real dvdl, prescorr, enercorr, dvdlcorr;
163 gmx_bool bEnergyOutOfBounds;
164 const char *tpid_leg[2] = {"direct", "reweighted"};
166 /* Since there is no upper limit to the insertion energies,
167 * we need to set an upper limit for the distribution output.
169 real bU_bin_limit = 50;
170 real bU_logV_bin_limit = bU_bin_limit + 10;
174 top = gmx_mtop_generate_local_top(top_global, inputrec);
176 groups = &top_global->groups;
178 bCavity = (inputrec->eI == eiTPIC);
181 ptr = getenv("GMX_TPIC_MASSES");
188 /* Read (multiple) masses from env var GMX_TPIC_MASSES,
189 * The center of mass of the last atoms is then used for TPIC.
192 while (sscanf(ptr, "%lf%n", &dbl, &i) > 0)
194 srenew(mass_cavity, nat_cavity+1);
195 mass_cavity[nat_cavity] = dbl;
196 fprintf(fplog, "mass[%d] = %f\n",
197 nat_cavity+1, mass_cavity[nat_cavity]);
203 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
209 init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
210 state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
211 /* We never need full pbc for TPI */
213 /* Determine the temperature for the Boltzmann weighting */
214 temp = inputrec->opts.ref_t[0];
217 for (i = 1; (i < inputrec->opts.ngtc); i++)
219 if (inputrec->opts.ref_t[i] != temp)
221 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
222 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
226 "\n The temperature for test particle insertion is %.3f K\n\n",
229 beta = 1.0/(BOLTZ*temp);
231 /* Number of insertions per frame */
232 nsteps = inputrec->nsteps;
234 /* Use the same neighborlist with more insertions points
235 * in a sphere of radius drmax around the initial point
237 /* This should be a proper mdp parameter */
238 drmax = inputrec->rtpi;
240 /* An environment variable can be set to dump all configurations
241 * to pdb with an insertion energy <= this value.
243 dump_pdb = getenv("GMX_TPI_DUMP");
247 sscanf(dump_pdb, "%lf", &dump_ener);
250 atoms2md(top_global, inputrec, 0, NULL, 0, top_global->natoms, mdatoms);
251 update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
254 init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
255 snew(f, top_global->natoms);
257 /* Print to log file */
258 runtime_start(runtime);
259 print_date_and_time(fplog, cr->nodeid,
260 "Started Test Particle Insertion", runtime);
261 wallcycle_start(wcycle, ewcRUN);
263 /* The last charge group is the group to be inserted */
264 cg_tp = top->cgs.nr - 1;
265 a_tp0 = top->cgs.index[cg_tp];
266 a_tp1 = top->cgs.index[cg_tp+1];
269 fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
271 if (a_tp1 - a_tp0 > 1 &&
272 (inputrec->rlist < inputrec->rcoulomb ||
273 inputrec->rlist < inputrec->rvdw))
275 gmx_fatal(FARGS, "Can not do TPI for multi-atom molecule with a twin-range cut-off");
277 snew(x_mol, a_tp1-a_tp0);
279 bDispCorr = (inputrec->eDispCorr != edispcNO);
281 for (i = a_tp0; i < a_tp1; i++)
283 /* Copy the coordinates of the molecule to be insterted */
284 copy_rvec(state->x[i], x_mol[i-a_tp0]);
285 /* Check if we need to print electrostatic energies */
286 bCharge |= (mdatoms->chargeA[i] != 0 ||
287 (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
289 bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype != eelRF_NEC);
291 calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state->x, fr->cg_cm);
294 if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
296 fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
297 fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
302 /* Center the molecule to be inserted at zero */
303 for (i = 0; i < a_tp1-a_tp0; i++)
305 rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
311 fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
312 a_tp1-a_tp0, bCharge ? "with" : "without");
314 fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
315 nsteps, opt2fn("-rerun", nfile, fnm));
320 if (inputrec->nstlist > 1)
322 if (drmax == 0 && a_tp1-a_tp0 == 1)
324 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);
328 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
336 fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
340 ngid = groups->grps[egcENER].nr;
341 gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
354 if (EEL_FULL(fr->eeltype))
359 snew(sum_UgembU, nener);
361 /* Initialize random generator */
362 tpi_rand = gmx_rng_init(inputrec->ld_seed);
366 fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
367 "TPI energies", "Time (ps)",
368 "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
369 xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
372 sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
373 leg[e++] = strdup(str);
374 sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
375 leg[e++] = strdup(str);
376 sprintf(str, "f. <e\\S-\\betaU\\N>");
377 leg[e++] = strdup(str);
378 sprintf(str, "f. V");
379 leg[e++] = strdup(str);
380 sprintf(str, "f. <Ue\\S-\\betaU\\N>");
381 leg[e++] = strdup(str);
382 for (i = 0; i < ngid; i++)
384 sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
385 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
386 leg[e++] = strdup(str);
390 sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
391 leg[e++] = strdup(str);
395 for (i = 0; i < ngid; i++)
397 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
398 *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
399 leg[e++] = strdup(str);
403 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
404 leg[e++] = strdup(str);
406 if (EEL_FULL(fr->eeltype))
408 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
409 leg[e++] = strdup(str);
412 xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
413 for (i = 0; i < 4+nener; i++)
427 bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
428 &rerun_fr, TRX_NEED_X);
431 if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
432 mdatoms->nr - (a_tp1 - a_tp0))
434 gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
435 "is not equal the number in the run input file (%d) "
436 "minus the number of atoms to insert (%d)\n",
437 rerun_fr.natoms, bCavity ? " minus one" : "",
438 mdatoms->nr, a_tp1-a_tp0);
441 refvolshift = log(det(rerun_fr.box));
444 /* Make sure we don't detect SSE overflow generated before this point */
445 gmx_mm_check_and_reset_overflow();
448 while (bNotLastFrame)
450 lambda = rerun_fr.lambda;
454 for (e = 0; e < nener; e++)
459 /* Copy the coordinates from the input trajectory */
460 for (i = 0; i < rerun_fr.natoms; i++)
462 copy_rvec(rerun_fr.x[i], state->x[i]);
464 copy_mat(rerun_fr.box, state->box);
469 bStateChanged = TRUE;
471 for (step = 0; step < nsteps; step++)
473 /* In parallel all nodes generate all random configurations.
474 * In that way the result is identical to a single cpu tpi run.
478 /* Random insertion in the whole volume */
479 bNS = (step % inputrec->nstlist == 0);
482 /* Generate a random position in the box */
483 x_init[XX] = gmx_rng_uniform_real(tpi_rand)*state->box[XX][XX];
484 x_init[YY] = gmx_rng_uniform_real(tpi_rand)*state->box[YY][YY];
485 x_init[ZZ] = gmx_rng_uniform_real(tpi_rand)*state->box[ZZ][ZZ];
487 if (inputrec->nstlist == 1)
489 copy_rvec(x_init, x_tp);
493 /* Generate coordinates within |dx|=drmax of x_init */
496 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
497 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
498 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
500 while (norm2(dx) > drmax*drmax);
501 rvec_add(x_init, dx, x_tp);
506 /* Random insertion around a cavity location
507 * given by the last coordinate of the trajectory.
513 /* Copy the location of the cavity */
514 copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
518 /* Determine the center of mass of the last molecule */
521 for (i = 0; i < nat_cavity; i++)
523 for (d = 0; d < DIM; d++)
526 mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
528 mass_tot += mass_cavity[i];
530 for (d = 0; d < DIM; d++)
532 x_init[d] /= mass_tot;
536 /* Generate coordinates within |dx|=drmax of x_init */
539 dx[XX] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
540 dx[YY] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
541 dx[ZZ] = (2*gmx_rng_uniform_real(tpi_rand) - 1)*drmax;
543 while (norm2(dx) > drmax*drmax);
544 rvec_add(x_init, dx, x_tp);
547 if (a_tp1 - a_tp0 == 1)
549 /* Insert a single atom, just copy the insertion location */
550 copy_rvec(x_tp, state->x[a_tp0]);
554 /* Copy the coordinates from the top file */
555 for (i = a_tp0; i < a_tp1; i++)
557 copy_rvec(x_mol[i-a_tp0], state->x[i]);
559 /* Rotate the molecule randomly */
560 rotate_conf(a_tp1-a_tp0, state->x+a_tp0, NULL,
561 2*M_PI*gmx_rng_uniform_real(tpi_rand),
562 2*M_PI*gmx_rng_uniform_real(tpi_rand),
563 2*M_PI*gmx_rng_uniform_real(tpi_rand));
564 /* Shift to the insertion location */
565 for (i = a_tp0; i < a_tp1; i++)
567 rvec_inc(state->x[i], x_tp);
571 /* Check if this insertion belongs to this node */
575 switch (inputrec->eI)
578 bOurStep = ((step / inputrec->nstlist) % nnodes == cr->nodeid);
581 bOurStep = (step % nnodes == cr->nodeid);
584 gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
589 /* Clear some matrix variables */
590 clear_mat(force_vir);
591 clear_mat(shake_vir);
595 /* Set the charge group center of mass of the test particle */
596 copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
598 /* Calc energy (no forces) on new positions.
599 * Since we only need the intermolecular energy
600 * and the RF exclusion terms of the inserted molecule occur
601 * within a single charge group we can pass NULL for the graph.
602 * This also avoids shifts that would move charge groups
605 * Some checks above ensure than we can not have
606 * twin-range interactions together with nstlist > 1,
607 * therefore we do not need to remember the LR energies.
609 /* Make do_force do a single node force calculation */
611 do_force(fplog, cr, inputrec,
612 step, nrnb, wcycle, top, top_global, &top_global->groups,
613 state->box, state->x, &state->hist,
614 f, force_vir, mdatoms, enerd, fcd,
616 NULL, fr, NULL, mu_tot, t, NULL, NULL, FALSE,
617 GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
618 (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS | GMX_FORCE_DO_LR : 0) |
619 (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
621 bStateChanged = FALSE;
624 /* Calculate long range corrections to pressure and energy */
625 calc_dispcorr(fplog, inputrec, fr, step, top_global->natoms, state->box,
626 lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
627 /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
628 enerd->term[F_DISPCORR] = enercorr;
629 enerd->term[F_EPOT] += enercorr;
630 enerd->term[F_PRES] += prescorr;
631 enerd->term[F_DVDL_VDW] += dvdlcorr;
633 epot = enerd->term[F_EPOT];
634 bEnergyOutOfBounds = FALSE;
636 /* With SSE the energy can overflow, check for this */
637 if (gmx_mm_check_and_reset_overflow())
641 fprintf(debug, "Found an SSE overflow, assuming the energy is out of bounds\n");
643 bEnergyOutOfBounds = TRUE;
646 /* If the compiler doesn't optimize this check away
647 * we catch the NAN energies.
648 * The epot>GMX_REAL_MAX check catches inf values,
649 * which should nicely result in embU=0 through the exp below,
650 * but it does not hurt to check anyhow.
652 /* Non-bonded Interaction usually diverge at r=0.
653 * With tabulated interaction functions the first few entries
654 * should be capped in a consistent fashion between
655 * repulsion, dispersion and Coulomb to avoid accidental
656 * negative values in the total energy.
657 * The table generation code in tables.c does this.
658 * With user tbales the user should take care of this.
660 if (epot != epot || epot > GMX_REAL_MAX)
662 bEnergyOutOfBounds = TRUE;
664 if (bEnergyOutOfBounds)
668 fprintf(debug, "\n time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, step, epot);
674 embU = exp(-beta*epot);
676 /* Determine the weighted energy contributions of each energy group */
678 sum_UgembU[e++] += epot*embU;
681 for (i = 0; i < ngid; i++)
684 (enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] +
685 enerd->grpp.ener[egBHAMLR][GID(i, gid_tp, ngid)])*embU;
690 for (i = 0; i < ngid; i++)
693 (enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] +
694 enerd->grpp.ener[egLJLR][GID(i, gid_tp, ngid)])*embU;
699 sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
703 for (i = 0; i < ngid; i++)
706 (enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] +
707 enerd->grpp.ener[egCOULLR][GID(i, gid_tp, ngid)])*embU;
711 sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
713 if (EEL_FULL(fr->eeltype))
715 sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
720 if (embU == 0 || beta*epot > bU_bin_limit)
726 i = (int)((bU_logV_bin_limit
727 - (beta*epot - logV + refvolshift))*invbinw
735 realloc_bins(&bin, &nbin, i+10);
742 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
743 step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
746 if (dump_pdb && epot <= dump_ener)
748 sprintf(str, "t%g_step%d.pdb", t, step);
749 sprintf(str2, "t: %f step %d ener: %f", t, step, epot);
750 write_sto_conf_mtop(str, str2, top_global, state->x, state->v,
751 inputrec->ePBC, state->box);
758 /* When running in parallel sum the energies over the processes */
759 gmx_sumd(1, &sum_embU, cr);
760 gmx_sumd(nener, sum_UgembU, cr);
765 VembU_all += V*sum_embU/nsteps;
769 if (bVerbose || frame%10 == 0 || frame < 10)
771 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
772 -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
775 fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
777 VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
778 sum_embU == 0 ? 20/beta : -log(sum_embU/nsteps)/beta,
780 for (e = 0; e < nener; e++)
782 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
784 fprintf(fp_tpi, "\n");
788 bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
789 } /* End of the loop */
790 runtime_end(runtime);
796 gmx_fio_fclose(fp_tpi);
801 fprintf(fplog, "\n");
802 fprintf(fplog, " <V> = %12.5e nm^3\n", V_all/frame);
803 fprintf(fplog, " <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
806 /* Write the Boltzmann factor histogram */
809 /* When running in parallel sum the bins over the processes */
812 realloc_bins(&bin, &nbin, i);
813 gmx_sumd(nbin, bin, cr);
817 fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
818 "TPI energy distribution",
819 "\\betaU - log(V/<V>)", "count", oenv);
820 sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
821 xvgr_subtitle(fp_tpi, str, oenv);
822 xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
823 for (i = nbin-1; i > 0; i--)
825 bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
826 fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
829 bin[i]*exp(-bUlogV)*V_all/VembU_all);
831 gmx_fio_fclose(fp_tpi);
837 runtime->nsteps_done = frame*inputrec->nsteps;