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, double *p_k, double *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, double *p_k, double *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;
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;
367 real *numweighted_lamee, *logfrac;
369 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;
371 /* if we have equilibrated the weights, exit now */
377 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
379 dfhist->bEquil = TRUE;
380 /* zero out the visited states so we know how many equilibrated states we have
382 for (i = 0; i < nlim; i++)
384 dfhist->n_at_lam[i] = 0;
389 /* If we reached this far, we have not equilibrated yet, keep on
390 going resetting the weights */
392 if (EWL(expand->elamstats))
394 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
396 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
397 dfhist->wl_histo[fep_state] += 1.0;
399 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
403 /* first increment count */
404 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
405 for (i = 0; i < nlim; i++)
407 dfhist->wl_histo[i] += (real)p_k[i];
410 /* then increment weights (uses count) */
412 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
414 for (i = 0; i < nlim; i++)
416 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
418 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
423 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
424 dfhist->sum_weights[i] -= log(di);
430 zero_sum_weights = dfhist->sum_weights[0];
431 for (i = 0; i < nlim; i++)
433 dfhist->sum_weights[i] -= zero_sum_weights;
437 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
440 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
441 maxc = 2*expand->c_range+1;
444 snew(lam_variance, nlim);
446 snew(omegap_array, maxc);
447 snew(weightsp_array, maxc);
448 snew(varp_array, maxc);
449 snew(dwp_array, maxc);
451 snew(omegam_array, maxc);
452 snew(weightsm_array, maxc);
453 snew(varm_array, maxc);
454 snew(dwm_array, maxc);
456 /* unpack the current lambdas -- we will only update 2 of these */
458 for (i = 0; i < nlim-1; i++)
459 { /* only through the second to last */
460 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
461 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
464 /* accumulate running averages */
465 for (nval = 0; nval < maxc; nval++)
467 /* constants for later use */
468 cnval = (real)(nval-expand->c_range);
469 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
472 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
473 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
475 de_function = 1.0/(1.0+de);
477 else if (expand->elamstats == elamstatsMETROPOLIS)
485 de_function = 1.0/de;
488 dfhist->accum_m[fep_state][nval] += de_function;
489 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
492 if (fep_state < nlim-1)
494 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
495 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
497 de_function = 1.0/(1.0+de);
499 else if (expand->elamstats == elamstatsMETROPOLIS)
507 de_function = 1.0/de;
510 dfhist->accum_p[fep_state][nval] += de_function;
511 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
514 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
516 n0 = dfhist->n_at_lam[fep_state];
519 nm1 = dfhist->n_at_lam[fep_state-1];
525 if (fep_state < nlim-1)
527 np1 = dfhist->n_at_lam[fep_state+1];
534 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
535 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;
539 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
540 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
541 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
542 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
545 if ((fep_state > 0 ) && (nm1 > 0))
547 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
548 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
551 if ((fep_state < nlim-1) && (np1 > 0))
553 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
554 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
568 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
572 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
574 if ((n0 > 0) && (nm1 > 0))
576 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
577 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
581 if (fep_state < nlim-1)
585 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
589 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
591 if ((n0 > 0) && (np1 > 0))
593 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
594 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
600 omegam_array[nval] = omega_m1_0;
604 omegam_array[nval] = 0;
606 weightsm_array[nval] = clam_weightsm;
607 varm_array[nval] = clam_varm;
610 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
614 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
619 omegap_array[nval] = omega_p1_0;
623 omegap_array[nval] = 0;
625 weightsp_array[nval] = clam_weightsp;
626 varp_array[nval] = clam_varp;
627 if ((np1 > 0) && (n0 > 0))
629 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
633 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
638 /* find the C's closest to the old weights value */
640 min_nvalm = FindMinimum(dwm_array, maxc);
641 omega_m1_0 = omegam_array[min_nvalm];
642 clam_weightsm = weightsm_array[min_nvalm];
643 clam_varm = varm_array[min_nvalm];
645 min_nvalp = FindMinimum(dwp_array, maxc);
646 omega_p1_0 = omegap_array[min_nvalp];
647 clam_weightsp = weightsp_array[min_nvalp];
648 clam_varp = varp_array[min_nvalp];
650 clam_osum = omega_m1_0 + omega_p1_0;
654 clam_minvar = 0.5*log(clam_osum);
659 lam_dg[fep_state-1] = clam_weightsm;
660 lam_variance[fep_state-1] = clam_varm;
663 if (fep_state < nlim-1)
665 lam_dg[fep_state] = clam_weightsp;
666 lam_variance[fep_state] = clam_varp;
669 if (expand->elamstats == elamstatsMINVAR)
671 bSufficientSamples = TRUE;
672 /* make sure they are all past a threshold */
673 for (i = 0; i < nlim; i++)
675 if (dfhist->n_at_lam[i] < expand->minvarmin)
677 bSufficientSamples = FALSE;
680 if (bSufficientSamples)
682 dfhist->sum_minvar[fep_state] = clam_minvar;
685 for (i = 0; i < nlim; i++)
687 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
689 expand->minvar_const = clam_minvar;
690 dfhist->sum_minvar[fep_state] = 0.0;
694 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
699 /* we need to rezero minvar now, since it could change at fep_state = 0 */
700 dfhist->sum_dg[0] = 0.0;
701 dfhist->sum_variance[0] = 0.0;
702 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
704 for (i = 1; i < nlim; i++)
706 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
707 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
708 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
715 sfree(weightsm_array);
720 sfree(weightsp_array);
727 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k, gmx_rng_t rng)
729 /* Choose new lambda value, and update transition matrix */
731 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
732 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
734 double *propose, *accept, *remainder;
737 gmx_bool bRestricted;
739 starting_fep_state = fep_state;
740 lamnew = fep_state; /* so that there is a default setting -- stays the same */
742 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
744 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
746 /* Use a marching method to run through the lambdas and get preliminary free energy data,
747 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
749 /* if we have enough at this lambda, move on to the next one */
751 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
753 lamnew = fep_state+1;
754 if (lamnew == nlim) /* whoops, stepped too far! */
769 snew(remainder, nlim);
771 for (i = 0; i < expand->lmc_repeats; i++)
774 for (ifep = 0; ifep < nlim; ifep++)
780 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
783 /* use the Gibbs sampler, with restricted range */
784 if (expand->gibbsdeltalam < 0)
792 minfep = fep_state - expand->gibbsdeltalam;
793 maxfep = fep_state + expand->gibbsdeltalam;
804 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
806 if (expand->elmcmove == elmcmoveGIBBS)
808 for (ifep = minfep; ifep <= maxfep; ifep++)
810 propose[ifep] = p_k[ifep];
814 r1 = gmx_rng_uniform_real(rng);
815 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
817 if (r1 <= p_k[lamnew])
824 else if (expand->elmcmove == elmcmoveMETGIBBS)
827 /* Metropolized Gibbs sampling */
828 for (ifep = minfep; ifep <= maxfep; ifep++)
830 remainder[ifep] = 1 - p_k[ifep];
833 /* find the proposal probabilities */
835 if (remainder[fep_state] == 0)
837 /* only the current state has any probability */
838 /* we have to stay at the current state */
843 for (ifep = minfep; ifep <= maxfep; ifep++)
845 if (ifep != fep_state)
847 propose[ifep] = p_k[ifep]/remainder[fep_state];
855 r1 = gmx_rng_uniform_real(rng);
856 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
858 pnorm = p_k[lamtrial]/remainder[fep_state];
859 if (lamtrial != fep_state)
869 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
871 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
872 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
873 if (trialprob < tprob)
877 r2 = gmx_rng_uniform_real(rng);
888 /* now figure out the acceptance probability for each */
889 for (ifep = minfep; ifep <= maxfep; ifep++)
892 if (remainder[ifep] != 0)
894 trialprob = (remainder[fep_state])/(remainder[ifep]);
898 trialprob = 1.0; /* this state is the only choice! */
900 if (trialprob < tprob)
904 /* probability for fep_state=0, but that's fine, it's never proposed! */
905 accept[ifep] = tprob;
911 /* it's possible some rounding is failing */
912 if (gmx_within_tol(remainder[fep_state],0,50*GMX_DOUBLE_EPS))
914 /* numerical rounding error -- no state other than the original has weight */
919 /* probably not a numerical issue */
921 int nerror = 200+(maxfep-minfep+1)*60;
923 snew(errorstr, nerror);
924 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
925 of sum weights. Generated detailed info for failure */
926 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);
927 for (ifep = minfep; ifep <= maxfep; ifep++)
929 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
931 gmx_fatal(FARGS, errorstr);
935 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
937 /* use the metropolis sampler with trial +/- 1 */
938 r1 = gmx_rng_uniform_real(rng);
943 lamtrial = fep_state;
947 lamtrial = fep_state-1;
952 if (fep_state == nlim-1)
954 lamtrial = fep_state;
958 lamtrial = fep_state+1;
962 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
963 if (expand->elmcmove == elmcmoveMETROPOLIS)
967 if (trialprob < tprob)
971 propose[fep_state] = 0;
972 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
973 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
974 accept[lamtrial] = tprob;
977 else if (expand->elmcmove == elmcmoveBARKER)
979 tprob = 1.0/(1.0+exp(-de));
981 propose[fep_state] = (1-tprob);
982 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
983 accept[fep_state] = 1.0;
984 accept[lamtrial] = 1.0;
987 r2 = gmx_rng_uniform_real(rng);
998 for (ifep = 0; ifep < nlim; ifep++)
1000 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1001 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1006 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1015 /* print out the weights to the log, along with current state */
1016 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1017 int nlam, int frequency, gmx_large_int_t step)
1019 int nlim, i, ifep, jfep;
1020 real dw, dg, dv, dm, Tprint;
1022 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1023 gmx_bool bSimTemp = FALSE;
1025 nlim = fep->n_lambda;
1026 if (simtemp != NULL)
1031 if (mod(step, frequency) == 0)
1033 fprintf(outfile, " MC-lambda information\n");
1034 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1036 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1038 fprintf(outfile, " N");
1039 for (i = 0; i < efptNR; i++)
1041 if (fep->separate_dvdl[i])
1043 fprintf(outfile, "%7s", print_names[i]);
1045 else if ((i == efptTEMPERATURE) && bSimTemp)
1047 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1050 fprintf(outfile, " Count ");
1051 if (expand->elamstats == elamstatsMINVAR)
1053 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1057 fprintf(outfile, "G(in kT) dG(in kT)\n");
1059 for (ifep = 0; ifep < nlim; ifep++)
1070 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1071 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1072 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1073 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1076 fprintf(outfile, "%3d", (ifep+1));
1077 for (i = 0; i < efptNR; i++)
1079 if (fep->separate_dvdl[i])
1081 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1083 else if (i == efptTEMPERATURE && bSimTemp)
1085 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1088 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1090 if (expand->elamstats == elamstatsWL)
1092 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1096 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1099 else /* we have equilibrated weights */
1101 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1103 if (expand->elamstats == elamstatsMINVAR)
1105 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1109 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1113 fprintf(outfile, " <<\n");
1117 fprintf(outfile, " \n");
1120 fprintf(outfile, "\n");
1122 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1124 fprintf(outfile, " Transition Matrix\n");
1125 for (ifep = 0; ifep < nlim; ifep++)
1127 fprintf(outfile, "%12d", (ifep+1));
1129 fprintf(outfile, "\n");
1130 for (ifep = 0; ifep < nlim; ifep++)
1132 for (jfep = 0; jfep < nlim; jfep++)
1134 if (dfhist->n_at_lam[ifep] > 0)
1136 if (expand->bSymmetrizedTMatrix)
1138 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1142 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1149 fprintf(outfile, "%12.8f", Tprint);
1151 fprintf(outfile, "%3d\n", (ifep+1));
1154 fprintf(outfile, " Empirical Transition Matrix\n");
1155 for (ifep = 0; ifep < nlim; ifep++)
1157 fprintf(outfile, "%12d", (ifep+1));
1159 fprintf(outfile, "\n");
1160 for (ifep = 0; ifep < nlim; ifep++)
1162 for (jfep = 0; jfep < nlim; jfep++)
1164 if (dfhist->n_at_lam[ifep] > 0)
1166 if (expand->bSymmetrizedTMatrix)
1168 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1172 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1179 fprintf(outfile, "%12.8f", Tprint);
1181 fprintf(outfile, "%3d\n", (ifep+1));
1187 extern void get_mc_state(gmx_rng_t rng, t_state *state)
1189 gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
1192 extern void set_mc_state(gmx_rng_t rng, t_state *state)
1194 gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
1197 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1198 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
1199 gmx_large_int_t step, gmx_rng_t mcrng,
1200 rvec *v, t_mdatoms *mdatoms)
1202 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1204 int i, nlam, nlim, lamnew, totalsamples;
1205 real oneovert, maxscaled = 0, maxweighted = 0;
1208 double *temperature_lambdas;
1209 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1211 expand = ir->expandedvals;
1212 simtemp = ir->simtempvals;
1213 nlim = ir->fepvals->n_lambda;
1214 nlam = state->fep_state;
1216 snew(scaled_lamee, nlim);
1217 snew(weighted_lamee, nlim);
1218 snew(pfep_lamee, nlim);
1221 if (expand->bInit_weights) /* if initialized weights, we need to fill them in */
1223 dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
1224 for (i = 0; i < nlim; i++)
1226 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
1227 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
1229 expand->bInit_weights = FALSE;
1232 /* update the count at the current lambda*/
1233 dfhist->n_at_lam[nlam]++;
1235 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1236 pressure controlled.*/
1239 where does this PV term go?
1240 for (i=0;i<nlim;i++)
1242 fep_lamee[i] += pVTerm;
1246 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1247 /* we don't need to include the pressure term, since the volume is the same between the two.
1248 is there some term we are neglecting, however? */
1250 if (ir->efep != efepNO)
1252 for (i = 0; i < nlim; i++)
1256 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1257 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1258 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
1262 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1263 /* mc_temp is currently set to the system reft unless otherwise defined */
1266 /* save these energies for printing, so they don't get overwritten by the next step */
1267 /* they aren't overwritten in the non-free energy case, but we always print with these
1275 for (i = 0; i < nlim; i++)
1277 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
1282 for (i = 0; i < nlim; i++)
1284 pfep_lamee[i] = scaled_lamee[i];
1286 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1289 maxscaled = scaled_lamee[i];
1290 maxweighted = weighted_lamee[i];
1294 if (scaled_lamee[i] > maxscaled)
1296 maxscaled = scaled_lamee[i];
1298 if (weighted_lamee[i] > maxweighted)
1300 maxweighted = weighted_lamee[i];
1305 for (i = 0; i < nlim; i++)
1307 scaled_lamee[i] -= maxscaled;
1308 weighted_lamee[i] -= maxweighted;
1311 /* update weights - we decide whether or not to actually do this inside */
1313 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, nlam, scaled_lamee, weighted_lamee, step);
1314 if (bDoneEquilibrating)
1318 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1322 lamnew = ChooseNewLambda(nlim, expand, dfhist, nlam, weighted_lamee, p_k, mcrng);
1323 /* if using simulated tempering, we need to adjust the temperatures */
1324 if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
1329 int nstart, nend, gt;
1331 snew(buf_ngtc, ir->opts.ngtc);
1333 for (i = 0; i < ir->opts.ngtc; i++)
1335 if (ir->opts.ref_t[i] > 0)
1337 told = ir->opts.ref_t[i];
1338 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1339 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1343 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1345 nstart = mdatoms->start;
1346 nend = nstart + mdatoms->homenr;
1347 for (n = nstart; n < nend; n++)
1352 gt = mdatoms->cTC[n];
1354 for (d = 0; d < DIM; d++)
1356 v[n][d] *= buf_ngtc[gt];
1360 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1362 /* we need to recalculate the masses if the temperature has changed */
1363 init_npt_masses(ir, state, MassQ, FALSE);
1364 for (i = 0; i < state->nnhpres; i++)
1366 for (j = 0; j < ir->opts.nhchainlength; j++)
1368 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1371 for (i = 0; i < ir->opts.ngtc; i++)
1373 for (j = 0; j < ir->opts.nhchainlength; j++)
1375 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1382 /* now check on the Wang-Landau updating critera */
1384 if (EWL(expand->elamstats))
1386 bSwitchtoOneOverT = FALSE;
1387 if (expand->bWLoneovert)
1390 for (i = 0; i < nlim; i++)
1392 totalsamples += dfhist->n_at_lam[i];
1394 oneovert = (1.0*nlim)/totalsamples;
1395 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1396 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1397 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1398 (dfhist->wl_delta < expand->init_wl_delta))
1400 bSwitchtoOneOverT = TRUE;
1403 if (bSwitchtoOneOverT)
1405 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1409 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1412 for (i = 0; i < nlim; i++)
1414 dfhist->wl_histo[i] = 0;
1416 dfhist->wl_delta *= expand->wl_scale;
1419 fprintf(log, "\nStep %d: weights are now:", (int)step);
1420 for (i = 0; i < nlim; i++)
1422 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1430 sfree(scaled_lamee);
1431 sfree(weighted_lamee);