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/bondf.h"
54 #include "gromacs/legacyheaders/pme.h"
55 #include "gromacs/legacyheaders/disre.h"
56 #include "gromacs/legacyheaders/orires.h"
57 #include "gromacs/legacyheaders/network.h"
58 #include "gromacs/legacyheaders/calcmu.h"
59 #include "gromacs/legacyheaders/constr.h"
60 #include "gromacs/random/random.h"
61 #include "gromacs/legacyheaders/domdec.h"
62 #include "gromacs/legacyheaders/macros.h"
64 #include "gromacs/fileio/confio.h"
65 #include "gromacs/fileio/gmxfio.h"
66 #include "gromacs/fileio/trnio.h"
67 #include "gromacs/fileio/xtcio.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxmpi.h"
72 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
75 dfhist->wl_delta = expand->init_wl_delta;
76 for (i = 0; i < nlim; i++)
78 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
79 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
83 /* Eventually should contain all the functions needed to initialize expanded ensemble
84 before the md loop starts */
85 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
89 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
93 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
100 maxene = ene[minfep];
101 /* find the maximum value */
102 for (i = minfep; i <= maxfep; i++)
109 /* find the denominator */
110 for (i = minfep; i <= maxfep; i++)
112 *pks += exp(ene[i]-maxene);
115 for (i = minfep; i <= maxfep; i++)
117 p_k[i] = exp(ene[i]-maxene) / *pks;
121 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
130 for (i = 0; i < nlim; i++)
134 /* add the delta, since we need to make sure it's greater than zero, and
135 we need a non-arbitrary number? */
136 nene[i] = ene[i] + log(nvals[i]+delta);
140 nene[i] = ene[i] + log(nvals[i]);
144 /* find the maximum value */
146 for (i = 0; i < nlim; i++)
148 if (nene[i] > maxene)
154 /* subtract off the maximum, avoiding overflow */
155 for (i = 0; i < nlim; i++)
160 /* find the denominator */
161 for (i = 0; i < nlim; i++)
163 *pks += exp(nene[i]);
167 for (i = 0; i < nlim; i++)
169 p_k[i] = exp(nene[i]) / *pks;
174 real do_logsum(int N, real *a_n)
178 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
183 /* compute maximum argument to exp(.) */
186 for (i = 1; i < N; i++)
188 maxarg = max(maxarg, a_n[i]);
191 /* compute sum of exp(a_n - maxarg) */
193 for (i = 0; i < N; i++)
195 sum = sum + exp(a_n[i] - maxarg);
198 /* compute log sum */
199 logsum = log(sum) + maxarg;
203 int FindMinimum(real *min_metric, int N)
210 min_val = min_metric[0];
212 for (nval = 0; nval < N; nval++)
214 if (min_metric[nval] < min_val)
216 min_val = min_metric[nval];
223 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
231 for (i = 0; i < nhisto; i++)
238 /* no samples! is bad!*/
242 nmean /= (real)nhisto;
245 for (i = 0; i < nhisto; i++)
247 /* make sure that all points are in the ratio < x < 1/ratio range */
248 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
257 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
261 gmx_bool bDoneEquilibrating = TRUE;
264 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
266 /* calculate the total number of samples */
267 switch (expand->elmceq)
270 /* We have not equilibrated, and won't, ever. */
273 /* we have equilibrated -- we're done */
276 /* first, check if we are equilibrating by steps, if we're still under */
277 if (step < expand->equil_steps)
279 bDoneEquilibrating = FALSE;
284 for (i = 0; i < nlim; i++)
286 totalsamples += dfhist->n_at_lam[i];
288 if (totalsamples < expand->equil_samples)
290 bDoneEquilibrating = FALSE;
294 for (i = 0; i < nlim; i++)
296 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
299 bDoneEquilibrating = FALSE;
305 if (EWL(expand->elamstats)) /* This check is in readir as well, but
308 if (dfhist->wl_delta > expand->equil_wl_delta)
310 bDoneEquilibrating = FALSE;
315 /* we can use the flatness as a judge of good weights, as long as
316 we're not doing minvar, or Wang-Landau.
317 But turn off for now until we figure out exactly how we do this.
320 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
322 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
323 floats for this histogram function. */
326 snew(modhisto, nlim);
327 for (i = 0; i < nlim; i++)
329 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
331 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
335 bDoneEquilibrating = FALSE;
339 bDoneEquilibrating = TRUE;
341 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
343 if (expand->lmc_forced_nstart > 0)
345 for (i = 0; i < nlim; i++)
347 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
350 bDoneEquilibrating = FALSE;
355 return bDoneEquilibrating;
358 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
359 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
361 real maxdiff = 0.000000001;
362 gmx_bool bSufficientSamples;
363 int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
364 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
365 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
366 real de, de_function, dr, denom, maxdr;
367 real min_val, cnval, zero_sum_weights;
368 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
369 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
370 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
373 real *numweighted_lamee, *logfrac;
375 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;
377 /* if we have equilibrated the weights, exit now */
383 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
385 dfhist->bEquil = TRUE;
386 /* zero out the visited states so we know how many equilibrated states we have
388 for (i = 0; i < nlim; i++)
390 dfhist->n_at_lam[i] = 0;
395 /* If we reached this far, we have not equilibrated yet, keep on
396 going resetting the weights */
398 if (EWL(expand->elamstats))
400 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
402 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
403 dfhist->wl_histo[fep_state] += 1.0;
405 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
409 /* first increment count */
410 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
411 for (i = 0; i < nlim; i++)
413 dfhist->wl_histo[i] += (real)p_k[i];
416 /* then increment weights (uses count) */
418 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
420 for (i = 0; i < nlim; i++)
422 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
424 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
429 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
430 dfhist->sum_weights[i] -= log(di);
436 zero_sum_weights = dfhist->sum_weights[0];
437 for (i = 0; i < nlim; i++)
439 dfhist->sum_weights[i] -= zero_sum_weights;
443 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
446 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
447 maxc = 2*expand->c_range+1;
450 snew(lam_variance, nlim);
452 snew(omegap_array, maxc);
453 snew(weightsp_array, maxc);
454 snew(varp_array, maxc);
455 snew(dwp_array, maxc);
457 snew(omegam_array, maxc);
458 snew(weightsm_array, maxc);
459 snew(varm_array, maxc);
460 snew(dwm_array, maxc);
462 /* unpack the current lambdas -- we will only update 2 of these */
464 for (i = 0; i < nlim-1; i++)
465 { /* only through the second to last */
466 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
467 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
470 /* accumulate running averages */
471 for (nval = 0; nval < maxc; nval++)
473 /* constants for later use */
474 cnval = (real)(nval-expand->c_range);
475 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
478 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
479 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
481 de_function = 1.0/(1.0+de);
483 else if (expand->elamstats == elamstatsMETROPOLIS)
491 de_function = 1.0/de;
494 dfhist->accum_m[fep_state][nval] += de_function;
495 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
498 if (fep_state < nlim-1)
500 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
501 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
503 de_function = 1.0/(1.0+de);
505 else if (expand->elamstats == elamstatsMETROPOLIS)
513 de_function = 1.0/de;
516 dfhist->accum_p[fep_state][nval] += de_function;
517 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
520 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
522 n0 = dfhist->n_at_lam[fep_state];
525 nm1 = dfhist->n_at_lam[fep_state-1];
531 if (fep_state < nlim-1)
533 np1 = dfhist->n_at_lam[fep_state+1];
540 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
541 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;
545 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
546 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
547 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
548 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
551 if ((fep_state > 0 ) && (nm1 > 0))
553 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
554 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
557 if ((fep_state < nlim-1) && (np1 > 0))
559 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
560 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
574 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
578 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
580 if ((n0 > 0) && (nm1 > 0))
582 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
583 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
587 if (fep_state < nlim-1)
591 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
595 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
597 if ((n0 > 0) && (np1 > 0))
599 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
600 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
606 omegam_array[nval] = omega_m1_0;
610 omegam_array[nval] = 0;
612 weightsm_array[nval] = clam_weightsm;
613 varm_array[nval] = clam_varm;
616 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
620 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
625 omegap_array[nval] = omega_p1_0;
629 omegap_array[nval] = 0;
631 weightsp_array[nval] = clam_weightsp;
632 varp_array[nval] = clam_varp;
633 if ((np1 > 0) && (n0 > 0))
635 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
639 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
644 /* find the C's closest to the old weights value */
646 min_nvalm = FindMinimum(dwm_array, maxc);
647 omega_m1_0 = omegam_array[min_nvalm];
648 clam_weightsm = weightsm_array[min_nvalm];
649 clam_varm = varm_array[min_nvalm];
651 min_nvalp = FindMinimum(dwp_array, maxc);
652 omega_p1_0 = omegap_array[min_nvalp];
653 clam_weightsp = weightsp_array[min_nvalp];
654 clam_varp = varp_array[min_nvalp];
656 clam_osum = omega_m1_0 + omega_p1_0;
660 clam_minvar = 0.5*log(clam_osum);
665 lam_dg[fep_state-1] = clam_weightsm;
666 lam_variance[fep_state-1] = clam_varm;
669 if (fep_state < nlim-1)
671 lam_dg[fep_state] = clam_weightsp;
672 lam_variance[fep_state] = clam_varp;
675 if (expand->elamstats == elamstatsMINVAR)
677 bSufficientSamples = TRUE;
678 /* make sure they are all past a threshold */
679 for (i = 0; i < nlim; i++)
681 if (dfhist->n_at_lam[i] < expand->minvarmin)
683 bSufficientSamples = FALSE;
686 if (bSufficientSamples)
688 dfhist->sum_minvar[fep_state] = clam_minvar;
691 for (i = 0; i < nlim; i++)
693 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
695 expand->minvar_const = clam_minvar;
696 dfhist->sum_minvar[fep_state] = 0.0;
700 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
705 /* we need to rezero minvar now, since it could change at fep_state = 0 */
706 dfhist->sum_dg[0] = 0.0;
707 dfhist->sum_variance[0] = 0.0;
708 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
710 for (i = 1; i < nlim; i++)
712 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
713 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
714 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
721 sfree(weightsm_array);
726 sfree(weightsp_array);
733 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
734 gmx_int64_t seed, gmx_int64_t step)
736 /* Choose new lambda value, and update transition matrix */
738 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
739 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
741 double *propose, *accept, *remainder;
744 gmx_bool bRestricted;
746 starting_fep_state = fep_state;
747 lamnew = fep_state; /* so that there is a default setting -- stays the same */
749 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
751 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
753 /* Use a marching method to run through the lambdas and get preliminary free energy data,
754 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
756 /* if we have enough at this lambda, move on to the next one */
758 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
760 lamnew = fep_state+1;
761 if (lamnew == nlim) /* whoops, stepped too far! */
776 snew(remainder, nlim);
778 for (i = 0; i < expand->lmc_repeats; i++)
782 gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
784 for (ifep = 0; ifep < nlim; ifep++)
790 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
793 /* use the Gibbs sampler, with restricted range */
794 if (expand->gibbsdeltalam < 0)
802 minfep = fep_state - expand->gibbsdeltalam;
803 maxfep = fep_state + expand->gibbsdeltalam;
814 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
816 if (expand->elmcmove == elmcmoveGIBBS)
818 for (ifep = minfep; ifep <= maxfep; ifep++)
820 propose[ifep] = p_k[ifep];
825 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
827 if (r1 <= p_k[lamnew])
834 else if (expand->elmcmove == elmcmoveMETGIBBS)
837 /* Metropolized Gibbs sampling */
838 for (ifep = minfep; ifep <= maxfep; ifep++)
840 remainder[ifep] = 1 - p_k[ifep];
843 /* find the proposal probabilities */
845 if (remainder[fep_state] == 0)
847 /* only the current state has any probability */
848 /* we have to stay at the current state */
853 for (ifep = minfep; ifep <= maxfep; ifep++)
855 if (ifep != fep_state)
857 propose[ifep] = p_k[ifep]/remainder[fep_state];
866 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
868 pnorm = p_k[lamtrial]/remainder[fep_state];
869 if (lamtrial != fep_state)
879 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
881 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
882 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
883 if (trialprob < tprob)
898 /* now figure out the acceptance probability for each */
899 for (ifep = minfep; ifep <= maxfep; ifep++)
902 if (remainder[ifep] != 0)
904 trialprob = (remainder[fep_state])/(remainder[ifep]);
908 trialprob = 1.0; /* this state is the only choice! */
910 if (trialprob < tprob)
914 /* probability for fep_state=0, but that's fine, it's never proposed! */
915 accept[ifep] = tprob;
921 /* it's possible some rounding is failing */
922 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
924 /* numerical rounding error -- no state other than the original has weight */
929 /* probably not a numerical issue */
931 int nerror = 200+(maxfep-minfep+1)*60;
933 snew(errorstr, nerror);
934 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
935 of sum weights. Generated detailed info for failure */
936 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);
937 for (ifep = minfep; ifep <= maxfep; ifep++)
939 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
941 gmx_fatal(FARGS, errorstr);
945 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
947 /* use the metropolis sampler with trial +/- 1 */
953 lamtrial = fep_state;
957 lamtrial = fep_state-1;
962 if (fep_state == nlim-1)
964 lamtrial = fep_state;
968 lamtrial = fep_state+1;
972 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
973 if (expand->elmcmove == elmcmoveMETROPOLIS)
977 if (trialprob < tprob)
981 propose[fep_state] = 0;
982 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
983 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
984 accept[lamtrial] = tprob;
987 else if (expand->elmcmove == elmcmoveBARKER)
989 tprob = 1.0/(1.0+exp(-de));
991 propose[fep_state] = (1-tprob);
992 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
993 accept[fep_state] = 1.0;
994 accept[lamtrial] = 1.0;
1008 for (ifep = 0; ifep < nlim; ifep++)
1010 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1011 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1016 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1025 /* print out the weights to the log, along with current state */
1026 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1027 int fep_state, int frequency, gmx_int64_t step)
1029 int nlim, i, ifep, jfep;
1030 real dw, dg, dv, dm, Tprint;
1032 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1033 gmx_bool bSimTemp = FALSE;
1035 nlim = fep->n_lambda;
1036 if (simtemp != NULL)
1041 if (mod(step, frequency) == 0)
1043 fprintf(outfile, " MC-lambda information\n");
1044 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1046 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1048 fprintf(outfile, " N");
1049 for (i = 0; i < efptNR; i++)
1051 if (fep->separate_dvdl[i])
1053 fprintf(outfile, "%7s", print_names[i]);
1055 else if ((i == efptTEMPERATURE) && bSimTemp)
1057 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1060 fprintf(outfile, " Count ");
1061 if (expand->elamstats == elamstatsMINVAR)
1063 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1067 fprintf(outfile, "G(in kT) dG(in kT)\n");
1069 for (ifep = 0; ifep < nlim; ifep++)
1080 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1081 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1082 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1083 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1086 fprintf(outfile, "%3d", (ifep+1));
1087 for (i = 0; i < efptNR; i++)
1089 if (fep->separate_dvdl[i])
1091 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1093 else if (i == efptTEMPERATURE && bSimTemp)
1095 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1098 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1100 if (expand->elamstats == elamstatsWL)
1102 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1106 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1109 else /* we have equilibrated weights */
1111 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1113 if (expand->elamstats == elamstatsMINVAR)
1115 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1119 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1121 if (ifep == fep_state)
1123 fprintf(outfile, " <<\n");
1127 fprintf(outfile, " \n");
1130 fprintf(outfile, "\n");
1132 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1134 fprintf(outfile, " Transition Matrix\n");
1135 for (ifep = 0; ifep < nlim; ifep++)
1137 fprintf(outfile, "%12d", (ifep+1));
1139 fprintf(outfile, "\n");
1140 for (ifep = 0; ifep < nlim; ifep++)
1142 for (jfep = 0; jfep < nlim; jfep++)
1144 if (dfhist->n_at_lam[ifep] > 0)
1146 if (expand->bSymmetrizedTMatrix)
1148 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1152 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1159 fprintf(outfile, "%12.8f", Tprint);
1161 fprintf(outfile, "%3d\n", (ifep+1));
1164 fprintf(outfile, " Empirical Transition Matrix\n");
1165 for (ifep = 0; ifep < nlim; ifep++)
1167 fprintf(outfile, "%12d", (ifep+1));
1169 fprintf(outfile, "\n");
1170 for (ifep = 0; ifep < nlim; ifep++)
1172 for (jfep = 0; jfep < nlim; jfep++)
1174 if (dfhist->n_at_lam[ifep] > 0)
1176 if (expand->bSymmetrizedTMatrix)
1178 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1182 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1189 fprintf(outfile, "%12.8f", Tprint);
1191 fprintf(outfile, "%3d\n", (ifep+1));
1197 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1198 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1200 rvec *v, t_mdatoms *mdatoms)
1201 /* Note that the state variable is only needed for simulated tempering, not
1202 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1204 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1206 int i, nlim, lamnew, totalsamples;
1207 real oneovert, maxscaled = 0, maxweighted = 0;
1210 double *temperature_lambdas;
1211 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1213 expand = ir->expandedvals;
1214 simtemp = ir->simtempvals;
1215 nlim = ir->fepvals->n_lambda;
1217 snew(scaled_lamee, nlim);
1218 snew(weighted_lamee, nlim);
1219 snew(pfep_lamee, nlim);
1222 /* update the count at the current lambda*/
1223 dfhist->n_at_lam[fep_state]++;
1225 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1226 pressure controlled.*/
1229 where does this PV term go?
1230 for (i=0;i<nlim;i++)
1232 fep_lamee[i] += pVTerm;
1236 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1237 /* we don't need to include the pressure term, since the volume is the same between the two.
1238 is there some term we are neglecting, however? */
1240 if (ir->efep != efepNO)
1242 for (i = 0; i < nlim; i++)
1246 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1247 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1248 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1252 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1253 /* mc_temp is currently set to the system reft unless otherwise defined */
1256 /* save these energies for printing, so they don't get overwritten by the next step */
1257 /* they aren't overwritten in the non-free energy case, but we always print with these
1265 for (i = 0; i < nlim; i++)
1267 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1272 for (i = 0; i < nlim; i++)
1274 pfep_lamee[i] = scaled_lamee[i];
1276 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1279 maxscaled = scaled_lamee[i];
1280 maxweighted = weighted_lamee[i];
1284 if (scaled_lamee[i] > maxscaled)
1286 maxscaled = scaled_lamee[i];
1288 if (weighted_lamee[i] > maxweighted)
1290 maxweighted = weighted_lamee[i];
1295 for (i = 0; i < nlim; i++)
1297 scaled_lamee[i] -= maxscaled;
1298 weighted_lamee[i] -= maxweighted;
1301 /* update weights - we decide whether or not to actually do this inside */
1303 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1304 if (bDoneEquilibrating)
1308 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1312 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1313 ir->expandedvals->lmc_seed, step);
1314 /* if using simulated tempering, we need to adjust the temperatures */
1315 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1320 int nstart, nend, gt;
1322 snew(buf_ngtc, ir->opts.ngtc);
1324 for (i = 0; i < ir->opts.ngtc; i++)
1326 if (ir->opts.ref_t[i] > 0)
1328 told = ir->opts.ref_t[i];
1329 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1330 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1334 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1337 nend = mdatoms->homenr;
1338 for (n = nstart; n < nend; n++)
1343 gt = mdatoms->cTC[n];
1345 for (d = 0; d < DIM; d++)
1347 v[n][d] *= buf_ngtc[gt];
1351 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1353 /* we need to recalculate the masses if the temperature has changed */
1354 init_npt_masses(ir, state, MassQ, FALSE);
1355 for (i = 0; i < state->nnhpres; i++)
1357 for (j = 0; j < ir->opts.nhchainlength; j++)
1359 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1362 for (i = 0; i < ir->opts.ngtc; i++)
1364 for (j = 0; j < ir->opts.nhchainlength; j++)
1366 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1373 /* now check on the Wang-Landau updating critera */
1375 if (EWL(expand->elamstats))
1377 bSwitchtoOneOverT = FALSE;
1378 if (expand->bWLoneovert)
1381 for (i = 0; i < nlim; i++)
1383 totalsamples += dfhist->n_at_lam[i];
1385 oneovert = (1.0*nlim)/totalsamples;
1386 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1387 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1388 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1389 (dfhist->wl_delta < expand->init_wl_delta))
1391 bSwitchtoOneOverT = TRUE;
1394 if (bSwitchtoOneOverT)
1396 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1400 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1403 for (i = 0; i < nlim; i++)
1405 dfhist->wl_histo[i] = 0;
1407 dfhist->wl_delta *= expand->wl_scale;
1410 fprintf(log, "\nStep %d: weights are now:", (int)step);
1411 for (i = 0; i < nlim; i++)
1413 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1421 sfree(scaled_lamee);
1422 sfree(weighted_lamee);