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.
43 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/fileio/confio.h"
47 #include "chargegroup.h"
48 #include "gromacs/math/vec.h"
52 #include "gromacs/math/units.h"
62 #include "gromacs/random/random.h"
66 #include "gromacs/fileio/confio.h"
67 #include "gromacs/fileio/gmxfio.h"
68 #include "gromacs/fileio/trnio.h"
69 #include "gromacs/fileio/xtcio.h"
70 #include "gromacs/timing/wallcycle.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxmpi.h"
74 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
77 dfhist->wl_delta = expand->init_wl_delta;
78 for (i = 0; i < nlim; i++)
80 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
81 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
85 /* Eventually should contain all the functions needed to initialize expanded ensemble
86 before the md loop starts */
87 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
91 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
95 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
102 maxene = ene[minfep];
103 /* find the maximum value */
104 for (i = minfep; i <= maxfep; i++)
111 /* find the denominator */
112 for (i = minfep; i <= maxfep; i++)
114 *pks += exp(ene[i]-maxene);
117 for (i = minfep; i <= maxfep; i++)
119 p_k[i] = exp(ene[i]-maxene) / *pks;
123 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
132 for (i = 0; i < nlim; i++)
136 /* add the delta, since we need to make sure it's greater than zero, and
137 we need a non-arbitrary number? */
138 nene[i] = ene[i] + log(nvals[i]+delta);
142 nene[i] = ene[i] + log(nvals[i]);
146 /* find the maximum value */
148 for (i = 0; i < nlim; i++)
150 if (nene[i] > maxene)
156 /* subtract off the maximum, avoiding overflow */
157 for (i = 0; i < nlim; i++)
162 /* find the denominator */
163 for (i = 0; i < nlim; i++)
165 *pks += exp(nene[i]);
169 for (i = 0; i < nlim; i++)
171 p_k[i] = exp(nene[i]) / *pks;
176 real do_logsum(int N, real *a_n)
180 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
185 /* compute maximum argument to exp(.) */
188 for (i = 1; i < N; i++)
190 maxarg = max(maxarg, a_n[i]);
193 /* compute sum of exp(a_n - maxarg) */
195 for (i = 0; i < N; i++)
197 sum = sum + exp(a_n[i] - maxarg);
200 /* compute log sum */
201 logsum = log(sum) + maxarg;
205 int FindMinimum(real *min_metric, int N)
212 min_val = min_metric[0];
214 for (nval = 0; nval < N; nval++)
216 if (min_metric[nval] < min_val)
218 min_val = min_metric[nval];
225 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
233 for (i = 0; i < nhisto; i++)
240 /* no samples! is bad!*/
244 nmean /= (real)nhisto;
247 for (i = 0; i < nhisto; i++)
249 /* make sure that all points are in the ratio < x < 1/ratio range */
250 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
259 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
263 gmx_bool bDoneEquilibrating = TRUE;
266 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
268 /* calculate the total number of samples */
269 switch (expand->elmceq)
272 /* We have not equilibrated, and won't, ever. */
275 /* we have equilibrated -- we're done */
278 /* first, check if we are equilibrating by steps, if we're still under */
279 if (step < expand->equil_steps)
281 bDoneEquilibrating = FALSE;
286 for (i = 0; i < nlim; i++)
288 totalsamples += dfhist->n_at_lam[i];
290 if (totalsamples < expand->equil_samples)
292 bDoneEquilibrating = FALSE;
296 for (i = 0; i < nlim; i++)
298 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
301 bDoneEquilibrating = FALSE;
307 if (EWL(expand->elamstats)) /* This check is in readir as well, but
310 if (dfhist->wl_delta > expand->equil_wl_delta)
312 bDoneEquilibrating = FALSE;
317 /* we can use the flatness as a judge of good weights, as long as
318 we're not doing minvar, or Wang-Landau.
319 But turn off for now until we figure out exactly how we do this.
322 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
324 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
325 floats for this histogram function. */
328 snew(modhisto, nlim);
329 for (i = 0; i < nlim; i++)
331 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
333 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
337 bDoneEquilibrating = FALSE;
341 bDoneEquilibrating = TRUE;
343 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
345 if (expand->lmc_forced_nstart > 0)
347 for (i = 0; i < nlim; i++)
349 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
352 bDoneEquilibrating = FALSE;
357 return bDoneEquilibrating;
360 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
361 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
363 real maxdiff = 0.000000001;
364 gmx_bool bSufficientSamples;
365 int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
366 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
367 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
368 real de, de_function, dr, denom, maxdr;
369 real min_val, cnval, zero_sum_weights;
370 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
371 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
372 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
375 real *numweighted_lamee, *logfrac;
377 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;
379 /* if we have equilibrated the weights, exit now */
385 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
387 dfhist->bEquil = TRUE;
388 /* zero out the visited states so we know how many equilibrated states we have
390 for (i = 0; i < nlim; i++)
392 dfhist->n_at_lam[i] = 0;
397 /* If we reached this far, we have not equilibrated yet, keep on
398 going resetting the weights */
400 if (EWL(expand->elamstats))
402 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
404 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
405 dfhist->wl_histo[fep_state] += 1.0;
407 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
411 /* first increment count */
412 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
413 for (i = 0; i < nlim; i++)
415 dfhist->wl_histo[i] += (real)p_k[i];
418 /* then increment weights (uses count) */
420 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
422 for (i = 0; i < nlim; i++)
424 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
426 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
431 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
432 dfhist->sum_weights[i] -= log(di);
438 zero_sum_weights = dfhist->sum_weights[0];
439 for (i = 0; i < nlim; i++)
441 dfhist->sum_weights[i] -= zero_sum_weights;
445 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
448 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
449 maxc = 2*expand->c_range+1;
452 snew(lam_variance, nlim);
454 snew(omegap_array, maxc);
455 snew(weightsp_array, maxc);
456 snew(varp_array, maxc);
457 snew(dwp_array, maxc);
459 snew(omegam_array, maxc);
460 snew(weightsm_array, maxc);
461 snew(varm_array, maxc);
462 snew(dwm_array, maxc);
464 /* unpack the current lambdas -- we will only update 2 of these */
466 for (i = 0; i < nlim-1; i++)
467 { /* only through the second to last */
468 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
469 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
472 /* accumulate running averages */
473 for (nval = 0; nval < maxc; nval++)
475 /* constants for later use */
476 cnval = (real)(nval-expand->c_range);
477 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
480 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
481 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
483 de_function = 1.0/(1.0+de);
485 else if (expand->elamstats == elamstatsMETROPOLIS)
493 de_function = 1.0/de;
496 dfhist->accum_m[fep_state][nval] += de_function;
497 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
500 if (fep_state < nlim-1)
502 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
503 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
505 de_function = 1.0/(1.0+de);
507 else if (expand->elamstats == elamstatsMETROPOLIS)
515 de_function = 1.0/de;
518 dfhist->accum_p[fep_state][nval] += de_function;
519 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
522 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
524 n0 = dfhist->n_at_lam[fep_state];
527 nm1 = dfhist->n_at_lam[fep_state-1];
533 if (fep_state < nlim-1)
535 np1 = dfhist->n_at_lam[fep_state+1];
542 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
543 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;
547 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
548 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
549 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
550 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
553 if ((fep_state > 0 ) && (nm1 > 0))
555 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
556 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
559 if ((fep_state < nlim-1) && (np1 > 0))
561 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
562 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
576 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
580 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
582 if ((n0 > 0) && (nm1 > 0))
584 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
585 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
589 if (fep_state < nlim-1)
593 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
597 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
599 if ((n0 > 0) && (np1 > 0))
601 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
602 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
608 omegam_array[nval] = omega_m1_0;
612 omegam_array[nval] = 0;
614 weightsm_array[nval] = clam_weightsm;
615 varm_array[nval] = clam_varm;
618 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
622 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
627 omegap_array[nval] = omega_p1_0;
631 omegap_array[nval] = 0;
633 weightsp_array[nval] = clam_weightsp;
634 varp_array[nval] = clam_varp;
635 if ((np1 > 0) && (n0 > 0))
637 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
641 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
646 /* find the C's closest to the old weights value */
648 min_nvalm = FindMinimum(dwm_array, maxc);
649 omega_m1_0 = omegam_array[min_nvalm];
650 clam_weightsm = weightsm_array[min_nvalm];
651 clam_varm = varm_array[min_nvalm];
653 min_nvalp = FindMinimum(dwp_array, maxc);
654 omega_p1_0 = omegap_array[min_nvalp];
655 clam_weightsp = weightsp_array[min_nvalp];
656 clam_varp = varp_array[min_nvalp];
658 clam_osum = omega_m1_0 + omega_p1_0;
662 clam_minvar = 0.5*log(clam_osum);
667 lam_dg[fep_state-1] = clam_weightsm;
668 lam_variance[fep_state-1] = clam_varm;
671 if (fep_state < nlim-1)
673 lam_dg[fep_state] = clam_weightsp;
674 lam_variance[fep_state] = clam_varp;
677 if (expand->elamstats == elamstatsMINVAR)
679 bSufficientSamples = TRUE;
680 /* make sure they are all past a threshold */
681 for (i = 0; i < nlim; i++)
683 if (dfhist->n_at_lam[i] < expand->minvarmin)
685 bSufficientSamples = FALSE;
688 if (bSufficientSamples)
690 dfhist->sum_minvar[fep_state] = clam_minvar;
693 for (i = 0; i < nlim; i++)
695 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
697 expand->minvar_const = clam_minvar;
698 dfhist->sum_minvar[fep_state] = 0.0;
702 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
707 /* we need to rezero minvar now, since it could change at fep_state = 0 */
708 dfhist->sum_dg[0] = 0.0;
709 dfhist->sum_variance[0] = 0.0;
710 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
712 for (i = 1; i < nlim; i++)
714 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
715 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
716 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
723 sfree(weightsm_array);
728 sfree(weightsp_array);
735 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
736 gmx_int64_t seed, gmx_int64_t step)
738 /* Choose new lambda value, and update transition matrix */
740 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
741 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
743 double *propose, *accept, *remainder;
746 gmx_bool bRestricted;
748 starting_fep_state = fep_state;
749 lamnew = fep_state; /* so that there is a default setting -- stays the same */
751 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
753 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
755 /* Use a marching method to run through the lambdas and get preliminary free energy data,
756 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
758 /* if we have enough at this lambda, move on to the next one */
760 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
762 lamnew = fep_state+1;
763 if (lamnew == nlim) /* whoops, stepped too far! */
778 snew(remainder, nlim);
780 for (i = 0; i < expand->lmc_repeats; i++)
784 gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
786 for (ifep = 0; ifep < nlim; ifep++)
792 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
795 /* use the Gibbs sampler, with restricted range */
796 if (expand->gibbsdeltalam < 0)
804 minfep = fep_state - expand->gibbsdeltalam;
805 maxfep = fep_state + expand->gibbsdeltalam;
816 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
818 if (expand->elmcmove == elmcmoveGIBBS)
820 for (ifep = minfep; ifep <= maxfep; ifep++)
822 propose[ifep] = p_k[ifep];
827 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
829 if (r1 <= p_k[lamnew])
836 else if (expand->elmcmove == elmcmoveMETGIBBS)
839 /* Metropolized Gibbs sampling */
840 for (ifep = minfep; ifep <= maxfep; ifep++)
842 remainder[ifep] = 1 - p_k[ifep];
845 /* find the proposal probabilities */
847 if (remainder[fep_state] == 0)
849 /* only the current state has any probability */
850 /* we have to stay at the current state */
855 for (ifep = minfep; ifep <= maxfep; ifep++)
857 if (ifep != fep_state)
859 propose[ifep] = p_k[ifep]/remainder[fep_state];
868 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
870 pnorm = p_k[lamtrial]/remainder[fep_state];
871 if (lamtrial != fep_state)
881 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
883 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
884 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
885 if (trialprob < tprob)
900 /* now figure out the acceptance probability for each */
901 for (ifep = minfep; ifep <= maxfep; ifep++)
904 if (remainder[ifep] != 0)
906 trialprob = (remainder[fep_state])/(remainder[ifep]);
910 trialprob = 1.0; /* this state is the only choice! */
912 if (trialprob < tprob)
916 /* probability for fep_state=0, but that's fine, it's never proposed! */
917 accept[ifep] = tprob;
923 /* it's possible some rounding is failing */
924 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
926 /* numerical rounding error -- no state other than the original has weight */
931 /* probably not a numerical issue */
933 int nerror = 200+(maxfep-minfep+1)*60;
935 snew(errorstr, nerror);
936 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
937 of sum weights. Generated detailed info for failure */
938 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);
939 for (ifep = minfep; ifep <= maxfep; ifep++)
941 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
943 gmx_fatal(FARGS, errorstr);
947 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
949 /* use the metropolis sampler with trial +/- 1 */
955 lamtrial = fep_state;
959 lamtrial = fep_state-1;
964 if (fep_state == nlim-1)
966 lamtrial = fep_state;
970 lamtrial = fep_state+1;
974 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
975 if (expand->elmcmove == elmcmoveMETROPOLIS)
979 if (trialprob < tprob)
983 propose[fep_state] = 0;
984 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
985 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
986 accept[lamtrial] = tprob;
989 else if (expand->elmcmove == elmcmoveBARKER)
991 tprob = 1.0/(1.0+exp(-de));
993 propose[fep_state] = (1-tprob);
994 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
995 accept[fep_state] = 1.0;
996 accept[lamtrial] = 1.0;
1010 for (ifep = 0; ifep < nlim; ifep++)
1012 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1013 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1018 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1027 /* print out the weights to the log, along with current state */
1028 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1029 int fep_state, int frequency, gmx_int64_t step)
1031 int nlim, i, ifep, jfep;
1032 real dw, dg, dv, dm, Tprint;
1034 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1035 gmx_bool bSimTemp = FALSE;
1037 nlim = fep->n_lambda;
1038 if (simtemp != NULL)
1043 if (mod(step, frequency) == 0)
1045 fprintf(outfile, " MC-lambda information\n");
1046 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1048 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1050 fprintf(outfile, " N");
1051 for (i = 0; i < efptNR; i++)
1053 if (fep->separate_dvdl[i])
1055 fprintf(outfile, "%7s", print_names[i]);
1057 else if ((i == efptTEMPERATURE) && bSimTemp)
1059 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1062 fprintf(outfile, " Count ");
1063 if (expand->elamstats == elamstatsMINVAR)
1065 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1069 fprintf(outfile, "G(in kT) dG(in kT)\n");
1071 for (ifep = 0; ifep < nlim; ifep++)
1082 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1083 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1084 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1085 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1088 fprintf(outfile, "%3d", (ifep+1));
1089 for (i = 0; i < efptNR; i++)
1091 if (fep->separate_dvdl[i])
1093 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1095 else if (i == efptTEMPERATURE && bSimTemp)
1097 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1100 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1102 if (expand->elamstats == elamstatsWL)
1104 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1108 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1111 else /* we have equilibrated weights */
1113 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1115 if (expand->elamstats == elamstatsMINVAR)
1117 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1121 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1123 if (ifep == fep_state)
1125 fprintf(outfile, " <<\n");
1129 fprintf(outfile, " \n");
1132 fprintf(outfile, "\n");
1134 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1136 fprintf(outfile, " Transition Matrix\n");
1137 for (ifep = 0; ifep < nlim; ifep++)
1139 fprintf(outfile, "%12d", (ifep+1));
1141 fprintf(outfile, "\n");
1142 for (ifep = 0; ifep < nlim; ifep++)
1144 for (jfep = 0; jfep < nlim; jfep++)
1146 if (dfhist->n_at_lam[ifep] > 0)
1148 if (expand->bSymmetrizedTMatrix)
1150 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1154 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1161 fprintf(outfile, "%12.8f", Tprint);
1163 fprintf(outfile, "%3d\n", (ifep+1));
1166 fprintf(outfile, " Empirical Transition Matrix\n");
1167 for (ifep = 0; ifep < nlim; ifep++)
1169 fprintf(outfile, "%12d", (ifep+1));
1171 fprintf(outfile, "\n");
1172 for (ifep = 0; ifep < nlim; ifep++)
1174 for (jfep = 0; jfep < nlim; jfep++)
1176 if (dfhist->n_at_lam[ifep] > 0)
1178 if (expand->bSymmetrizedTMatrix)
1180 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1184 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1191 fprintf(outfile, "%12.8f", Tprint);
1193 fprintf(outfile, "%3d\n", (ifep+1));
1199 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1200 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1202 rvec *v, t_mdatoms *mdatoms)
1203 /* Note that the state variable is only needed for simulated tempering, not
1204 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1206 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1208 int i, nlim, lamnew, totalsamples;
1209 real oneovert, maxscaled = 0, maxweighted = 0;
1212 double *temperature_lambdas;
1213 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1215 expand = ir->expandedvals;
1216 simtemp = ir->simtempvals;
1217 nlim = ir->fepvals->n_lambda;
1219 snew(scaled_lamee, nlim);
1220 snew(weighted_lamee, nlim);
1221 snew(pfep_lamee, nlim);
1224 /* update the count at the current lambda*/
1225 dfhist->n_at_lam[fep_state]++;
1227 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1228 pressure controlled.*/
1231 where does this PV term go?
1232 for (i=0;i<nlim;i++)
1234 fep_lamee[i] += pVTerm;
1238 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1239 /* we don't need to include the pressure term, since the volume is the same between the two.
1240 is there some term we are neglecting, however? */
1242 if (ir->efep != efepNO)
1244 for (i = 0; i < nlim; i++)
1248 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1249 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1250 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1254 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1255 /* mc_temp is currently set to the system reft unless otherwise defined */
1258 /* save these energies for printing, so they don't get overwritten by the next step */
1259 /* they aren't overwritten in the non-free energy case, but we always print with these
1267 for (i = 0; i < nlim; i++)
1269 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1274 for (i = 0; i < nlim; i++)
1276 pfep_lamee[i] = scaled_lamee[i];
1278 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1281 maxscaled = scaled_lamee[i];
1282 maxweighted = weighted_lamee[i];
1286 if (scaled_lamee[i] > maxscaled)
1288 maxscaled = scaled_lamee[i];
1290 if (weighted_lamee[i] > maxweighted)
1292 maxweighted = weighted_lamee[i];
1297 for (i = 0; i < nlim; i++)
1299 scaled_lamee[i] -= maxscaled;
1300 weighted_lamee[i] -= maxweighted;
1303 /* update weights - we decide whether or not to actually do this inside */
1305 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1306 if (bDoneEquilibrating)
1310 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1314 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1315 ir->expandedvals->lmc_seed, step);
1316 /* if using simulated tempering, we need to adjust the temperatures */
1317 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1322 int nstart, nend, gt;
1324 snew(buf_ngtc, ir->opts.ngtc);
1326 for (i = 0; i < ir->opts.ngtc; i++)
1328 if (ir->opts.ref_t[i] > 0)
1330 told = ir->opts.ref_t[i];
1331 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1332 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1336 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1339 nend = mdatoms->homenr;
1340 for (n = nstart; n < nend; n++)
1345 gt = mdatoms->cTC[n];
1347 for (d = 0; d < DIM; d++)
1349 v[n][d] *= buf_ngtc[gt];
1353 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1355 /* we need to recalculate the masses if the temperature has changed */
1356 init_npt_masses(ir, state, MassQ, FALSE);
1357 for (i = 0; i < state->nnhpres; i++)
1359 for (j = 0; j < ir->opts.nhchainlength; j++)
1361 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1364 for (i = 0; i < ir->opts.ngtc; i++)
1366 for (j = 0; j < ir->opts.nhchainlength; j++)
1368 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1375 /* now check on the Wang-Landau updating critera */
1377 if (EWL(expand->elamstats))
1379 bSwitchtoOneOverT = FALSE;
1380 if (expand->bWLoneovert)
1383 for (i = 0; i < nlim; i++)
1385 totalsamples += dfhist->n_at_lam[i];
1387 oneovert = (1.0*nlim)/totalsamples;
1388 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1389 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1390 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1391 (dfhist->wl_delta < expand->init_wl_delta))
1393 bSwitchtoOneOverT = TRUE;
1396 if (bSwitchtoOneOverT)
1398 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1402 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1405 for (i = 0; i < nlim; i++)
1407 dfhist->wl_histo[i] = 0;
1409 dfhist->wl_delta *= expand->wl_scale;
1412 fprintf(log, "\nStep %d: weights are now:", (int)step);
1413 for (i = 0; i < nlim; i++)
1415 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1423 sfree(scaled_lamee);
1424 sfree(weighted_lamee);