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.
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/utility/smalloc.h"
42 #include "gromacs/legacyheaders/names.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/legacyheaders/txtdump.h"
45 #include "gromacs/legacyheaders/chargegroup.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/legacyheaders/mdrun.h"
49 #include "gromacs/legacyheaders/update.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/legacyheaders/mdatoms.h"
52 #include "gromacs/legacyheaders/force.h"
53 #include "gromacs/legacyheaders/pme.h"
54 #include "gromacs/legacyheaders/disre.h"
55 #include "gromacs/legacyheaders/orires.h"
56 #include "gromacs/legacyheaders/network.h"
57 #include "gromacs/legacyheaders/calcmu.h"
58 #include "gromacs/legacyheaders/constr.h"
59 #include "gromacs/random/random.h"
60 #include "gromacs/legacyheaders/domdec.h"
61 #include "gromacs/legacyheaders/macros.h"
63 #include "gromacs/fileio/confio.h"
64 #include "gromacs/fileio/gmxfio.h"
65 #include "gromacs/fileio/trnio.h"
66 #include "gromacs/fileio/xtcio.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxmpi.h"
71 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
74 dfhist->wl_delta = expand->init_wl_delta;
75 for (i = 0; i < nlim; i++)
77 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
78 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
82 /* Eventually should contain all the functions needed to initialize expanded ensemble
83 before the md loop starts */
84 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
88 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
92 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
100 /* find the maximum value */
101 for (i = minfep; i <= maxfep; i++)
108 /* find the denominator */
109 for (i = minfep; i <= maxfep; i++)
111 *pks += exp(ene[i]-maxene);
114 for (i = minfep; i <= maxfep; i++)
116 p_k[i] = exp(ene[i]-maxene) / *pks;
120 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
129 for (i = 0; i < nlim; i++)
133 /* add the delta, since we need to make sure it's greater than zero, and
134 we need a non-arbitrary number? */
135 nene[i] = ene[i] + log(nvals[i]+delta);
139 nene[i] = ene[i] + log(nvals[i]);
143 /* find the maximum value */
145 for (i = 0; i < nlim; i++)
147 if (nene[i] > maxene)
153 /* subtract off the maximum, avoiding overflow */
154 for (i = 0; i < nlim; i++)
159 /* find the denominator */
160 for (i = 0; i < nlim; i++)
162 *pks += exp(nene[i]);
166 for (i = 0; i < nlim; i++)
168 p_k[i] = exp(nene[i]) / *pks;
173 real do_logsum(int N, real *a_n)
177 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
182 /* compute maximum argument to exp(.) */
185 for (i = 1; i < N; i++)
187 maxarg = max(maxarg, a_n[i]);
190 /* compute sum of exp(a_n - maxarg) */
192 for (i = 0; i < N; i++)
194 sum = sum + exp(a_n[i] - maxarg);
197 /* compute log sum */
198 logsum = log(sum) + maxarg;
202 int FindMinimum(real *min_metric, int N)
209 min_val = min_metric[0];
211 for (nval = 0; nval < N; nval++)
213 if (min_metric[nval] < min_val)
215 min_val = min_metric[nval];
222 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
230 for (i = 0; i < nhisto; i++)
237 /* no samples! is bad!*/
241 nmean /= (real)nhisto;
244 for (i = 0; i < nhisto; i++)
246 /* make sure that all points are in the ratio < x < 1/ratio range */
247 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
256 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
260 gmx_bool bDoneEquilibrating = TRUE;
263 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
265 /* calculate the total number of samples */
266 switch (expand->elmceq)
269 /* We have not equilibrated, and won't, ever. */
272 /* we have equilibrated -- we're done */
275 /* first, check if we are equilibrating by steps, if we're still under */
276 if (step < expand->equil_steps)
278 bDoneEquilibrating = FALSE;
283 for (i = 0; i < nlim; i++)
285 totalsamples += dfhist->n_at_lam[i];
287 if (totalsamples < expand->equil_samples)
289 bDoneEquilibrating = FALSE;
293 for (i = 0; i < nlim; i++)
295 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
298 bDoneEquilibrating = FALSE;
304 if (EWL(expand->elamstats)) /* This check is in readir as well, but
307 if (dfhist->wl_delta > expand->equil_wl_delta)
309 bDoneEquilibrating = FALSE;
314 /* we can use the flatness as a judge of good weights, as long as
315 we're not doing minvar, or Wang-Landau.
316 But turn off for now until we figure out exactly how we do this.
319 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
321 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
322 floats for this histogram function. */
325 snew(modhisto, nlim);
326 for (i = 0; i < nlim; i++)
328 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
330 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
334 bDoneEquilibrating = FALSE;
338 bDoneEquilibrating = TRUE;
340 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
342 if (expand->lmc_forced_nstart > 0)
344 for (i = 0; i < nlim; i++)
346 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
349 bDoneEquilibrating = FALSE;
354 return bDoneEquilibrating;
357 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
358 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
360 real maxdiff = 0.000000001;
361 gmx_bool bSufficientSamples;
362 int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
363 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
364 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
365 real de, de_function, dr, denom, maxdr;
366 real min_val, cnval, zero_sum_weights;
367 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
368 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
369 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
372 real *numweighted_lamee, *logfrac;
374 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;
376 /* if we have equilibrated the weights, exit now */
382 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
384 dfhist->bEquil = TRUE;
385 /* zero out the visited states so we know how many equilibrated states we have
387 for (i = 0; i < nlim; i++)
389 dfhist->n_at_lam[i] = 0;
394 /* If we reached this far, we have not equilibrated yet, keep on
395 going resetting the weights */
397 if (EWL(expand->elamstats))
399 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
401 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
402 dfhist->wl_histo[fep_state] += 1.0;
404 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
408 /* first increment count */
409 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
410 for (i = 0; i < nlim; i++)
412 dfhist->wl_histo[i] += (real)p_k[i];
415 /* then increment weights (uses count) */
417 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
419 for (i = 0; i < nlim; i++)
421 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
423 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
428 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
429 dfhist->sum_weights[i] -= log(di);
435 zero_sum_weights = dfhist->sum_weights[0];
436 for (i = 0; i < nlim; i++)
438 dfhist->sum_weights[i] -= zero_sum_weights;
442 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
445 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
446 maxc = 2*expand->c_range+1;
449 snew(lam_variance, nlim);
451 snew(omegap_array, maxc);
452 snew(weightsp_array, maxc);
453 snew(varp_array, maxc);
454 snew(dwp_array, maxc);
456 snew(omegam_array, maxc);
457 snew(weightsm_array, maxc);
458 snew(varm_array, maxc);
459 snew(dwm_array, maxc);
461 /* unpack the current lambdas -- we will only update 2 of these */
463 for (i = 0; i < nlim-1; i++)
464 { /* only through the second to last */
465 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
466 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
469 /* accumulate running averages */
470 for (nval = 0; nval < maxc; nval++)
472 /* constants for later use */
473 cnval = (real)(nval-expand->c_range);
474 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
477 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
478 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
480 de_function = 1.0/(1.0+de);
482 else if (expand->elamstats == elamstatsMETROPOLIS)
490 de_function = 1.0/de;
493 dfhist->accum_m[fep_state][nval] += de_function;
494 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
497 if (fep_state < nlim-1)
499 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
500 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
502 de_function = 1.0/(1.0+de);
504 else if (expand->elamstats == elamstatsMETROPOLIS)
512 de_function = 1.0/de;
515 dfhist->accum_p[fep_state][nval] += de_function;
516 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
519 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
521 n0 = dfhist->n_at_lam[fep_state];
524 nm1 = dfhist->n_at_lam[fep_state-1];
530 if (fep_state < nlim-1)
532 np1 = dfhist->n_at_lam[fep_state+1];
539 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
540 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;
544 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
545 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
546 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
547 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
550 if ((fep_state > 0 ) && (nm1 > 0))
552 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
553 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
556 if ((fep_state < nlim-1) && (np1 > 0))
558 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
559 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
573 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
577 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
579 if ((n0 > 0) && (nm1 > 0))
581 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
582 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
586 if (fep_state < nlim-1)
590 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
594 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
596 if ((n0 > 0) && (np1 > 0))
598 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
599 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
605 omegam_array[nval] = omega_m1_0;
609 omegam_array[nval] = 0;
611 weightsm_array[nval] = clam_weightsm;
612 varm_array[nval] = clam_varm;
615 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
619 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
624 omegap_array[nval] = omega_p1_0;
628 omegap_array[nval] = 0;
630 weightsp_array[nval] = clam_weightsp;
631 varp_array[nval] = clam_varp;
632 if ((np1 > 0) && (n0 > 0))
634 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
638 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
643 /* find the C's closest to the old weights value */
645 min_nvalm = FindMinimum(dwm_array, maxc);
646 omega_m1_0 = omegam_array[min_nvalm];
647 clam_weightsm = weightsm_array[min_nvalm];
648 clam_varm = varm_array[min_nvalm];
650 min_nvalp = FindMinimum(dwp_array, maxc);
651 omega_p1_0 = omegap_array[min_nvalp];
652 clam_weightsp = weightsp_array[min_nvalp];
653 clam_varp = varp_array[min_nvalp];
655 clam_osum = omega_m1_0 + omega_p1_0;
659 clam_minvar = 0.5*log(clam_osum);
664 lam_dg[fep_state-1] = clam_weightsm;
665 lam_variance[fep_state-1] = clam_varm;
668 if (fep_state < nlim-1)
670 lam_dg[fep_state] = clam_weightsp;
671 lam_variance[fep_state] = clam_varp;
674 if (expand->elamstats == elamstatsMINVAR)
676 bSufficientSamples = TRUE;
677 /* make sure they are all past a threshold */
678 for (i = 0; i < nlim; i++)
680 if (dfhist->n_at_lam[i] < expand->minvarmin)
682 bSufficientSamples = FALSE;
685 if (bSufficientSamples)
687 dfhist->sum_minvar[fep_state] = clam_minvar;
690 for (i = 0; i < nlim; i++)
692 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
694 expand->minvar_const = clam_minvar;
695 dfhist->sum_minvar[fep_state] = 0.0;
699 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
704 /* we need to rezero minvar now, since it could change at fep_state = 0 */
705 dfhist->sum_dg[0] = 0.0;
706 dfhist->sum_variance[0] = 0.0;
707 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
709 for (i = 1; i < nlim; i++)
711 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
712 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
713 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
720 sfree(weightsm_array);
725 sfree(weightsp_array);
732 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
733 gmx_int64_t seed, gmx_int64_t step)
735 /* Choose new lambda value, and update transition matrix */
737 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
738 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
740 double *propose, *accept, *remainder;
743 gmx_bool bRestricted;
745 starting_fep_state = fep_state;
746 lamnew = fep_state; /* so that there is a default setting -- stays the same */
748 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
750 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
752 /* Use a marching method to run through the lambdas and get preliminary free energy data,
753 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
755 /* if we have enough at this lambda, move on to the next one */
757 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
759 lamnew = fep_state+1;
760 if (lamnew == nlim) /* whoops, stepped too far! */
775 snew(remainder, nlim);
777 for (i = 0; i < expand->lmc_repeats; i++)
781 gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
783 for (ifep = 0; ifep < nlim; ifep++)
789 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
792 /* use the Gibbs sampler, with restricted range */
793 if (expand->gibbsdeltalam < 0)
801 minfep = fep_state - expand->gibbsdeltalam;
802 maxfep = fep_state + expand->gibbsdeltalam;
813 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
815 if (expand->elmcmove == elmcmoveGIBBS)
817 for (ifep = minfep; ifep <= maxfep; ifep++)
819 propose[ifep] = p_k[ifep];
824 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
826 if (r1 <= p_k[lamnew])
833 else if (expand->elmcmove == elmcmoveMETGIBBS)
836 /* Metropolized Gibbs sampling */
837 for (ifep = minfep; ifep <= maxfep; ifep++)
839 remainder[ifep] = 1 - p_k[ifep];
842 /* find the proposal probabilities */
844 if (remainder[fep_state] == 0)
846 /* only the current state has any probability */
847 /* we have to stay at the current state */
852 for (ifep = minfep; ifep <= maxfep; ifep++)
854 if (ifep != fep_state)
856 propose[ifep] = p_k[ifep]/remainder[fep_state];
865 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
867 pnorm = p_k[lamtrial]/remainder[fep_state];
868 if (lamtrial != fep_state)
878 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
880 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
881 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
882 if (trialprob < tprob)
897 /* now figure out the acceptance probability for each */
898 for (ifep = minfep; ifep <= maxfep; ifep++)
901 if (remainder[ifep] != 0)
903 trialprob = (remainder[fep_state])/(remainder[ifep]);
907 trialprob = 1.0; /* this state is the only choice! */
909 if (trialprob < tprob)
913 /* probability for fep_state=0, but that's fine, it's never proposed! */
914 accept[ifep] = tprob;
920 /* it's possible some rounding is failing */
921 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
923 /* numerical rounding error -- no state other than the original has weight */
928 /* probably not a numerical issue */
930 int nerror = 200+(maxfep-minfep+1)*60;
932 snew(errorstr, nerror);
933 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
934 of sum weights. Generated detailed info for failure */
935 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);
936 for (ifep = minfep; ifep <= maxfep; ifep++)
938 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
940 gmx_fatal(FARGS, errorstr);
944 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
946 /* use the metropolis sampler with trial +/- 1 */
952 lamtrial = fep_state;
956 lamtrial = fep_state-1;
961 if (fep_state == nlim-1)
963 lamtrial = fep_state;
967 lamtrial = fep_state+1;
971 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
972 if (expand->elmcmove == elmcmoveMETROPOLIS)
976 if (trialprob < tprob)
980 propose[fep_state] = 0;
981 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
982 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
983 accept[lamtrial] = tprob;
986 else if (expand->elmcmove == elmcmoveBARKER)
988 tprob = 1.0/(1.0+exp(-de));
990 propose[fep_state] = (1-tprob);
991 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
992 accept[fep_state] = 1.0;
993 accept[lamtrial] = 1.0;
1007 for (ifep = 0; ifep < nlim; ifep++)
1009 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1010 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1015 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1024 /* print out the weights to the log, along with current state */
1025 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1026 int fep_state, int frequency, gmx_int64_t step)
1028 int nlim, i, ifep, jfep;
1029 real dw, dg, dv, dm, Tprint;
1031 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1032 gmx_bool bSimTemp = FALSE;
1034 nlim = fep->n_lambda;
1035 if (simtemp != NULL)
1040 if (mod(step, frequency) == 0)
1042 fprintf(outfile, " MC-lambda information\n");
1043 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1045 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1047 fprintf(outfile, " N");
1048 for (i = 0; i < efptNR; i++)
1050 if (fep->separate_dvdl[i])
1052 fprintf(outfile, "%7s", print_names[i]);
1054 else if ((i == efptTEMPERATURE) && bSimTemp)
1056 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1059 fprintf(outfile, " Count ");
1060 if (expand->elamstats == elamstatsMINVAR)
1062 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1066 fprintf(outfile, "G(in kT) dG(in kT)\n");
1068 for (ifep = 0; ifep < nlim; ifep++)
1079 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1080 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1081 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1082 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1085 fprintf(outfile, "%3d", (ifep+1));
1086 for (i = 0; i < efptNR; i++)
1088 if (fep->separate_dvdl[i])
1090 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1092 else if (i == efptTEMPERATURE && bSimTemp)
1094 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1097 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1099 if (expand->elamstats == elamstatsWL)
1101 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1105 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1108 else /* we have equilibrated weights */
1110 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1112 if (expand->elamstats == elamstatsMINVAR)
1114 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1118 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1120 if (ifep == fep_state)
1122 fprintf(outfile, " <<\n");
1126 fprintf(outfile, " \n");
1129 fprintf(outfile, "\n");
1131 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1133 fprintf(outfile, " Transition Matrix\n");
1134 for (ifep = 0; ifep < nlim; ifep++)
1136 fprintf(outfile, "%12d", (ifep+1));
1138 fprintf(outfile, "\n");
1139 for (ifep = 0; ifep < nlim; ifep++)
1141 for (jfep = 0; jfep < nlim; jfep++)
1143 if (dfhist->n_at_lam[ifep] > 0)
1145 if (expand->bSymmetrizedTMatrix)
1147 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1151 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1158 fprintf(outfile, "%12.8f", Tprint);
1160 fprintf(outfile, "%3d\n", (ifep+1));
1163 fprintf(outfile, " Empirical Transition Matrix\n");
1164 for (ifep = 0; ifep < nlim; ifep++)
1166 fprintf(outfile, "%12d", (ifep+1));
1168 fprintf(outfile, "\n");
1169 for (ifep = 0; ifep < nlim; ifep++)
1171 for (jfep = 0; jfep < nlim; jfep++)
1173 if (dfhist->n_at_lam[ifep] > 0)
1175 if (expand->bSymmetrizedTMatrix)
1177 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1181 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1188 fprintf(outfile, "%12.8f", Tprint);
1190 fprintf(outfile, "%3d\n", (ifep+1));
1196 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1197 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1199 rvec *v, t_mdatoms *mdatoms)
1200 /* Note that the state variable is only needed for simulated tempering, not
1201 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1203 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1205 int i, nlim, lamnew, totalsamples;
1206 real oneovert, maxscaled = 0, maxweighted = 0;
1209 double *temperature_lambdas;
1210 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1212 expand = ir->expandedvals;
1213 simtemp = ir->simtempvals;
1214 nlim = ir->fepvals->n_lambda;
1216 snew(scaled_lamee, nlim);
1217 snew(weighted_lamee, nlim);
1218 snew(pfep_lamee, nlim);
1221 /* update the count at the current lambda*/
1222 dfhist->n_at_lam[fep_state]++;
1224 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1225 pressure controlled.*/
1228 where does this PV term go?
1229 for (i=0;i<nlim;i++)
1231 fep_lamee[i] += pVTerm;
1235 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1236 /* we don't need to include the pressure term, since the volume is the same between the two.
1237 is there some term we are neglecting, however? */
1239 if (ir->efep != efepNO)
1241 for (i = 0; i < nlim; i++)
1245 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1246 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1247 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1251 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1252 /* mc_temp is currently set to the system reft unless otherwise defined */
1255 /* save these energies for printing, so they don't get overwritten by the next step */
1256 /* they aren't overwritten in the non-free energy case, but we always print with these
1264 for (i = 0; i < nlim; i++)
1266 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1271 for (i = 0; i < nlim; i++)
1273 pfep_lamee[i] = scaled_lamee[i];
1275 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1278 maxscaled = scaled_lamee[i];
1279 maxweighted = weighted_lamee[i];
1283 if (scaled_lamee[i] > maxscaled)
1285 maxscaled = scaled_lamee[i];
1287 if (weighted_lamee[i] > maxweighted)
1289 maxweighted = weighted_lamee[i];
1294 for (i = 0; i < nlim; i++)
1296 scaled_lamee[i] -= maxscaled;
1297 weighted_lamee[i] -= maxweighted;
1300 /* update weights - we decide whether or not to actually do this inside */
1302 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1303 if (bDoneEquilibrating)
1307 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1311 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1312 ir->expandedvals->lmc_seed, step);
1313 /* if using simulated tempering, we need to adjust the temperatures */
1314 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1319 int nstart, nend, gt;
1321 snew(buf_ngtc, ir->opts.ngtc);
1323 for (i = 0; i < ir->opts.ngtc; i++)
1325 if (ir->opts.ref_t[i] > 0)
1327 told = ir->opts.ref_t[i];
1328 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1329 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1333 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1336 nend = mdatoms->homenr;
1337 for (n = nstart; n < nend; n++)
1342 gt = mdatoms->cTC[n];
1344 for (d = 0; d < DIM; d++)
1346 v[n][d] *= buf_ngtc[gt];
1350 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1352 /* we need to recalculate the masses if the temperature has changed */
1353 init_npt_masses(ir, state, MassQ, FALSE);
1354 for (i = 0; i < state->nnhpres; i++)
1356 for (j = 0; j < ir->opts.nhchainlength; j++)
1358 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1361 for (i = 0; i < ir->opts.ngtc; i++)
1363 for (j = 0; j < ir->opts.nhchainlength; j++)
1365 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1372 /* now check on the Wang-Landau updating critera */
1374 if (EWL(expand->elamstats))
1376 bSwitchtoOneOverT = FALSE;
1377 if (expand->bWLoneovert)
1380 for (i = 0; i < nlim; i++)
1382 totalsamples += dfhist->n_at_lam[i];
1384 oneovert = (1.0*nlim)/totalsamples;
1385 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1386 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1387 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1388 (dfhist->wl_delta < expand->init_wl_delta))
1390 bSwitchtoOneOverT = TRUE;
1393 if (bSwitchtoOneOverT)
1395 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1399 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1402 for (i = 0; i < nlim; i++)
1404 dfhist->wl_histo[i] = 0;
1406 dfhist->wl_delta *= expand->wl_scale;
1409 fprintf(log, "\nStep %d: weights are now:", (int)step);
1410 for (i = 0; i < nlim; i++)
1412 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1420 sfree(scaled_lamee);
1421 sfree(weighted_lamee);