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"
80 #include "gmx_random.h"
83 #include "gmx_wallcycle.h"
93 void GenerateGibbsProbabilities(real *ene, real *p_k, real *pks, int minfep, int maxfep)
100 maxene = ene[minfep];
101 /* find the maximum value */
102 for (i = minfep; i <= maxfep; i++)
109 /* find the denominator */
110 for (i = minfep; i <= maxfep; i++)
112 *pks += exp(ene[i]-maxene);
115 for (i = minfep; i <= maxfep; i++)
117 p_k[i] = exp(ene[i]-maxene) / *pks;
121 void GenerateWeightedGibbsProbabilities(real *ene, real *p_k, real *pks, int nlim, real *nvals, real delta)
130 for (i = 0; i < nlim; i++)
134 /* add the delta, since we need to make sure it's greater than zero, and
135 we need a non-arbitrary number? */
136 nene[i] = ene[i] + log(nvals[i]+delta);
140 nene[i] = ene[i] + log(nvals[i]);
144 /* find the maximum value */
146 for (i = 0; i < nlim; i++)
148 if (nene[i] > maxene)
154 /* subtract off the maximum, avoiding overflow */
155 for (i = 0; i < nlim; i++)
160 /* find the denominator */
161 for (i = 0; i < nlim; i++)
163 *pks += exp(nene[i]);
167 for (i = 0; i < nlim; i++)
169 p_k[i] = exp(nene[i]) / *pks;
174 real do_logsum(int N, real *a_n)
178 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
183 /* compute maximum argument to exp(.) */
186 for (i = 1; i < N; i++)
188 maxarg = max(maxarg, a_n[i]);
191 /* compute sum of exp(a_n - maxarg) */
193 for (i = 0; i < N; i++)
195 sum = sum + exp(a_n[i] - maxarg);
198 /* compute log sum */
199 logsum = log(sum) + maxarg;
203 int FindMinimum(real *min_metric, int N)
210 min_val = min_metric[0];
212 for (nval = 0; nval < N; nval++)
214 if (min_metric[nval] < min_val)
216 min_val = min_metric[nval];
223 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
231 for (i = 0; i < nhisto; i++)
238 /* no samples! is bad!*/
242 nmean /= (real)nhisto;
245 for (i = 0; i < nhisto; i++)
247 /* make sure that all points are in the ratio < x < 1/ratio range */
248 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
257 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_large_int_t step)
261 gmx_bool bDoneEquilibrating = TRUE;
264 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
266 /* calculate the total number of samples */
267 switch (expand->elmceq)
270 /* We have not equilibrated, and won't, ever. */
273 /* we have equilibrated -- we're done */
276 /* first, check if we are equilibrating by steps, if we're still under */
277 if (step < expand->equil_steps)
279 bDoneEquilibrating = FALSE;
284 for (i = 0; i < nlim; i++)
286 totalsamples += dfhist->n_at_lam[i];
288 if (totalsamples < expand->equil_samples)
290 bDoneEquilibrating = FALSE;
294 for (i = 0; i < nlim; i++)
296 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
299 bDoneEquilibrating = FALSE;
305 if (EWL(expand->elamstats)) /* This check is in readir as well, but
308 if (dfhist->wl_delta > expand->equil_wl_delta)
310 bDoneEquilibrating = FALSE;
315 /* we can use the flatness as a judge of good weights, as long as
316 we're not doing minvar, or Wang-Landau.
317 But turn off for now until we figure out exactly how we do this.
320 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
322 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
323 floats for this histogram function. */
326 snew(modhisto, nlim);
327 for (i = 0; i < nlim; i++)
329 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
331 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
335 bDoneEquilibrating = FALSE;
339 bDoneEquilibrating = TRUE;
341 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
343 if (expand->lmc_forced_nstart > 0)
345 for (i = 0; i < nlim; i++)
347 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
350 bDoneEquilibrating = FALSE;
355 return bDoneEquilibrating;
358 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
359 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_large_int_t step)
361 real maxdiff = 0.000000001;
362 gmx_bool bSufficientSamples;
363 int i, k, n, nz, indexi, indexk, min_n, max_n, nlam, totali;
364 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
365 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
366 real de, de_function, dr, denom, maxdr, pks = 0;
367 real min_val, cnval, zero_sum_weights;
368 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
369 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
370 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg, *p_k;
371 real *numweighted_lamee, *logfrac;
373 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;
375 /* if we have equilibrated the weights, exit now */
381 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
383 dfhist->bEquil = TRUE;
384 /* zero out the visited states so we know how many equilibrated states we have
386 for (i = 0; i < nlim; i++)
388 dfhist->n_at_lam[i] = 0;
393 /* If we reached this far, we have not equilibrated yet, keep on
394 going resetting the weights */
396 if (EWL(expand->elamstats))
398 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
400 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
401 dfhist->wl_histo[fep_state] += 1.0;
403 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
407 /* first increment count */
408 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
409 for (i = 0; i < nlim; i++)
411 dfhist->wl_histo[i] += p_k[i];
414 /* then increment weights (uses count) */
416 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
418 for (i = 0; i < nlim; i++)
420 dfhist->sum_weights[i] -= dfhist->wl_delta*p_k[i];
422 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
427 di = 1+dfhist->wl_delta*p_k[i];
428 dfhist->sum_weights[i] -= log(di);
434 zero_sum_weights = dfhist->sum_weights[0];
435 for (i = 0; i < nlim; i++)
437 dfhist->sum_weights[i] -= zero_sum_weights;
441 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
444 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
445 maxc = 2*expand->c_range+1;
448 snew(lam_variance, nlim);
450 snew(omegap_array, maxc);
451 snew(weightsp_array, maxc);
452 snew(varp_array, maxc);
453 snew(dwp_array, maxc);
455 snew(omegam_array, maxc);
456 snew(weightsm_array, maxc);
457 snew(varm_array, maxc);
458 snew(dwm_array, maxc);
460 /* unpack the current lambdas -- we will only update 2 of these */
462 for (i = 0; i < nlim-1; i++)
463 { /* only through the second to last */
464 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
465 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
468 /* accumulate running averages */
469 for (nval = 0; nval < maxc; nval++)
471 /* constants for later use */
472 cnval = (real)(nval-expand->c_range);
473 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
476 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
477 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
479 de_function = 1.0/(1.0+de);
481 else if (expand->elamstats == elamstatsMETROPOLIS)
489 de_function = 1.0/de;
492 dfhist->accum_m[fep_state][nval] += de_function;
493 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
496 if (fep_state < nlim-1)
498 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
499 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
501 de_function = 1.0/(1.0+de);
503 else if (expand->elamstats == elamstatsMETROPOLIS)
511 de_function = 1.0/de;
514 dfhist->accum_p[fep_state][nval] += de_function;
515 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
518 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
520 n0 = dfhist->n_at_lam[fep_state];
523 nm1 = dfhist->n_at_lam[fep_state-1];
529 if (fep_state < nlim-1)
531 np1 = dfhist->n_at_lam[fep_state+1];
538 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
539 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;
543 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
544 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
545 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
546 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
549 if ((fep_state > 0 ) && (nm1 > 0))
551 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
552 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
555 if ((fep_state < nlim-1) && (np1 > 0))
557 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
558 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
572 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
576 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
578 if ((n0 > 0) && (nm1 > 0))
580 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
581 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
585 if (fep_state < nlim-1)
589 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
593 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
595 if ((n0 > 0) && (np1 > 0))
597 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
598 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
604 omegam_array[nval] = omega_m1_0;
608 omegam_array[nval] = 0;
610 weightsm_array[nval] = clam_weightsm;
611 varm_array[nval] = clam_varm;
614 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
618 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
623 omegap_array[nval] = omega_p1_0;
627 omegap_array[nval] = 0;
629 weightsp_array[nval] = clam_weightsp;
630 varp_array[nval] = clam_varp;
631 if ((np1 > 0) && (n0 > 0))
633 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
637 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
642 /* find the C's closest to the old weights value */
644 min_nvalm = FindMinimum(dwm_array, maxc);
645 omega_m1_0 = omegam_array[min_nvalm];
646 clam_weightsm = weightsm_array[min_nvalm];
647 clam_varm = varm_array[min_nvalm];
649 min_nvalp = FindMinimum(dwp_array, maxc);
650 omega_p1_0 = omegap_array[min_nvalp];
651 clam_weightsp = weightsp_array[min_nvalp];
652 clam_varp = varp_array[min_nvalp];
654 clam_osum = omega_m1_0 + omega_p1_0;
658 clam_minvar = 0.5*log(clam_osum);
663 lam_dg[fep_state-1] = clam_weightsm;
664 lam_variance[fep_state-1] = clam_varm;
667 if (fep_state < nlim-1)
669 lam_dg[fep_state] = clam_weightsp;
670 lam_variance[fep_state] = clam_varp;
673 if (expand->elamstats == elamstatsMINVAR)
675 bSufficientSamples = TRUE;
676 /* make sure they are all past a threshold */
677 for (i = 0; i < nlim; i++)
679 if (dfhist->n_at_lam[i] < expand->minvarmin)
681 bSufficientSamples = FALSE;
684 if (bSufficientSamples)
686 dfhist->sum_minvar[fep_state] = clam_minvar;
689 for (i = 0; i < nlim; i++)
691 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
693 expand->minvar_const = clam_minvar;
694 dfhist->sum_minvar[fep_state] = 0.0;
698 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
703 /* we need to rezero minvar now, since it could change at fep_state = 0 */
704 dfhist->sum_dg[0] = 0.0;
705 dfhist->sum_variance[0] = 0.0;
706 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
708 for (i = 1; i < nlim; i++)
710 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
711 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
712 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
719 sfree(weightsm_array);
724 sfree(weightsp_array);
731 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)
733 /* Choose new lambda value, and update transition matrix */
735 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
736 real r1, r2, pks, de_old, de_new, de, trialprob, tprob = 0;
738 real *propose, *accept, *remainder;
740 gmx_bool bRestricted;
742 starting_fep_state = fep_state;
743 lamnew = fep_state; /* so that there is a default setting -- stays the same */
745 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
747 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
749 /* Use a marching method to run through the lambdas and get preliminary free energy data,
750 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
752 /* if we have enough at this lambda, move on to the next one */
754 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
756 lamnew = fep_state+1;
757 if (lamnew == nlim) /* whoops, stepped too far! */
772 snew(remainder, nlim);
774 for (i = 0; i < expand->lmc_repeats; i++)
777 for (ifep = 0; ifep < nlim; ifep++)
783 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
786 /* use the Gibbs sampler, with restricted range */
787 if (expand->gibbsdeltalam < 0)
795 minfep = fep_state - expand->gibbsdeltalam;
796 maxfep = fep_state + expand->gibbsdeltalam;
807 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
809 if (expand->elmcmove == elmcmoveGIBBS)
811 for (ifep = minfep; ifep <= maxfep; ifep++)
813 propose[ifep] = p_k[ifep];
817 r1 = gmx_rng_uniform_real(rng);
818 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
820 if (r1 <= p_k[lamnew])
827 else if (expand->elmcmove == elmcmoveMETGIBBS)
830 /* Metropolized Gibbs sampling */
831 for (ifep = minfep; ifep <= maxfep; ifep++)
833 remainder[ifep] = 1 - p_k[ifep];
836 /* find the proposal probabilities */
838 if (remainder[fep_state] == 0)
840 /* only the current state has any probability */
841 /* we have to stay at the current state */
846 for (ifep = minfep; ifep <= maxfep; ifep++)
848 if (ifep != fep_state)
850 propose[ifep] = p_k[ifep]/remainder[fep_state];
858 r1 = gmx_rng_uniform_real(rng);
859 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
861 pnorm = p_k[lamtrial]/remainder[fep_state];
862 if (lamtrial != fep_state)
872 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
874 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
875 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
876 if (trialprob < tprob)
880 r2 = gmx_rng_uniform_real(rng);
891 /* now figure out the acceptance probability for each */
892 for (ifep = minfep; ifep <= maxfep; ifep++)
895 if (remainder[ifep] != 0)
897 trialprob = (remainder[fep_state])/(remainder[ifep]);
901 trialprob = 1.0; /* this state is the only choice! */
903 if (trialprob < tprob)
907 /* probability for fep_state=0, but that's fine, it's never proposed! */
908 accept[ifep] = tprob;
914 /* it's possible some rounding is failing */
915 if (remainder[fep_state] < 2.0e-15)
917 /* probably numerical rounding error -- no state other than the original has weight */
922 /* probably not a numerical issue */
924 int nerror = 200+(maxfep-minfep+1)*60;
926 snew(errorstr, nerror);
927 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
928 of sum weights. Generated detailed info for failure */
929 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);
930 for (ifep = minfep; ifep <= maxfep; ifep++)
932 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
934 gmx_fatal(FARGS, errorstr);
938 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
940 /* use the metropolis sampler with trial +/- 1 */
941 r1 = gmx_rng_uniform_real(rng);
946 lamtrial = fep_state;
950 lamtrial = fep_state-1;
955 if (fep_state == nlim-1)
957 lamtrial = fep_state;
961 lamtrial = fep_state+1;
965 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
966 if (expand->elmcmove == elmcmoveMETROPOLIS)
970 if (trialprob < tprob)
974 propose[fep_state] = 0;
975 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
976 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
977 accept[lamtrial] = tprob;
980 else if (expand->elmcmove == elmcmoveBARKER)
982 tprob = 1.0/(1.0+exp(-de));
984 propose[fep_state] = (1-tprob);
985 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
986 accept[fep_state] = 1.0;
987 accept[lamtrial] = 1.0;
990 r2 = gmx_rng_uniform_real(rng);
1001 for (ifep = 0; ifep < nlim; ifep++)
1003 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1004 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1009 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1018 /* print out the weights to the log, along with current state */
1019 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1020 int nlam, int frequency, gmx_large_int_t step)
1022 int nlim, i, ifep, jfep;
1023 real dw, dg, dv, dm, Tprint;
1025 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1026 gmx_bool bSimTemp = FALSE;
1028 nlim = fep->n_lambda;
1029 if (simtemp != NULL)
1034 if (mod(step, frequency) == 0)
1036 fprintf(outfile, " MC-lambda information\n");
1037 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1039 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1041 fprintf(outfile, " N");
1042 for (i = 0; i < efptNR; i++)
1044 if (fep->separate_dvdl[i])
1046 fprintf(outfile, "%7s", print_names[i]);
1048 else if ((i == efptTEMPERATURE) && bSimTemp)
1050 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1053 fprintf(outfile, " Count ");
1054 if (expand->elamstats == elamstatsMINVAR)
1056 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1060 fprintf(outfile, "G(in kT) dG(in kT)\n");
1062 for (ifep = 0; ifep < nlim; ifep++)
1073 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1074 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1075 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1076 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1079 fprintf(outfile, "%3d", (ifep+1));
1080 for (i = 0; i < efptNR; i++)
1082 if (fep->separate_dvdl[i])
1084 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1086 else if (i == efptTEMPERATURE && bSimTemp)
1088 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1091 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1093 if (expand->elamstats == elamstatsWL)
1095 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1099 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1102 else /* we have equilibrated weights */
1104 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1106 if (expand->elamstats == elamstatsMINVAR)
1108 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1112 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1116 fprintf(outfile, " <<\n");
1120 fprintf(outfile, " \n");
1123 fprintf(outfile, "\n");
1125 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1127 fprintf(outfile, " Transition Matrix\n");
1128 for (ifep = 0; ifep < nlim; ifep++)
1130 fprintf(outfile, "%12d", (ifep+1));
1132 fprintf(outfile, "\n");
1133 for (ifep = 0; ifep < nlim; ifep++)
1135 for (jfep = 0; jfep < nlim; jfep++)
1137 if (dfhist->n_at_lam[ifep] > 0)
1139 if (expand->bSymmetrizedTMatrix)
1141 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1145 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1152 fprintf(outfile, "%12.8f", Tprint);
1154 fprintf(outfile, "%3d\n", (ifep+1));
1157 fprintf(outfile, " Empirical Transition Matrix\n");
1158 for (ifep = 0; ifep < nlim; ifep++)
1160 fprintf(outfile, "%12d", (ifep+1));
1162 fprintf(outfile, "\n");
1163 for (ifep = 0; ifep < nlim; ifep++)
1165 for (jfep = 0; jfep < nlim; jfep++)
1167 if (dfhist->n_at_lam[ifep] > 0)
1169 if (expand->bSymmetrizedTMatrix)
1171 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1175 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1182 fprintf(outfile, "%12.8f", Tprint);
1184 fprintf(outfile, "%3d\n", (ifep+1));
1190 extern void get_mc_state(gmx_rng_t rng, t_state *state)
1192 gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
1195 extern void set_mc_state(gmx_rng_t rng, t_state *state)
1197 gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
1200 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1201 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
1202 gmx_large_int_t step, gmx_rng_t mcrng,
1203 rvec *v, t_mdatoms *mdatoms)
1205 real *pfep_lamee, *p_k, *scaled_lamee, *weighted_lamee;
1206 int i, nlam, nlim, lamnew, totalsamples;
1207 real oneovert, maxscaled = 0, maxweighted = 0;
1210 double *temperature_lambdas;
1211 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1213 expand = ir->expandedvals;
1214 simtemp = ir->simtempvals;
1215 nlim = ir->fepvals->n_lambda;
1216 nlam = state->fep_state;
1218 snew(scaled_lamee, nlim);
1219 snew(weighted_lamee, nlim);
1220 snew(pfep_lamee, nlim);
1223 if (expand->bInit_weights) /* if initialized weights, we need to fill them in */
1225 dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
1226 for (i = 0; i < nlim; i++)
1228 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
1229 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
1231 expand->bInit_weights = FALSE;
1234 /* update the count at the current lambda*/
1235 dfhist->n_at_lam[nlam]++;
1237 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1238 pressure controlled.*/
1241 where does this PV term go?
1242 for (i=0;i<nlim;i++)
1244 fep_lamee[i] += pVTerm;
1248 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1249 /* we don't need to include the pressure term, since the volume is the same between the two.
1250 is there some term we are neglecting, however? */
1252 if (ir->efep != efepNO)
1254 for (i = 0; i < nlim; i++)
1258 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1259 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1260 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
1264 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1265 /* mc_temp is currently set to the system reft unless otherwise defined */
1268 /* save these energies for printing, so they don't get overwritten by the next step */
1269 /* they aren't overwritten in the non-free energy case, but we always print with these
1277 for (i = 0; i < nlim; i++)
1279 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
1284 for (i = 0; i < nlim; i++)
1286 pfep_lamee[i] = scaled_lamee[i];
1288 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1291 maxscaled = scaled_lamee[i];
1292 maxweighted = weighted_lamee[i];
1296 if (scaled_lamee[i] > maxscaled)
1298 maxscaled = scaled_lamee[i];
1300 if (weighted_lamee[i] > maxweighted)
1302 maxweighted = weighted_lamee[i];
1307 for (i = 0; i < nlim; i++)
1309 scaled_lamee[i] -= maxscaled;
1310 weighted_lamee[i] -= maxweighted;
1313 /* update weights - we decide whether or not to actually do this inside */
1315 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, nlam, scaled_lamee, weighted_lamee, step);
1316 if (bDoneEquilibrating)
1320 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1324 lamnew = ChooseNewLambda(log, nlim, expand, dfhist, nlam, weighted_lamee, p_k, mcrng);
1325 /* if using simulated tempering, we need to adjust the temperatures */
1326 if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
1331 int nstart, nend, gt;
1333 snew(buf_ngtc, ir->opts.ngtc);
1335 for (i = 0; i < ir->opts.ngtc; i++)
1337 if (ir->opts.ref_t[i] > 0)
1339 told = ir->opts.ref_t[i];
1340 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1341 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1345 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1347 nstart = mdatoms->start;
1348 nend = nstart + mdatoms->homenr;
1349 for (n = nstart; n < nend; n++)
1354 gt = mdatoms->cTC[n];
1356 for (d = 0; d < DIM; d++)
1358 v[n][d] *= buf_ngtc[gt];
1362 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1364 /* we need to recalculate the masses if the temperature has changed */
1365 init_npt_masses(ir, state, MassQ, FALSE);
1366 for (i = 0; i < state->nnhpres; i++)
1368 for (j = 0; j < ir->opts.nhchainlength; j++)
1370 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1373 for (i = 0; i < ir->opts.ngtc; i++)
1375 for (j = 0; j < ir->opts.nhchainlength; j++)
1377 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1384 /* now check on the Wang-Landau updating critera */
1386 if (EWL(expand->elamstats))
1388 bSwitchtoOneOverT = FALSE;
1389 if (expand->bWLoneovert)
1392 for (i = 0; i < nlim; i++)
1394 totalsamples += dfhist->n_at_lam[i];
1396 oneovert = (1.0*nlim)/totalsamples;
1397 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1398 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1399 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1400 (dfhist->wl_delta < expand->init_wl_delta))
1402 bSwitchtoOneOverT = TRUE;
1405 if (bSwitchtoOneOverT)
1407 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1411 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1414 for (i = 0; i < nlim; i++)
1416 dfhist->wl_histo[i] = 0;
1418 dfhist->wl_delta *= expand->wl_scale;
1421 fprintf(log, "\nStep %d: weights are now:", (int)step);
1422 for (i = 0; i < nlim; i++)
1424 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1431 sfree(scaled_lamee);
1432 sfree(weighted_lamee);