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
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROwing Monsters And Cloning Shrimps
40 #include <catamount/dclock.h>
46 #ifdef HAVE_SYS_TIME_H
59 #include "chargegroup.h"
79 #include "gmx_random.h"
82 #include "gmx_wallcycle.h"
85 #include "gromacs/utility/gmxmpi.h"
87 void GenerateGibbsProbabilities(real *ene, real *p_k, real *pks, int minfep, int maxfep)
95 /* find the maximum value */
96 for (i = minfep; i <= maxfep; i++)
103 /* find the denominator */
104 for (i = minfep; i <= maxfep; i++)
106 *pks += exp(ene[i]-maxene);
109 for (i = minfep; i <= maxfep; i++)
111 p_k[i] = exp(ene[i]-maxene) / *pks;
115 void GenerateWeightedGibbsProbabilities(real *ene, real *p_k, real *pks, int nlim, real *nvals, real delta)
124 for (i = 0; i < nlim; i++)
128 /* add the delta, since we need to make sure it's greater than zero, and
129 we need a non-arbitrary number? */
130 nene[i] = ene[i] + log(nvals[i]+delta);
134 nene[i] = ene[i] + log(nvals[i]);
138 /* find the maximum value */
140 for (i = 0; i < nlim; i++)
142 if (nene[i] > maxene)
148 /* subtract off the maximum, avoiding overflow */
149 for (i = 0; i < nlim; i++)
154 /* find the denominator */
155 for (i = 0; i < nlim; i++)
157 *pks += exp(nene[i]);
161 for (i = 0; i < nlim; i++)
163 p_k[i] = exp(nene[i]) / *pks;
168 real do_logsum(int N, real *a_n)
172 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
177 /* compute maximum argument to exp(.) */
180 for (i = 1; i < N; i++)
182 maxarg = max(maxarg, a_n[i]);
185 /* compute sum of exp(a_n - maxarg) */
187 for (i = 0; i < N; i++)
189 sum = sum + exp(a_n[i] - maxarg);
192 /* compute log sum */
193 logsum = log(sum) + maxarg;
197 int FindMinimum(real *min_metric, int N)
204 min_val = min_metric[0];
206 for (nval = 0; nval < N; nval++)
208 if (min_metric[nval] < min_val)
210 min_val = min_metric[nval];
217 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
225 for (i = 0; i < nhisto; i++)
232 /* no samples! is bad!*/
236 nmean /= (real)nhisto;
239 for (i = 0; i < nhisto; i++)
241 /* make sure that all points are in the ratio < x < 1/ratio range */
242 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
251 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_large_int_t step)
255 gmx_bool bDoneEquilibrating = TRUE;
258 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
260 /* calculate the total number of samples */
261 switch (expand->elmceq)
264 /* We have not equilibrated, and won't, ever. */
267 /* we have equilibrated -- we're done */
270 /* first, check if we are equilibrating by steps, if we're still under */
271 if (step < expand->equil_steps)
273 bDoneEquilibrating = FALSE;
278 for (i = 0; i < nlim; i++)
280 totalsamples += dfhist->n_at_lam[i];
282 if (totalsamples < expand->equil_samples)
284 bDoneEquilibrating = FALSE;
288 for (i = 0; i < nlim; i++)
290 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
293 bDoneEquilibrating = FALSE;
299 if (EWL(expand->elamstats)) /* This check is in readir as well, but
302 if (dfhist->wl_delta > expand->equil_wl_delta)
304 bDoneEquilibrating = FALSE;
309 /* we can use the flatness as a judge of good weights, as long as
310 we're not doing minvar, or Wang-Landau.
311 But turn off for now until we figure out exactly how we do this.
314 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
316 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
317 floats for this histogram function. */
320 snew(modhisto, nlim);
321 for (i = 0; i < nlim; i++)
323 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
325 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
329 bDoneEquilibrating = FALSE;
333 bDoneEquilibrating = TRUE;
335 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
337 if (expand->lmc_forced_nstart > 0)
339 for (i = 0; i < nlim; i++)
341 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
344 bDoneEquilibrating = FALSE;
349 return bDoneEquilibrating;
352 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
353 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_large_int_t step)
355 real maxdiff = 0.000000001;
356 gmx_bool bSufficientSamples;
357 int i, k, n, nz, indexi, indexk, min_n, max_n, nlam, totali;
358 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
359 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
360 real de, de_function, dr, denom, maxdr, pks = 0;
361 real min_val, cnval, zero_sum_weights;
362 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
363 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
364 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg, *p_k;
365 real *numweighted_lamee, *logfrac;
367 real chi_m1_0, chi_p1_0, chi_m2_0, chi_p2_0, chi_p1_m1, chi_p2_m1, chi_m1_p1, chi_m2_p1;
369 /* if we have equilibrated the weights, exit now */
375 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
377 dfhist->bEquil = TRUE;
378 /* zero out the visited states so we know how many equilibrated states we have
380 for (i = 0; i < nlim; i++)
382 dfhist->n_at_lam[i] = 0;
387 /* If we reached this far, we have not equilibrated yet, keep on
388 going resetting the weights */
390 if (EWL(expand->elamstats))
392 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
394 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
395 dfhist->wl_histo[fep_state] += 1.0;
397 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
401 /* first increment count */
402 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
403 for (i = 0; i < nlim; i++)
405 dfhist->wl_histo[i] += p_k[i];
408 /* then increment weights (uses count) */
410 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
412 for (i = 0; i < nlim; i++)
414 dfhist->sum_weights[i] -= dfhist->wl_delta*p_k[i];
416 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
421 di = 1+dfhist->wl_delta*p_k[i];
422 dfhist->sum_weights[i] -= log(di);
428 zero_sum_weights = dfhist->sum_weights[0];
429 for (i = 0; i < nlim; i++)
431 dfhist->sum_weights[i] -= zero_sum_weights;
435 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
438 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
439 maxc = 2*expand->c_range+1;
442 snew(lam_variance, nlim);
444 snew(omegap_array, maxc);
445 snew(weightsp_array, maxc);
446 snew(varp_array, maxc);
447 snew(dwp_array, maxc);
449 snew(omegam_array, maxc);
450 snew(weightsm_array, maxc);
451 snew(varm_array, maxc);
452 snew(dwm_array, maxc);
454 /* unpack the current lambdas -- we will only update 2 of these */
456 for (i = 0; i < nlim-1; i++)
457 { /* only through the second to last */
458 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
459 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
462 /* accumulate running averages */
463 for (nval = 0; nval < maxc; nval++)
465 /* constants for later use */
466 cnval = (real)(nval-expand->c_range);
467 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
470 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
471 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
473 de_function = 1.0/(1.0+de);
475 else if (expand->elamstats == elamstatsMETROPOLIS)
483 de_function = 1.0/de;
486 dfhist->accum_m[fep_state][nval] += de_function;
487 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
490 if (fep_state < nlim-1)
492 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
493 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
495 de_function = 1.0/(1.0+de);
497 else if (expand->elamstats == elamstatsMETROPOLIS)
505 de_function = 1.0/de;
508 dfhist->accum_p[fep_state][nval] += de_function;
509 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
512 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
514 n0 = dfhist->n_at_lam[fep_state];
517 nm1 = dfhist->n_at_lam[fep_state-1];
523 if (fep_state < nlim-1)
525 np1 = dfhist->n_at_lam[fep_state+1];
532 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
533 chi_m1_0 = chi_p1_0 = chi_m2_0 = chi_p2_0 = chi_p1_m1 = chi_p2_m1 = chi_m1_p1 = chi_m2_p1 = 0;
537 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
538 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
539 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
540 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
543 if ((fep_state > 0 ) && (nm1 > 0))
545 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
546 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
549 if ((fep_state < nlim-1) && (np1 > 0))
551 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
552 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
566 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
570 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
572 if ((n0 > 0) && (nm1 > 0))
574 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
575 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
579 if (fep_state < nlim-1)
583 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
587 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
589 if ((n0 > 0) && (np1 > 0))
591 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
592 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
598 omegam_array[nval] = omega_m1_0;
602 omegam_array[nval] = 0;
604 weightsm_array[nval] = clam_weightsm;
605 varm_array[nval] = clam_varm;
608 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
612 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
617 omegap_array[nval] = omega_p1_0;
621 omegap_array[nval] = 0;
623 weightsp_array[nval] = clam_weightsp;
624 varp_array[nval] = clam_varp;
625 if ((np1 > 0) && (n0 > 0))
627 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
631 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
636 /* find the C's closest to the old weights value */
638 min_nvalm = FindMinimum(dwm_array, maxc);
639 omega_m1_0 = omegam_array[min_nvalm];
640 clam_weightsm = weightsm_array[min_nvalm];
641 clam_varm = varm_array[min_nvalm];
643 min_nvalp = FindMinimum(dwp_array, maxc);
644 omega_p1_0 = omegap_array[min_nvalp];
645 clam_weightsp = weightsp_array[min_nvalp];
646 clam_varp = varp_array[min_nvalp];
648 clam_osum = omega_m1_0 + omega_p1_0;
652 clam_minvar = 0.5*log(clam_osum);
657 lam_dg[fep_state-1] = clam_weightsm;
658 lam_variance[fep_state-1] = clam_varm;
661 if (fep_state < nlim-1)
663 lam_dg[fep_state] = clam_weightsp;
664 lam_variance[fep_state] = clam_varp;
667 if (expand->elamstats == elamstatsMINVAR)
669 bSufficientSamples = TRUE;
670 /* make sure they are all past a threshold */
671 for (i = 0; i < nlim; i++)
673 if (dfhist->n_at_lam[i] < expand->minvarmin)
675 bSufficientSamples = FALSE;
678 if (bSufficientSamples)
680 dfhist->sum_minvar[fep_state] = clam_minvar;
683 for (i = 0; i < nlim; i++)
685 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
687 expand->minvar_const = clam_minvar;
688 dfhist->sum_minvar[fep_state] = 0.0;
692 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
697 /* we need to rezero minvar now, since it could change at fep_state = 0 */
698 dfhist->sum_dg[0] = 0.0;
699 dfhist->sum_variance[0] = 0.0;
700 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
702 for (i = 1; i < nlim; i++)
704 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
705 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
706 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
713 sfree(weightsm_array);
718 sfree(weightsp_array);
725 static int ChooseNewLambda(FILE *log, int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, real *p_k, gmx_rng_t rng)
727 /* Choose new lambda value, and update transition matrix */
729 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
730 real r1, r2, pks, de_old, de_new, de, trialprob, tprob = 0;
732 real *propose, *accept, *remainder;
734 gmx_bool bRestricted;
736 starting_fep_state = fep_state;
737 lamnew = fep_state; /* so that there is a default setting -- stays the same */
739 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
741 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
743 /* Use a marching method to run through the lambdas and get preliminary free energy data,
744 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
746 /* if we have enough at this lambda, move on to the next one */
748 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
750 lamnew = fep_state+1;
751 if (lamnew == nlim) /* whoops, stepped too far! */
766 snew(remainder, nlim);
768 for (i = 0; i < expand->lmc_repeats; i++)
771 for (ifep = 0; ifep < nlim; ifep++)
777 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
780 /* use the Gibbs sampler, with restricted range */
781 if (expand->gibbsdeltalam < 0)
789 minfep = fep_state - expand->gibbsdeltalam;
790 maxfep = fep_state + expand->gibbsdeltalam;
801 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
803 if (expand->elmcmove == elmcmoveGIBBS)
805 for (ifep = minfep; ifep <= maxfep; ifep++)
807 propose[ifep] = p_k[ifep];
811 r1 = gmx_rng_uniform_real(rng);
812 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
814 if (r1 <= p_k[lamnew])
821 else if (expand->elmcmove == elmcmoveMETGIBBS)
824 /* Metropolized Gibbs sampling */
825 for (ifep = minfep; ifep <= maxfep; ifep++)
827 remainder[ifep] = 1 - p_k[ifep];
830 /* find the proposal probabilities */
832 if (remainder[fep_state] == 0)
834 /* only the current state has any probability */
835 /* we have to stay at the current state */
840 for (ifep = minfep; ifep <= maxfep; ifep++)
842 if (ifep != fep_state)
844 propose[ifep] = p_k[ifep]/remainder[fep_state];
852 r1 = gmx_rng_uniform_real(rng);
853 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
855 pnorm = p_k[lamtrial]/remainder[fep_state];
856 if (lamtrial != fep_state)
866 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
868 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
869 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
870 if (trialprob < tprob)
874 r2 = gmx_rng_uniform_real(rng);
885 /* now figure out the acceptance probability for each */
886 for (ifep = minfep; ifep <= maxfep; ifep++)
889 if (remainder[ifep] != 0)
891 trialprob = (remainder[fep_state])/(remainder[ifep]);
895 trialprob = 1.0; /* this state is the only choice! */
897 if (trialprob < tprob)
901 /* probability for fep_state=0, but that's fine, it's never proposed! */
902 accept[ifep] = tprob;
908 /* it's possible some rounding is failing */
909 if (remainder[fep_state] < 2.0e-15)
911 /* probably numerical rounding error -- no state other than the original has weight */
916 /* probably not a numerical issue */
918 int nerror = 200+(maxfep-minfep+1)*60;
920 snew(errorstr, nerror);
921 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
922 of sum weights. Generated detailed info for failure */
923 loc += sprintf(errorstr, "Something wrong in choosing new lambda state with a Gibbs move -- probably underflow in weight determination.\nDenominator is: %3d%17.10e\n i dE numerator weights\n", 0, pks);
924 for (ifep = minfep; ifep <= maxfep; ifep++)
926 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
928 gmx_fatal(FARGS, errorstr);
932 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
934 /* use the metropolis sampler with trial +/- 1 */
935 r1 = gmx_rng_uniform_real(rng);
940 lamtrial = fep_state;
944 lamtrial = fep_state-1;
949 if (fep_state == nlim-1)
951 lamtrial = fep_state;
955 lamtrial = fep_state+1;
959 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
960 if (expand->elmcmove == elmcmoveMETROPOLIS)
964 if (trialprob < tprob)
968 propose[fep_state] = 0;
969 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
970 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
971 accept[lamtrial] = tprob;
974 else if (expand->elmcmove == elmcmoveBARKER)
976 tprob = 1.0/(1.0+exp(-de));
978 propose[fep_state] = (1-tprob);
979 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
980 accept[fep_state] = 1.0;
981 accept[lamtrial] = 1.0;
984 r2 = gmx_rng_uniform_real(rng);
995 for (ifep = 0; ifep < nlim; ifep++)
997 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
998 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1003 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1012 /* print out the weights to the log, along with current state */
1013 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1014 int nlam, int frequency, gmx_large_int_t step)
1016 int nlim, i, ifep, jfep;
1017 real dw, dg, dv, dm, Tprint;
1019 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1020 gmx_bool bSimTemp = FALSE;
1022 nlim = fep->n_lambda;
1023 if (simtemp != NULL)
1028 if (mod(step, frequency) == 0)
1030 fprintf(outfile, " MC-lambda information\n");
1031 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1033 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1035 fprintf(outfile, " N");
1036 for (i = 0; i < efptNR; i++)
1038 if (fep->separate_dvdl[i])
1040 fprintf(outfile, "%7s", print_names[i]);
1042 else if ((i == efptTEMPERATURE) && bSimTemp)
1044 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1047 fprintf(outfile, " Count ");
1048 if (expand->elamstats == elamstatsMINVAR)
1050 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1054 fprintf(outfile, "G(in kT) dG(in kT)\n");
1056 for (ifep = 0; ifep < nlim; ifep++)
1067 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1068 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1069 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1070 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1073 fprintf(outfile, "%3d", (ifep+1));
1074 for (i = 0; i < efptNR; i++)
1076 if (fep->separate_dvdl[i])
1078 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1080 else if (i == efptTEMPERATURE && bSimTemp)
1082 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1085 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1087 if (expand->elamstats == elamstatsWL)
1089 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1093 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1096 else /* we have equilibrated weights */
1098 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1100 if (expand->elamstats == elamstatsMINVAR)
1102 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1106 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1110 fprintf(outfile, " <<\n");
1114 fprintf(outfile, " \n");
1117 fprintf(outfile, "\n");
1119 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1121 fprintf(outfile, " Transition Matrix\n");
1122 for (ifep = 0; ifep < nlim; ifep++)
1124 fprintf(outfile, "%12d", (ifep+1));
1126 fprintf(outfile, "\n");
1127 for (ifep = 0; ifep < nlim; ifep++)
1129 for (jfep = 0; jfep < nlim; jfep++)
1131 if (dfhist->n_at_lam[ifep] > 0)
1133 if (expand->bSymmetrizedTMatrix)
1135 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1139 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1146 fprintf(outfile, "%12.8f", Tprint);
1148 fprintf(outfile, "%3d\n", (ifep+1));
1151 fprintf(outfile, " Empirical Transition Matrix\n");
1152 for (ifep = 0; ifep < nlim; ifep++)
1154 fprintf(outfile, "%12d", (ifep+1));
1156 fprintf(outfile, "\n");
1157 for (ifep = 0; ifep < nlim; ifep++)
1159 for (jfep = 0; jfep < nlim; jfep++)
1161 if (dfhist->n_at_lam[ifep] > 0)
1163 if (expand->bSymmetrizedTMatrix)
1165 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1169 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1176 fprintf(outfile, "%12.8f", Tprint);
1178 fprintf(outfile, "%3d\n", (ifep+1));
1184 extern void get_mc_state(gmx_rng_t rng, t_state *state)
1186 gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
1189 extern void set_mc_state(gmx_rng_t rng, t_state *state)
1191 gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
1194 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1195 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
1196 gmx_large_int_t step, gmx_rng_t mcrng,
1197 rvec *v, t_mdatoms *mdatoms)
1199 real *pfep_lamee, *p_k, *scaled_lamee, *weighted_lamee;
1200 int i, nlam, nlim, lamnew, totalsamples;
1201 real oneovert, maxscaled = 0, maxweighted = 0;
1204 double *temperature_lambdas;
1205 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1207 expand = ir->expandedvals;
1208 simtemp = ir->simtempvals;
1209 nlim = ir->fepvals->n_lambda;
1210 nlam = state->fep_state;
1212 snew(scaled_lamee, nlim);
1213 snew(weighted_lamee, nlim);
1214 snew(pfep_lamee, nlim);
1217 if (expand->bInit_weights) /* if initialized weights, we need to fill them in */
1219 dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
1220 for (i = 0; i < nlim; i++)
1222 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
1223 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
1225 expand->bInit_weights = FALSE;
1228 /* update the count at the current lambda*/
1229 dfhist->n_at_lam[nlam]++;
1231 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1232 pressure controlled.*/
1235 where does this PV term go?
1236 for (i=0;i<nlim;i++)
1238 fep_lamee[i] += pVTerm;
1242 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1243 /* we don't need to include the pressure term, since the volume is the same between the two.
1244 is there some term we are neglecting, however? */
1246 if (ir->efep != efepNO)
1248 for (i = 0; i < nlim; i++)
1252 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1253 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1254 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
1258 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1259 /* mc_temp is currently set to the system reft unless otherwise defined */
1262 /* save these energies for printing, so they don't get overwritten by the next step */
1263 /* they aren't overwritten in the non-free energy case, but we always print with these
1271 for (i = 0; i < nlim; i++)
1273 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
1278 for (i = 0; i < nlim; i++)
1280 pfep_lamee[i] = scaled_lamee[i];
1282 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1285 maxscaled = scaled_lamee[i];
1286 maxweighted = weighted_lamee[i];
1290 if (scaled_lamee[i] > maxscaled)
1292 maxscaled = scaled_lamee[i];
1294 if (weighted_lamee[i] > maxweighted)
1296 maxweighted = weighted_lamee[i];
1301 for (i = 0; i < nlim; i++)
1303 scaled_lamee[i] -= maxscaled;
1304 weighted_lamee[i] -= maxweighted;
1307 /* update weights - we decide whether or not to actually do this inside */
1309 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, nlam, scaled_lamee, weighted_lamee, step);
1310 if (bDoneEquilibrating)
1314 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1318 lamnew = ChooseNewLambda(log, nlim, expand, dfhist, nlam, weighted_lamee, p_k, mcrng);
1319 /* if using simulated tempering, we need to adjust the temperatures */
1320 if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
1325 int nstart, nend, gt;
1327 snew(buf_ngtc, ir->opts.ngtc);
1329 for (i = 0; i < ir->opts.ngtc; i++)
1331 if (ir->opts.ref_t[i] > 0)
1333 told = ir->opts.ref_t[i];
1334 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1335 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1339 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1341 nstart = mdatoms->start;
1342 nend = nstart + mdatoms->homenr;
1343 for (n = nstart; n < nend; n++)
1348 gt = mdatoms->cTC[n];
1350 for (d = 0; d < DIM; d++)
1352 v[n][d] *= buf_ngtc[gt];
1356 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1358 /* we need to recalculate the masses if the temperature has changed */
1359 init_npt_masses(ir, state, MassQ, FALSE);
1360 for (i = 0; i < state->nnhpres; i++)
1362 for (j = 0; j < ir->opts.nhchainlength; j++)
1364 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1367 for (i = 0; i < ir->opts.ngtc; i++)
1369 for (j = 0; j < ir->opts.nhchainlength; j++)
1371 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1378 /* now check on the Wang-Landau updating critera */
1380 if (EWL(expand->elamstats))
1382 bSwitchtoOneOverT = FALSE;
1383 if (expand->bWLoneovert)
1386 for (i = 0; i < nlim; i++)
1388 totalsamples += dfhist->n_at_lam[i];
1390 oneovert = (1.0*nlim)/totalsamples;
1391 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1392 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1393 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1394 (dfhist->wl_delta < expand->init_wl_delta))
1396 bSwitchtoOneOverT = TRUE;
1399 if (bSwitchtoOneOverT)
1401 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1405 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1408 for (i = 0; i < nlim; i++)
1410 dfhist->wl_histo[i] = 0;
1412 dfhist->wl_delta *= expand->wl_scale;
1415 fprintf(log, "\nStep %d: weights are now:", (int)step);
1416 for (i = 0; i < nlim; i++)
1418 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1426 sfree(scaled_lamee);
1427 sfree(weighted_lamee);