2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016, 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/domdec/domdec.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/fileio/xtcio.h"
47 #include "gromacs/gmxlib/chargegroup.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/listed-forces/disre.h"
51 #include "gromacs/listed-forces/orires.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/calcmu.h"
56 #include "gromacs/mdlib/constr.h"
57 #include "gromacs/mdlib/force.h"
58 #include "gromacs/mdlib/mdrun.h"
59 #include "gromacs/mdlib/update.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/random/threefry.h"
63 #include "gromacs/random/uniformrealdistribution.h"
64 #include "gromacs/timing/wallcycle.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/gmxmpi.h"
67 #include "gromacs/utility/smalloc.h"
69 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
72 dfhist->wl_delta = expand->init_wl_delta;
73 for (i = 0; i < nlim; i++)
75 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
76 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
80 /* Eventually should contain all the functions needed to initialize expanded ensemble
81 before the md loop starts */
82 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
86 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
90 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
98 /* find the maximum value */
99 for (i = minfep; i <= maxfep; i++)
106 /* find the denominator */
107 for (i = minfep; i <= maxfep; i++)
109 *pks += std::exp(ene[i]-maxene);
112 for (i = minfep; i <= maxfep; i++)
114 p_k[i] = std::exp(ene[i]-maxene) / *pks;
118 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
127 for (i = 0; i < nlim; i++)
131 /* add the delta, since we need to make sure it's greater than zero, and
132 we need a non-arbitrary number? */
133 nene[i] = ene[i] + std::log(nvals[i]+delta);
137 nene[i] = ene[i] + std::log(nvals[i]);
141 /* find the maximum value */
143 for (i = 0; i < nlim; i++)
145 if (nene[i] > maxene)
151 /* subtract off the maximum, avoiding overflow */
152 for (i = 0; i < nlim; i++)
157 /* find the denominator */
158 for (i = 0; i < nlim; i++)
160 *pks += std::exp(nene[i]);
164 for (i = 0; i < nlim; i++)
166 p_k[i] = std::exp(nene[i]) / *pks;
171 real do_logsum(int N, real *a_n)
175 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
180 /* compute maximum argument to exp(.) */
183 for (i = 1; i < N; i++)
185 maxarg = std::max(maxarg, a_n[i]);
188 /* compute sum of exp(a_n - maxarg) */
190 for (i = 0; i < N; i++)
192 sum = sum + std::exp(a_n[i] - maxarg);
195 /* compute log sum */
196 logsum = std::log(sum) + maxarg;
200 int FindMinimum(real *min_metric, int N)
207 min_val = min_metric[0];
209 for (nval = 0; nval < N; nval++)
211 if (min_metric[nval] < min_val)
213 min_val = min_metric[nval];
220 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
228 for (i = 0; i < nhisto; i++)
235 /* no samples! is bad!*/
239 nmean /= (real)nhisto;
242 for (i = 0; i < nhisto; i++)
244 /* make sure that all points are in the ratio < x < 1/ratio range */
245 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
254 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
258 gmx_bool bDoneEquilibrating = TRUE;
261 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
262 if (expand->lmc_forced_nstart > 0)
264 for (i = 0; i < nlim; i++)
266 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
269 bDoneEquilibrating = FALSE;
276 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
277 bDoneEquilibrating = TRUE;
279 /* calculate the total number of samples */
280 switch (expand->elmceq)
283 /* We have not equilibrated, and won't, ever. */
284 bDoneEquilibrating = FALSE;
287 /* we have equilibrated -- we're done */
288 bDoneEquilibrating = TRUE;
291 /* first, check if we are equilibrating by steps, if we're still under */
292 if (step < expand->equil_steps)
294 bDoneEquilibrating = FALSE;
299 for (i = 0; i < nlim; i++)
301 totalsamples += dfhist->n_at_lam[i];
303 if (totalsamples < expand->equil_samples)
305 bDoneEquilibrating = FALSE;
309 for (i = 0; i < nlim; i++)
311 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
314 bDoneEquilibrating = FALSE;
320 if (EWL(expand->elamstats)) /* This check is in readir as well, but
323 if (dfhist->wl_delta > expand->equil_wl_delta)
325 bDoneEquilibrating = FALSE;
330 /* we can use the flatness as a judge of good weights, as long as
331 we're not doing minvar, or Wang-Landau.
332 But turn off for now until we figure out exactly how we do this.
335 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
337 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
338 floats for this histogram function. */
341 snew(modhisto, nlim);
342 for (i = 0; i < nlim; i++)
344 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
346 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
350 bDoneEquilibrating = FALSE;
355 bDoneEquilibrating = TRUE;
359 return bDoneEquilibrating;
362 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
363 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
365 gmx_bool bSufficientSamples;
367 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
368 real omega_m1_0, omega_p1_0, clam_osum;
369 real de, de_function;
370 real cnval, zero_sum_weights;
371 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
372 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
373 real *lam_variance, *lam_dg;
376 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;
378 /* if we have equilibrated the weights, exit now */
384 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
386 dfhist->bEquil = TRUE;
387 /* zero out the visited states so we know how many equilibrated states we have
389 for (i = 0; i < nlim; i++)
391 dfhist->n_at_lam[i] = 0;
396 /* If we reached this far, we have not equilibrated yet, keep on
397 going resetting the weights */
399 if (EWL(expand->elamstats))
401 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
403 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
404 dfhist->wl_histo[fep_state] += 1.0;
406 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
410 /* first increment count */
411 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
412 for (i = 0; i < nlim; i++)
414 dfhist->wl_histo[i] += (real)p_k[i];
417 /* then increment weights (uses count) */
419 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
421 for (i = 0; i < nlim; i++)
423 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
425 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
430 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
431 dfhist->sum_weights[i] -= log(di);
437 zero_sum_weights = dfhist->sum_weights[0];
438 for (i = 0; i < nlim; i++)
440 dfhist->sum_weights[i] -= zero_sum_weights;
444 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
447 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
448 maxc = 2*expand->c_range+1;
451 snew(lam_variance, nlim);
453 snew(omegap_array, maxc);
454 snew(weightsp_array, maxc);
455 snew(varp_array, maxc);
456 snew(dwp_array, maxc);
458 snew(omegam_array, maxc);
459 snew(weightsm_array, maxc);
460 snew(varm_array, maxc);
461 snew(dwm_array, maxc);
463 /* unpack the current lambdas -- we will only update 2 of these */
465 for (i = 0; i < nlim-1; i++)
466 { /* only through the second to last */
467 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
468 lam_variance[i] = gmx::square(dfhist->sum_variance[i+1]) - gmx::square(dfhist->sum_variance[i]);
471 /* accumulate running averages */
472 for (nval = 0; nval < maxc; nval++)
474 /* constants for later use */
475 cnval = (real)(nval-expand->c_range);
476 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
479 de = std::exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
480 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
482 de_function = 1.0/(1.0+de);
484 else if (expand->elamstats == elamstatsMETROPOLIS)
492 de_function = 1.0/de;
495 dfhist->accum_m[fep_state][nval] += de_function;
496 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
499 if (fep_state < nlim-1)
501 de = std::exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
502 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
504 de_function = 1.0/(1.0+de);
506 else if (expand->elamstats == elamstatsMETROPOLIS)
514 de_function = 1.0/de;
517 dfhist->accum_p[fep_state][nval] += de_function;
518 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
521 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
523 n0 = dfhist->n_at_lam[fep_state];
526 nm1 = dfhist->n_at_lam[fep_state-1];
532 if (fep_state < nlim-1)
534 np1 = dfhist->n_at_lam[fep_state+1];
541 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
542 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;
546 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
547 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
548 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
549 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
552 if ((fep_state > 0 ) && (nm1 > 0))
554 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
555 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
558 if ((fep_state < nlim-1) && (np1 > 0))
560 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
561 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
575 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
578 real omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
579 clam_weightsm = (std::log(chi_m1_0) - std::log(chi_p1_m1)) + cnval;
580 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
585 if (fep_state < nlim-1)
589 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
592 real omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
593 clam_weightsp = (std::log(chi_m1_p1) - std::log(chi_p1_0)) + cnval;
594 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
601 omegam_array[nval] = omega_m1_0;
605 omegam_array[nval] = 0;
607 weightsm_array[nval] = clam_weightsm;
608 varm_array[nval] = clam_varm;
611 dwm_array[nval] = fabs( (cnval + std::log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
615 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
620 omegap_array[nval] = omega_p1_0;
624 omegap_array[nval] = 0;
626 weightsp_array[nval] = clam_weightsp;
627 varp_array[nval] = clam_varp;
628 if ((np1 > 0) && (n0 > 0))
630 dwp_array[nval] = fabs( (cnval + std::log((1.0*np1)/n0)) - lam_dg[fep_state] );
634 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
639 /* find the C's closest to the old weights value */
641 min_nvalm = FindMinimum(dwm_array, maxc);
642 omega_m1_0 = omegam_array[min_nvalm];
643 clam_weightsm = weightsm_array[min_nvalm];
644 clam_varm = varm_array[min_nvalm];
646 min_nvalp = FindMinimum(dwp_array, maxc);
647 omega_p1_0 = omegap_array[min_nvalp];
648 clam_weightsp = weightsp_array[min_nvalp];
649 clam_varp = varp_array[min_nvalp];
651 clam_osum = omega_m1_0 + omega_p1_0;
655 clam_minvar = 0.5*std::log(clam_osum);
660 lam_dg[fep_state-1] = clam_weightsm;
661 lam_variance[fep_state-1] = clam_varm;
664 if (fep_state < nlim-1)
666 lam_dg[fep_state] = clam_weightsp;
667 lam_variance[fep_state] = clam_varp;
670 if (expand->elamstats == elamstatsMINVAR)
672 bSufficientSamples = TRUE;
673 /* make sure they are all past a threshold */
674 for (i = 0; i < nlim; i++)
676 if (dfhist->n_at_lam[i] < expand->minvarmin)
678 bSufficientSamples = FALSE;
681 if (bSufficientSamples)
683 dfhist->sum_minvar[fep_state] = clam_minvar;
686 for (i = 0; i < nlim; i++)
688 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
690 expand->minvar_const = clam_minvar;
691 dfhist->sum_minvar[fep_state] = 0.0;
695 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
700 /* we need to rezero minvar now, since it could change at fep_state = 0 */
701 dfhist->sum_dg[0] = 0.0;
702 dfhist->sum_variance[0] = 0.0;
703 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
705 for (i = 1; i < nlim; i++)
707 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
708 dfhist->sum_variance[i] = std::sqrt(lam_variance[i-1] + gmx::square(dfhist->sum_variance[i-1]));
709 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
716 sfree(weightsm_array);
721 sfree(weightsp_array);
728 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
729 gmx_int64_t seed, gmx_int64_t step)
731 /* Choose new lambda value, and update transition matrix */
733 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
734 real r1, r2, de, trialprob, tprob = 0;
735 double *propose, *accept, *remainder;
738 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::ExpandedEnsemble); // We only draw once, so zero bits internal counter is fine
739 gmx::UniformRealDistribution<real> dist;
741 starting_fep_state = fep_state;
742 lamnew = fep_state; /* so that there is a default setting -- stays the same */
744 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
746 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
748 /* Use a marching method to run through the lambdas and get preliminary free energy data,
749 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
751 /* if we have enough at this lambda, move on to the next one */
753 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
755 lamnew = fep_state+1;
756 if (lamnew == nlim) /* whoops, stepped too far! */
771 snew(remainder, nlim);
773 for (i = 0; i < expand->lmc_repeats; i++)
775 rng.restart(step, i);
778 for (ifep = 0; ifep < nlim; ifep++)
784 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
786 /* use the Gibbs sampler, with restricted range */
787 if (expand->gibbsdeltalam < 0)
794 minfep = fep_state - expand->gibbsdeltalam;
795 maxfep = fep_state + expand->gibbsdeltalam;
806 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
808 if (expand->elmcmove == elmcmoveGIBBS)
810 for (ifep = minfep; ifep <= maxfep; ifep++)
812 propose[ifep] = p_k[ifep];
817 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
819 if (r1 <= p_k[lamnew])
826 else if (expand->elmcmove == elmcmoveMETGIBBS)
829 /* Metropolized Gibbs sampling */
830 for (ifep = minfep; ifep <= maxfep; ifep++)
832 remainder[ifep] = 1 - p_k[ifep];
835 /* find the proposal probabilities */
837 if (remainder[fep_state] == 0)
839 /* only the current state has any probability */
840 /* we have to stay at the current state */
845 for (ifep = minfep; ifep <= maxfep; ifep++)
847 if (ifep != fep_state)
849 propose[ifep] = p_k[ifep]/remainder[fep_state];
858 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
860 pnorm = p_k[lamtrial]/remainder[fep_state];
861 if (lamtrial != fep_state)
871 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
873 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
874 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
875 if (trialprob < tprob)
890 /* now figure out the acceptance probability for each */
891 for (ifep = minfep; ifep <= maxfep; ifep++)
894 if (remainder[ifep] != 0)
896 trialprob = (remainder[fep_state])/(remainder[ifep]);
900 trialprob = 1.0; /* this state is the only choice! */
902 if (trialprob < tprob)
906 /* probability for fep_state=0, but that's fine, it's never proposed! */
907 accept[ifep] = tprob;
913 /* it's possible some rounding is failing */
914 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
916 /* numerical rounding error -- no state other than the original has weight */
921 /* probably not a numerical issue */
923 int nerror = 200+(maxfep-minfep+1)*60;
925 snew(errorstr, nerror);
926 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
927 of sum weights. Generated detailed info for failure */
928 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);
929 for (ifep = minfep; ifep <= maxfep; ifep++)
931 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
933 gmx_fatal(FARGS, errorstr);
937 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
939 /* use the metropolis sampler with trial +/- 1 */
945 lamtrial = fep_state;
949 lamtrial = fep_state-1;
954 if (fep_state == nlim-1)
956 lamtrial = fep_state;
960 lamtrial = fep_state+1;
964 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
965 if (expand->elmcmove == elmcmoveMETROPOLIS)
968 trialprob = std::exp(de);
969 if (trialprob < tprob)
973 propose[fep_state] = 0;
974 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
975 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
976 accept[lamtrial] = tprob;
979 else if (expand->elmcmove == elmcmoveBARKER)
981 tprob = 1.0/(1.0+std::exp(-de));
983 propose[fep_state] = (1-tprob);
984 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
985 accept[fep_state] = 1.0;
986 accept[lamtrial] = 1.0;
1000 for (ifep = 0; ifep < nlim; ifep++)
1002 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1003 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1008 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1017 /* print out the weights to the log, along with current state */
1018 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1019 int fep_state, int frequency, gmx_int64_t step)
1021 int nlim, i, ifep, jfep;
1022 real dw, dg, dv, Tprint;
1023 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1024 gmx_bool bSimTemp = FALSE;
1026 nlim = fep->n_lambda;
1027 if (simtemp != NULL)
1032 if (step % frequency == 0)
1034 fprintf(outfile, " MC-lambda information\n");
1035 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1037 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1039 fprintf(outfile, " N");
1040 for (i = 0; i < efptNR; i++)
1042 if (fep->separate_dvdl[i])
1044 fprintf(outfile, "%7s", print_names[i]);
1046 else if ((i == efptTEMPERATURE) && bSimTemp)
1048 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1051 fprintf(outfile, " Count ");
1052 if (expand->elamstats == elamstatsMINVAR)
1054 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1058 fprintf(outfile, "G(in kT) dG(in kT)\n");
1060 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 = std::sqrt(gmx::square(dfhist->sum_variance[ifep+1]) - gmx::square(dfhist->sum_variance[ifep]));
1074 fprintf(outfile, "%3d", (ifep+1));
1075 for (i = 0; i < efptNR; i++)
1077 if (fep->separate_dvdl[i])
1079 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1081 else if (i == efptTEMPERATURE && bSimTemp)
1083 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1086 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1088 if (expand->elamstats == elamstatsWL)
1090 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1094 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1097 else /* we have equilibrated weights */
1099 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1101 if (expand->elamstats == elamstatsMINVAR)
1103 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1107 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1109 if (ifep == fep_state)
1111 fprintf(outfile, " <<\n");
1115 fprintf(outfile, " \n");
1118 fprintf(outfile, "\n");
1120 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1122 fprintf(outfile, " Transition Matrix\n");
1123 for (ifep = 0; ifep < nlim; ifep++)
1125 fprintf(outfile, "%12d", (ifep+1));
1127 fprintf(outfile, "\n");
1128 for (ifep = 0; ifep < nlim; ifep++)
1130 for (jfep = 0; jfep < nlim; jfep++)
1132 if (dfhist->n_at_lam[ifep] > 0)
1134 if (expand->bSymmetrizedTMatrix)
1136 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1140 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1147 fprintf(outfile, "%12.8f", Tprint);
1149 fprintf(outfile, "%3d\n", (ifep+1));
1152 fprintf(outfile, " Empirical Transition Matrix\n");
1153 for (ifep = 0; ifep < nlim; ifep++)
1155 fprintf(outfile, "%12d", (ifep+1));
1157 fprintf(outfile, "\n");
1158 for (ifep = 0; ifep < nlim; ifep++)
1160 for (jfep = 0; jfep < nlim; jfep++)
1162 if (dfhist->n_at_lam[ifep] > 0)
1164 if (expand->bSymmetrizedTMatrix)
1166 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1170 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1177 fprintf(outfile, "%12.8f", Tprint);
1179 fprintf(outfile, "%3d\n", (ifep+1));
1185 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1186 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1188 rvec *v, t_mdatoms *mdatoms)
1189 /* Note that the state variable is only needed for simulated tempering, not
1190 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1192 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1194 int i, nlim, lamnew, totalsamples;
1195 real oneovert, maxscaled = 0, maxweighted = 0;
1198 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1200 expand = ir->expandedvals;
1201 simtemp = ir->simtempvals;
1202 nlim = ir->fepvals->n_lambda;
1204 snew(scaled_lamee, nlim);
1205 snew(weighted_lamee, nlim);
1206 snew(pfep_lamee, nlim);
1209 /* update the count at the current lambda*/
1210 dfhist->n_at_lam[fep_state]++;
1212 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1213 pressure controlled.*/
1216 where does this PV term go?
1217 for (i=0;i<nlim;i++)
1219 fep_lamee[i] += pVTerm;
1223 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1224 /* we don't need to include the pressure term, since the volume is the same between the two.
1225 is there some term we are neglecting, however? */
1227 if (ir->efep != efepNO)
1229 for (i = 0; i < nlim; i++)
1233 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1234 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1235 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1239 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1240 /* mc_temp is currently set to the system reft unless otherwise defined */
1243 /* save these energies for printing, so they don't get overwritten by the next step */
1244 /* they aren't overwritten in the non-free energy case, but we always print with these
1252 for (i = 0; i < nlim; i++)
1254 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1259 for (i = 0; i < nlim; i++)
1261 pfep_lamee[i] = scaled_lamee[i];
1263 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1266 maxscaled = scaled_lamee[i];
1267 maxweighted = weighted_lamee[i];
1271 if (scaled_lamee[i] > maxscaled)
1273 maxscaled = scaled_lamee[i];
1275 if (weighted_lamee[i] > maxweighted)
1277 maxweighted = weighted_lamee[i];
1282 for (i = 0; i < nlim; i++)
1284 scaled_lamee[i] -= maxscaled;
1285 weighted_lamee[i] -= maxweighted;
1288 /* update weights - we decide whether or not to actually do this inside */
1290 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1291 if (bDoneEquilibrating)
1295 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1299 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1300 ir->expandedvals->lmc_seed, step);
1301 /* if using simulated tempering, we need to adjust the temperatures */
1302 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1307 int nstart, nend, gt;
1309 snew(buf_ngtc, ir->opts.ngtc);
1311 for (i = 0; i < ir->opts.ngtc; i++)
1313 if (ir->opts.ref_t[i] > 0)
1315 told = ir->opts.ref_t[i];
1316 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1317 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1321 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1324 nend = mdatoms->homenr;
1325 for (n = nstart; n < nend; n++)
1330 gt = mdatoms->cTC[n];
1332 for (d = 0; d < DIM; d++)
1334 v[n][d] *= buf_ngtc[gt];
1338 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1340 /* we need to recalculate the masses if the temperature has changed */
1341 init_npt_masses(ir, state, MassQ, FALSE);
1342 for (i = 0; i < state->nnhpres; i++)
1344 for (j = 0; j < ir->opts.nhchainlength; j++)
1346 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1349 for (i = 0; i < ir->opts.ngtc; i++)
1351 for (j = 0; j < ir->opts.nhchainlength; j++)
1353 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1360 /* now check on the Wang-Landau updating critera */
1362 if (EWL(expand->elamstats))
1364 bSwitchtoOneOverT = FALSE;
1365 if (expand->bWLoneovert)
1368 for (i = 0; i < nlim; i++)
1370 totalsamples += dfhist->n_at_lam[i];
1372 oneovert = (1.0*nlim)/totalsamples;
1373 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1374 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1375 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1376 (dfhist->wl_delta < expand->init_wl_delta))
1378 bSwitchtoOneOverT = TRUE;
1381 if (bSwitchtoOneOverT)
1383 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1387 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1390 for (i = 0; i < nlim; i++)
1392 dfhist->wl_histo[i] = 0;
1394 dfhist->wl_delta *= expand->wl_scale;
1397 fprintf(log, "\nStep %d: weights are now:", (int)step);
1398 for (i = 0; i < nlim; i++)
1400 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1408 sfree(scaled_lamee);
1409 sfree(weighted_lamee);