2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "chargegroup.h"
66 #include "gromacs/random/random.h"
71 #include "gromacs/fileio/confio.h"
72 #include "gromacs/fileio/gmxfio.h"
73 #include "gromacs/fileio/trnio.h"
74 #include "gromacs/fileio/xtcio.h"
75 #include "gromacs/timing/wallcycle.h"
76 #include "gromacs/utility/gmxmpi.h"
78 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
81 dfhist->wl_delta = expand->init_wl_delta;
82 for (i = 0; i < nlim; i++)
84 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
85 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
89 /* Eventually should contain all the functions needed to initialize expanded ensemble
90 before the md loop starts */
91 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, gmx_rng_t *mcrng, df_history_t *dfhist)
93 *mcrng = gmx_rng_init(ir->expandedvals->lmc_seed);
96 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
100 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
107 maxene = ene[minfep];
108 /* find the maximum value */
109 for (i = minfep; i <= maxfep; i++)
116 /* find the denominator */
117 for (i = minfep; i <= maxfep; i++)
119 *pks += exp(ene[i]-maxene);
122 for (i = minfep; i <= maxfep; i++)
124 p_k[i] = exp(ene[i]-maxene) / *pks;
128 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
137 for (i = 0; i < nlim; i++)
141 /* add the delta, since we need to make sure it's greater than zero, and
142 we need a non-arbitrary number? */
143 nene[i] = ene[i] + log(nvals[i]+delta);
147 nene[i] = ene[i] + log(nvals[i]);
151 /* find the maximum value */
153 for (i = 0; i < nlim; i++)
155 if (nene[i] > maxene)
161 /* subtract off the maximum, avoiding overflow */
162 for (i = 0; i < nlim; i++)
167 /* find the denominator */
168 for (i = 0; i < nlim; i++)
170 *pks += exp(nene[i]);
174 for (i = 0; i < nlim; i++)
176 p_k[i] = exp(nene[i]) / *pks;
181 real do_logsum(int N, real *a_n)
185 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
190 /* compute maximum argument to exp(.) */
193 for (i = 1; i < N; i++)
195 maxarg = max(maxarg, a_n[i]);
198 /* compute sum of exp(a_n - maxarg) */
200 for (i = 0; i < N; i++)
202 sum = sum + exp(a_n[i] - maxarg);
205 /* compute log sum */
206 logsum = log(sum) + maxarg;
210 int FindMinimum(real *min_metric, int N)
217 min_val = min_metric[0];
219 for (nval = 0; nval < N; nval++)
221 if (min_metric[nval] < min_val)
223 min_val = min_metric[nval];
230 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
238 for (i = 0; i < nhisto; i++)
245 /* no samples! is bad!*/
249 nmean /= (real)nhisto;
252 for (i = 0; i < nhisto; i++)
254 /* make sure that all points are in the ratio < x < 1/ratio range */
255 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
264 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
268 gmx_bool bDoneEquilibrating = TRUE;
271 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
273 /* calculate the total number of samples */
274 switch (expand->elmceq)
277 /* We have not equilibrated, and won't, ever. */
280 /* we have equilibrated -- we're done */
283 /* first, check if we are equilibrating by steps, if we're still under */
284 if (step < expand->equil_steps)
286 bDoneEquilibrating = FALSE;
291 for (i = 0; i < nlim; i++)
293 totalsamples += dfhist->n_at_lam[i];
295 if (totalsamples < expand->equil_samples)
297 bDoneEquilibrating = FALSE;
301 for (i = 0; i < nlim; i++)
303 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
306 bDoneEquilibrating = FALSE;
312 if (EWL(expand->elamstats)) /* This check is in readir as well, but
315 if (dfhist->wl_delta > expand->equil_wl_delta)
317 bDoneEquilibrating = FALSE;
322 /* we can use the flatness as a judge of good weights, as long as
323 we're not doing minvar, or Wang-Landau.
324 But turn off for now until we figure out exactly how we do this.
327 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
329 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
330 floats for this histogram function. */
333 snew(modhisto, nlim);
334 for (i = 0; i < nlim; i++)
336 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
338 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
342 bDoneEquilibrating = FALSE;
346 bDoneEquilibrating = TRUE;
348 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
350 if (expand->lmc_forced_nstart > 0)
352 for (i = 0; i < nlim; i++)
354 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
357 bDoneEquilibrating = FALSE;
362 return bDoneEquilibrating;
365 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
366 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
368 real maxdiff = 0.000000001;
369 gmx_bool bSufficientSamples;
370 int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
371 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
372 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
373 real de, de_function, dr, denom, maxdr;
374 real min_val, cnval, zero_sum_weights;
375 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
376 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
377 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
380 real *numweighted_lamee, *logfrac;
382 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;
384 /* if we have equilibrated the weights, exit now */
390 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
392 dfhist->bEquil = TRUE;
393 /* zero out the visited states so we know how many equilibrated states we have
395 for (i = 0; i < nlim; i++)
397 dfhist->n_at_lam[i] = 0;
402 /* If we reached this far, we have not equilibrated yet, keep on
403 going resetting the weights */
405 if (EWL(expand->elamstats))
407 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
409 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
410 dfhist->wl_histo[fep_state] += 1.0;
412 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
416 /* first increment count */
417 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
418 for (i = 0; i < nlim; i++)
420 dfhist->wl_histo[i] += (real)p_k[i];
423 /* then increment weights (uses count) */
425 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
427 for (i = 0; i < nlim; i++)
429 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
431 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
436 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
437 dfhist->sum_weights[i] -= log(di);
443 zero_sum_weights = dfhist->sum_weights[0];
444 for (i = 0; i < nlim; i++)
446 dfhist->sum_weights[i] -= zero_sum_weights;
450 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
453 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
454 maxc = 2*expand->c_range+1;
457 snew(lam_variance, nlim);
459 snew(omegap_array, maxc);
460 snew(weightsp_array, maxc);
461 snew(varp_array, maxc);
462 snew(dwp_array, maxc);
464 snew(omegam_array, maxc);
465 snew(weightsm_array, maxc);
466 snew(varm_array, maxc);
467 snew(dwm_array, maxc);
469 /* unpack the current lambdas -- we will only update 2 of these */
471 for (i = 0; i < nlim-1; i++)
472 { /* only through the second to last */
473 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
474 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
477 /* accumulate running averages */
478 for (nval = 0; nval < maxc; nval++)
480 /* constants for later use */
481 cnval = (real)(nval-expand->c_range);
482 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
485 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
486 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
488 de_function = 1.0/(1.0+de);
490 else if (expand->elamstats == elamstatsMETROPOLIS)
498 de_function = 1.0/de;
501 dfhist->accum_m[fep_state][nval] += de_function;
502 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
505 if (fep_state < nlim-1)
507 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
508 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
510 de_function = 1.0/(1.0+de);
512 else if (expand->elamstats == elamstatsMETROPOLIS)
520 de_function = 1.0/de;
523 dfhist->accum_p[fep_state][nval] += de_function;
524 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
527 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
529 n0 = dfhist->n_at_lam[fep_state];
532 nm1 = dfhist->n_at_lam[fep_state-1];
538 if (fep_state < nlim-1)
540 np1 = dfhist->n_at_lam[fep_state+1];
547 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
548 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;
552 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
553 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
554 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
555 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
558 if ((fep_state > 0 ) && (nm1 > 0))
560 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
561 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
564 if ((fep_state < nlim-1) && (np1 > 0))
566 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
567 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
581 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
585 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
587 if ((n0 > 0) && (nm1 > 0))
589 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
590 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
594 if (fep_state < nlim-1)
598 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
602 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
604 if ((n0 > 0) && (np1 > 0))
606 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
607 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
613 omegam_array[nval] = omega_m1_0;
617 omegam_array[nval] = 0;
619 weightsm_array[nval] = clam_weightsm;
620 varm_array[nval] = clam_varm;
623 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
627 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
632 omegap_array[nval] = omega_p1_0;
636 omegap_array[nval] = 0;
638 weightsp_array[nval] = clam_weightsp;
639 varp_array[nval] = clam_varp;
640 if ((np1 > 0) && (n0 > 0))
642 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
646 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
651 /* find the C's closest to the old weights value */
653 min_nvalm = FindMinimum(dwm_array, maxc);
654 omega_m1_0 = omegam_array[min_nvalm];
655 clam_weightsm = weightsm_array[min_nvalm];
656 clam_varm = varm_array[min_nvalm];
658 min_nvalp = FindMinimum(dwp_array, maxc);
659 omega_p1_0 = omegap_array[min_nvalp];
660 clam_weightsp = weightsp_array[min_nvalp];
661 clam_varp = varp_array[min_nvalp];
663 clam_osum = omega_m1_0 + omega_p1_0;
667 clam_minvar = 0.5*log(clam_osum);
672 lam_dg[fep_state-1] = clam_weightsm;
673 lam_variance[fep_state-1] = clam_varm;
676 if (fep_state < nlim-1)
678 lam_dg[fep_state] = clam_weightsp;
679 lam_variance[fep_state] = clam_varp;
682 if (expand->elamstats == elamstatsMINVAR)
684 bSufficientSamples = TRUE;
685 /* make sure they are all past a threshold */
686 for (i = 0; i < nlim; i++)
688 if (dfhist->n_at_lam[i] < expand->minvarmin)
690 bSufficientSamples = FALSE;
693 if (bSufficientSamples)
695 dfhist->sum_minvar[fep_state] = clam_minvar;
698 for (i = 0; i < nlim; i++)
700 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
702 expand->minvar_const = clam_minvar;
703 dfhist->sum_minvar[fep_state] = 0.0;
707 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
712 /* we need to rezero minvar now, since it could change at fep_state = 0 */
713 dfhist->sum_dg[0] = 0.0;
714 dfhist->sum_variance[0] = 0.0;
715 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
717 for (i = 1; i < nlim; i++)
719 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
720 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
721 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
728 sfree(weightsm_array);
733 sfree(weightsp_array);
740 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)
742 /* Choose new lambda value, and update transition matrix */
744 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
745 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
747 double *propose, *accept, *remainder;
750 gmx_bool bRestricted;
752 starting_fep_state = fep_state;
753 lamnew = fep_state; /* so that there is a default setting -- stays the same */
755 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
757 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
759 /* Use a marching method to run through the lambdas and get preliminary free energy data,
760 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
762 /* if we have enough at this lambda, move on to the next one */
764 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
766 lamnew = fep_state+1;
767 if (lamnew == nlim) /* whoops, stepped too far! */
782 snew(remainder, nlim);
784 for (i = 0; i < expand->lmc_repeats; i++)
787 for (ifep = 0; ifep < nlim; ifep++)
793 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
796 /* use the Gibbs sampler, with restricted range */
797 if (expand->gibbsdeltalam < 0)
805 minfep = fep_state - expand->gibbsdeltalam;
806 maxfep = fep_state + expand->gibbsdeltalam;
817 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
819 if (expand->elmcmove == elmcmoveGIBBS)
821 for (ifep = minfep; ifep <= maxfep; ifep++)
823 propose[ifep] = p_k[ifep];
827 r1 = gmx_rng_uniform_real(rng);
828 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
830 if (r1 <= p_k[lamnew])
837 else if (expand->elmcmove == elmcmoveMETGIBBS)
840 /* Metropolized Gibbs sampling */
841 for (ifep = minfep; ifep <= maxfep; ifep++)
843 remainder[ifep] = 1 - p_k[ifep];
846 /* find the proposal probabilities */
848 if (remainder[fep_state] == 0)
850 /* only the current state has any probability */
851 /* we have to stay at the current state */
856 for (ifep = minfep; ifep <= maxfep; ifep++)
858 if (ifep != fep_state)
860 propose[ifep] = p_k[ifep]/remainder[fep_state];
868 r1 = gmx_rng_uniform_real(rng);
869 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
871 pnorm = p_k[lamtrial]/remainder[fep_state];
872 if (lamtrial != fep_state)
882 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
884 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
885 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
886 if (trialprob < tprob)
890 r2 = gmx_rng_uniform_real(rng);
901 /* now figure out the acceptance probability for each */
902 for (ifep = minfep; ifep <= maxfep; ifep++)
905 if (remainder[ifep] != 0)
907 trialprob = (remainder[fep_state])/(remainder[ifep]);
911 trialprob = 1.0; /* this state is the only choice! */
913 if (trialprob < tprob)
917 /* probability for fep_state=0, but that's fine, it's never proposed! */
918 accept[ifep] = tprob;
924 /* it's possible some rounding is failing */
925 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
927 /* numerical rounding error -- no state other than the original has weight */
932 /* probably not a numerical issue */
934 int nerror = 200+(maxfep-minfep+1)*60;
936 snew(errorstr, nerror);
937 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
938 of sum weights. Generated detailed info for failure */
939 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);
940 for (ifep = minfep; ifep <= maxfep; ifep++)
942 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
944 gmx_fatal(FARGS, errorstr);
948 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
950 /* use the metropolis sampler with trial +/- 1 */
951 r1 = gmx_rng_uniform_real(rng);
956 lamtrial = fep_state;
960 lamtrial = fep_state-1;
965 if (fep_state == nlim-1)
967 lamtrial = fep_state;
971 lamtrial = fep_state+1;
975 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
976 if (expand->elmcmove == elmcmoveMETROPOLIS)
980 if (trialprob < tprob)
984 propose[fep_state] = 0;
985 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
986 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
987 accept[lamtrial] = tprob;
990 else if (expand->elmcmove == elmcmoveBARKER)
992 tprob = 1.0/(1.0+exp(-de));
994 propose[fep_state] = (1-tprob);
995 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
996 accept[fep_state] = 1.0;
997 accept[lamtrial] = 1.0;
1000 r2 = gmx_rng_uniform_real(rng);
1011 for (ifep = 0; ifep < nlim; ifep++)
1013 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1014 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1019 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1028 /* print out the weights to the log, along with current state */
1029 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1030 int fep_state, int frequency, gmx_int64_t step)
1032 int nlim, i, ifep, jfep;
1033 real dw, dg, dv, dm, Tprint;
1035 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1036 gmx_bool bSimTemp = FALSE;
1038 nlim = fep->n_lambda;
1039 if (simtemp != NULL)
1044 if (mod(step, frequency) == 0)
1046 fprintf(outfile, " MC-lambda information\n");
1047 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1049 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1051 fprintf(outfile, " N");
1052 for (i = 0; i < efptNR; i++)
1054 if (fep->separate_dvdl[i])
1056 fprintf(outfile, "%7s", print_names[i]);
1058 else if ((i == efptTEMPERATURE) && bSimTemp)
1060 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1063 fprintf(outfile, " Count ");
1064 if (expand->elamstats == elamstatsMINVAR)
1066 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1070 fprintf(outfile, "G(in kT) dG(in kT)\n");
1072 for (ifep = 0; ifep < nlim; ifep++)
1083 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1084 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1085 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1086 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1089 fprintf(outfile, "%3d", (ifep+1));
1090 for (i = 0; i < efptNR; i++)
1092 if (fep->separate_dvdl[i])
1094 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1096 else if (i == efptTEMPERATURE && bSimTemp)
1098 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1101 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1103 if (expand->elamstats == elamstatsWL)
1105 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1109 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1112 else /* we have equilibrated weights */
1114 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1116 if (expand->elamstats == elamstatsMINVAR)
1118 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1122 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1124 if (ifep == fep_state)
1126 fprintf(outfile, " <<\n");
1130 fprintf(outfile, " \n");
1133 fprintf(outfile, "\n");
1135 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1137 fprintf(outfile, " Transition Matrix\n");
1138 for (ifep = 0; ifep < nlim; ifep++)
1140 fprintf(outfile, "%12d", (ifep+1));
1142 fprintf(outfile, "\n");
1143 for (ifep = 0; ifep < nlim; ifep++)
1145 for (jfep = 0; jfep < nlim; jfep++)
1147 if (dfhist->n_at_lam[ifep] > 0)
1149 if (expand->bSymmetrizedTMatrix)
1151 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1155 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1162 fprintf(outfile, "%12.8f", Tprint);
1164 fprintf(outfile, "%3d\n", (ifep+1));
1167 fprintf(outfile, " Empirical Transition Matrix\n");
1168 for (ifep = 0; ifep < nlim; ifep++)
1170 fprintf(outfile, "%12d", (ifep+1));
1172 fprintf(outfile, "\n");
1173 for (ifep = 0; ifep < nlim; ifep++)
1175 for (jfep = 0; jfep < nlim; jfep++)
1177 if (dfhist->n_at_lam[ifep] > 0)
1179 if (expand->bSymmetrizedTMatrix)
1181 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1185 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1192 fprintf(outfile, "%12.8f", Tprint);
1194 fprintf(outfile, "%3d\n", (ifep+1));
1200 extern void get_mc_state(gmx_rng_t rng, t_state *state)
1202 gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
1205 extern void set_mc_state(gmx_rng_t rng, t_state *state)
1207 gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
1210 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1211 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1212 gmx_int64_t step, gmx_rng_t mcrng,
1213 rvec *v, t_mdatoms *mdatoms)
1214 /* Note that the state variable is only needed for simulated tempering, not
1215 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1217 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1219 int i, nlim, lamnew, totalsamples;
1220 real oneovert, maxscaled = 0, maxweighted = 0;
1223 double *temperature_lambdas;
1224 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1226 expand = ir->expandedvals;
1227 simtemp = ir->simtempvals;
1228 nlim = ir->fepvals->n_lambda;
1230 snew(scaled_lamee, nlim);
1231 snew(weighted_lamee, nlim);
1232 snew(pfep_lamee, nlim);
1235 /* update the count at the current lambda*/
1236 dfhist->n_at_lam[fep_state]++;
1238 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1239 pressure controlled.*/
1242 where does this PV term go?
1243 for (i=0;i<nlim;i++)
1245 fep_lamee[i] += pVTerm;
1249 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1250 /* we don't need to include the pressure term, since the volume is the same between the two.
1251 is there some term we are neglecting, however? */
1253 if (ir->efep != efepNO)
1255 for (i = 0; i < nlim; i++)
1259 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1260 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1261 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1265 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1266 /* mc_temp is currently set to the system reft unless otherwise defined */
1269 /* save these energies for printing, so they don't get overwritten by the next step */
1270 /* they aren't overwritten in the non-free energy case, but we always print with these
1278 for (i = 0; i < nlim; i++)
1280 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1285 for (i = 0; i < nlim; i++)
1287 pfep_lamee[i] = scaled_lamee[i];
1289 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1292 maxscaled = scaled_lamee[i];
1293 maxweighted = weighted_lamee[i];
1297 if (scaled_lamee[i] > maxscaled)
1299 maxscaled = scaled_lamee[i];
1301 if (weighted_lamee[i] > maxweighted)
1303 maxweighted = weighted_lamee[i];
1308 for (i = 0; i < nlim; i++)
1310 scaled_lamee[i] -= maxscaled;
1311 weighted_lamee[i] -= maxweighted;
1314 /* update weights - we decide whether or not to actually do this inside */
1316 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1317 if (bDoneEquilibrating)
1321 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1325 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k, mcrng);
1326 /* if using simulated tempering, we need to adjust the temperatures */
1327 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1332 int nstart, nend, gt;
1334 snew(buf_ngtc, ir->opts.ngtc);
1336 for (i = 0; i < ir->opts.ngtc; i++)
1338 if (ir->opts.ref_t[i] > 0)
1340 told = ir->opts.ref_t[i];
1341 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1342 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1346 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1348 nstart = mdatoms->start;
1349 nend = nstart + mdatoms->homenr;
1350 for (n = nstart; n < nend; n++)
1355 gt = mdatoms->cTC[n];
1357 for (d = 0; d < DIM; d++)
1359 v[n][d] *= buf_ngtc[gt];
1363 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1365 /* we need to recalculate the masses if the temperature has changed */
1366 init_npt_masses(ir, state, MassQ, FALSE);
1367 for (i = 0; i < state->nnhpres; i++)
1369 for (j = 0; j < ir->opts.nhchainlength; j++)
1371 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1374 for (i = 0; i < ir->opts.ngtc; i++)
1376 for (j = 0; j < ir->opts.nhchainlength; j++)
1378 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1385 /* now check on the Wang-Landau updating critera */
1387 if (EWL(expand->elamstats))
1389 bSwitchtoOneOverT = FALSE;
1390 if (expand->bWLoneovert)
1393 for (i = 0; i < nlim; i++)
1395 totalsamples += dfhist->n_at_lam[i];
1397 oneovert = (1.0*nlim)/totalsamples;
1398 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1399 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1400 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1401 (dfhist->wl_delta < expand->init_wl_delta))
1403 bSwitchtoOneOverT = TRUE;
1406 if (bSwitchtoOneOverT)
1408 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1412 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1415 for (i = 0; i < nlim; i++)
1417 dfhist->wl_histo[i] = 0;
1419 dfhist->wl_delta *= expand->wl_scale;
1422 fprintf(log, "\nStep %d: weights are now:", (int)step);
1423 for (i = 0; i < nlim; i++)
1425 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1433 sfree(scaled_lamee);
1434 sfree(weighted_lamee);