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.
42 #include "gromacs/utility/smalloc.h"
44 #include "gromacs/fileio/confio.h"
47 #include "chargegroup.h"
64 #include "gromacs/random/random.h"
68 #include "gromacs/fileio/confio.h"
69 #include "gromacs/fileio/gmxfio.h"
70 #include "gromacs/fileio/trnio.h"
71 #include "gromacs/fileio/xtcio.h"
72 #include "gromacs/timing/wallcycle.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/gmxmpi.h"
76 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
79 dfhist->wl_delta = expand->init_wl_delta;
80 for (i = 0; i < nlim; i++)
82 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
83 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
87 /* Eventually should contain all the functions needed to initialize expanded ensemble
88 before the md loop starts */
89 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
93 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
97 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
104 maxene = ene[minfep];
105 /* find the maximum value */
106 for (i = minfep; i <= maxfep; i++)
113 /* find the denominator */
114 for (i = minfep; i <= maxfep; i++)
116 *pks += exp(ene[i]-maxene);
119 for (i = minfep; i <= maxfep; i++)
121 p_k[i] = exp(ene[i]-maxene) / *pks;
125 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
134 for (i = 0; i < nlim; i++)
138 /* add the delta, since we need to make sure it's greater than zero, and
139 we need a non-arbitrary number? */
140 nene[i] = ene[i] + log(nvals[i]+delta);
144 nene[i] = ene[i] + log(nvals[i]);
148 /* find the maximum value */
150 for (i = 0; i < nlim; i++)
152 if (nene[i] > maxene)
158 /* subtract off the maximum, avoiding overflow */
159 for (i = 0; i < nlim; i++)
164 /* find the denominator */
165 for (i = 0; i < nlim; i++)
167 *pks += exp(nene[i]);
171 for (i = 0; i < nlim; i++)
173 p_k[i] = exp(nene[i]) / *pks;
178 real do_logsum(int N, real *a_n)
182 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
187 /* compute maximum argument to exp(.) */
190 for (i = 1; i < N; i++)
192 maxarg = max(maxarg, a_n[i]);
195 /* compute sum of exp(a_n - maxarg) */
197 for (i = 0; i < N; i++)
199 sum = sum + exp(a_n[i] - maxarg);
202 /* compute log sum */
203 logsum = log(sum) + maxarg;
207 int FindMinimum(real *min_metric, int N)
214 min_val = min_metric[0];
216 for (nval = 0; nval < N; nval++)
218 if (min_metric[nval] < min_val)
220 min_val = min_metric[nval];
227 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
235 for (i = 0; i < nhisto; i++)
242 /* no samples! is bad!*/
246 nmean /= (real)nhisto;
249 for (i = 0; i < nhisto; i++)
251 /* make sure that all points are in the ratio < x < 1/ratio range */
252 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
261 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
265 gmx_bool bDoneEquilibrating = TRUE;
268 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
270 /* calculate the total number of samples */
271 switch (expand->elmceq)
274 /* We have not equilibrated, and won't, ever. */
277 /* we have equilibrated -- we're done */
280 /* first, check if we are equilibrating by steps, if we're still under */
281 if (step < expand->equil_steps)
283 bDoneEquilibrating = FALSE;
288 for (i = 0; i < nlim; i++)
290 totalsamples += dfhist->n_at_lam[i];
292 if (totalsamples < expand->equil_samples)
294 bDoneEquilibrating = FALSE;
298 for (i = 0; i < nlim; i++)
300 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
303 bDoneEquilibrating = FALSE;
309 if (EWL(expand->elamstats)) /* This check is in readir as well, but
312 if (dfhist->wl_delta > expand->equil_wl_delta)
314 bDoneEquilibrating = FALSE;
319 /* we can use the flatness as a judge of good weights, as long as
320 we're not doing minvar, or Wang-Landau.
321 But turn off for now until we figure out exactly how we do this.
324 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
326 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
327 floats for this histogram function. */
330 snew(modhisto, nlim);
331 for (i = 0; i < nlim; i++)
333 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
335 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
339 bDoneEquilibrating = FALSE;
343 bDoneEquilibrating = TRUE;
345 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
347 if (expand->lmc_forced_nstart > 0)
349 for (i = 0; i < nlim; i++)
351 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
354 bDoneEquilibrating = FALSE;
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 real maxdiff = 0.000000001;
366 gmx_bool bSufficientSamples;
367 int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
368 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
369 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
370 real de, de_function, dr, denom, maxdr;
371 real min_val, cnval, zero_sum_weights;
372 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
373 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
374 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
377 real *numweighted_lamee, *logfrac;
379 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;
381 /* if we have equilibrated the weights, exit now */
387 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
389 dfhist->bEquil = TRUE;
390 /* zero out the visited states so we know how many equilibrated states we have
392 for (i = 0; i < nlim; i++)
394 dfhist->n_at_lam[i] = 0;
399 /* If we reached this far, we have not equilibrated yet, keep on
400 going resetting the weights */
402 if (EWL(expand->elamstats))
404 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
406 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
407 dfhist->wl_histo[fep_state] += 1.0;
409 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
413 /* first increment count */
414 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
415 for (i = 0; i < nlim; i++)
417 dfhist->wl_histo[i] += (real)p_k[i];
420 /* then increment weights (uses count) */
422 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
424 for (i = 0; i < nlim; i++)
426 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
428 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
433 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
434 dfhist->sum_weights[i] -= log(di);
440 zero_sum_weights = dfhist->sum_weights[0];
441 for (i = 0; i < nlim; i++)
443 dfhist->sum_weights[i] -= zero_sum_weights;
447 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
450 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
451 maxc = 2*expand->c_range+1;
454 snew(lam_variance, nlim);
456 snew(omegap_array, maxc);
457 snew(weightsp_array, maxc);
458 snew(varp_array, maxc);
459 snew(dwp_array, maxc);
461 snew(omegam_array, maxc);
462 snew(weightsm_array, maxc);
463 snew(varm_array, maxc);
464 snew(dwm_array, maxc);
466 /* unpack the current lambdas -- we will only update 2 of these */
468 for (i = 0; i < nlim-1; i++)
469 { /* only through the second to last */
470 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
471 lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
474 /* accumulate running averages */
475 for (nval = 0; nval < maxc; nval++)
477 /* constants for later use */
478 cnval = (real)(nval-expand->c_range);
479 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
482 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
483 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
485 de_function = 1.0/(1.0+de);
487 else if (expand->elamstats == elamstatsMETROPOLIS)
495 de_function = 1.0/de;
498 dfhist->accum_m[fep_state][nval] += de_function;
499 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
502 if (fep_state < nlim-1)
504 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
505 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
507 de_function = 1.0/(1.0+de);
509 else if (expand->elamstats == elamstatsMETROPOLIS)
517 de_function = 1.0/de;
520 dfhist->accum_p[fep_state][nval] += de_function;
521 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
524 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
526 n0 = dfhist->n_at_lam[fep_state];
529 nm1 = dfhist->n_at_lam[fep_state-1];
535 if (fep_state < nlim-1)
537 np1 = dfhist->n_at_lam[fep_state+1];
544 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
545 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;
549 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
550 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
551 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
552 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
555 if ((fep_state > 0 ) && (nm1 > 0))
557 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
558 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
561 if ((fep_state < nlim-1) && (np1 > 0))
563 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
564 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
578 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
582 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
584 if ((n0 > 0) && (nm1 > 0))
586 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
587 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
591 if (fep_state < nlim-1)
595 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
599 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
601 if ((n0 > 0) && (np1 > 0))
603 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
604 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
610 omegam_array[nval] = omega_m1_0;
614 omegam_array[nval] = 0;
616 weightsm_array[nval] = clam_weightsm;
617 varm_array[nval] = clam_varm;
620 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
624 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
629 omegap_array[nval] = omega_p1_0;
633 omegap_array[nval] = 0;
635 weightsp_array[nval] = clam_weightsp;
636 varp_array[nval] = clam_varp;
637 if ((np1 > 0) && (n0 > 0))
639 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
643 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
648 /* find the C's closest to the old weights value */
650 min_nvalm = FindMinimum(dwm_array, maxc);
651 omega_m1_0 = omegam_array[min_nvalm];
652 clam_weightsm = weightsm_array[min_nvalm];
653 clam_varm = varm_array[min_nvalm];
655 min_nvalp = FindMinimum(dwp_array, maxc);
656 omega_p1_0 = omegap_array[min_nvalp];
657 clam_weightsp = weightsp_array[min_nvalp];
658 clam_varp = varp_array[min_nvalp];
660 clam_osum = omega_m1_0 + omega_p1_0;
664 clam_minvar = 0.5*log(clam_osum);
669 lam_dg[fep_state-1] = clam_weightsm;
670 lam_variance[fep_state-1] = clam_varm;
673 if (fep_state < nlim-1)
675 lam_dg[fep_state] = clam_weightsp;
676 lam_variance[fep_state] = clam_varp;
679 if (expand->elamstats == elamstatsMINVAR)
681 bSufficientSamples = TRUE;
682 /* make sure they are all past a threshold */
683 for (i = 0; i < nlim; i++)
685 if (dfhist->n_at_lam[i] < expand->minvarmin)
687 bSufficientSamples = FALSE;
690 if (bSufficientSamples)
692 dfhist->sum_minvar[fep_state] = clam_minvar;
695 for (i = 0; i < nlim; i++)
697 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
699 expand->minvar_const = clam_minvar;
700 dfhist->sum_minvar[fep_state] = 0.0;
704 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
709 /* we need to rezero minvar now, since it could change at fep_state = 0 */
710 dfhist->sum_dg[0] = 0.0;
711 dfhist->sum_variance[0] = 0.0;
712 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
714 for (i = 1; i < nlim; i++)
716 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
717 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
718 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
725 sfree(weightsm_array);
730 sfree(weightsp_array);
737 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
738 gmx_int64_t seed, gmx_int64_t step)
740 /* Choose new lambda value, and update transition matrix */
742 int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
743 real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
745 double *propose, *accept, *remainder;
748 gmx_bool bRestricted;
750 starting_fep_state = fep_state;
751 lamnew = fep_state; /* so that there is a default setting -- stays the same */
753 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
755 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
757 /* Use a marching method to run through the lambdas and get preliminary free energy data,
758 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
760 /* if we have enough at this lambda, move on to the next one */
762 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
764 lamnew = fep_state+1;
765 if (lamnew == nlim) /* whoops, stepped too far! */
780 snew(remainder, nlim);
782 for (i = 0; i < expand->lmc_repeats; i++)
786 gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
788 for (ifep = 0; ifep < nlim; ifep++)
794 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
797 /* use the Gibbs sampler, with restricted range */
798 if (expand->gibbsdeltalam < 0)
806 minfep = fep_state - expand->gibbsdeltalam;
807 maxfep = fep_state + expand->gibbsdeltalam;
818 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
820 if (expand->elmcmove == elmcmoveGIBBS)
822 for (ifep = minfep; ifep <= maxfep; ifep++)
824 propose[ifep] = p_k[ifep];
829 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
831 if (r1 <= p_k[lamnew])
838 else if (expand->elmcmove == elmcmoveMETGIBBS)
841 /* Metropolized Gibbs sampling */
842 for (ifep = minfep; ifep <= maxfep; ifep++)
844 remainder[ifep] = 1 - p_k[ifep];
847 /* find the proposal probabilities */
849 if (remainder[fep_state] == 0)
851 /* only the current state has any probability */
852 /* we have to stay at the current state */
857 for (ifep = minfep; ifep <= maxfep; ifep++)
859 if (ifep != fep_state)
861 propose[ifep] = p_k[ifep]/remainder[fep_state];
870 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
872 pnorm = p_k[lamtrial]/remainder[fep_state];
873 if (lamtrial != fep_state)
883 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
885 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
886 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
887 if (trialprob < tprob)
902 /* now figure out the acceptance probability for each */
903 for (ifep = minfep; ifep <= maxfep; ifep++)
906 if (remainder[ifep] != 0)
908 trialprob = (remainder[fep_state])/(remainder[ifep]);
912 trialprob = 1.0; /* this state is the only choice! */
914 if (trialprob < tprob)
918 /* probability for fep_state=0, but that's fine, it's never proposed! */
919 accept[ifep] = tprob;
925 /* it's possible some rounding is failing */
926 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
928 /* numerical rounding error -- no state other than the original has weight */
933 /* probably not a numerical issue */
935 int nerror = 200+(maxfep-minfep+1)*60;
937 snew(errorstr, nerror);
938 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
939 of sum weights. Generated detailed info for failure */
940 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);
941 for (ifep = minfep; ifep <= maxfep; ifep++)
943 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
945 gmx_fatal(FARGS, errorstr);
949 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
951 /* use the metropolis sampler with trial +/- 1 */
957 lamtrial = fep_state;
961 lamtrial = fep_state-1;
966 if (fep_state == nlim-1)
968 lamtrial = fep_state;
972 lamtrial = fep_state+1;
976 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
977 if (expand->elmcmove == elmcmoveMETROPOLIS)
981 if (trialprob < tprob)
985 propose[fep_state] = 0;
986 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
987 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
988 accept[lamtrial] = tprob;
991 else if (expand->elmcmove == elmcmoveBARKER)
993 tprob = 1.0/(1.0+exp(-de));
995 propose[fep_state] = (1-tprob);
996 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
997 accept[fep_state] = 1.0;
998 accept[lamtrial] = 1.0;
1012 for (ifep = 0; ifep < nlim; ifep++)
1014 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1015 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1020 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1029 /* print out the weights to the log, along with current state */
1030 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1031 int fep_state, int frequency, gmx_int64_t step)
1033 int nlim, i, ifep, jfep;
1034 real dw, dg, dv, dm, Tprint;
1036 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1037 gmx_bool bSimTemp = FALSE;
1039 nlim = fep->n_lambda;
1040 if (simtemp != NULL)
1045 if (mod(step, frequency) == 0)
1047 fprintf(outfile, " MC-lambda information\n");
1048 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1050 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1052 fprintf(outfile, " N");
1053 for (i = 0; i < efptNR; i++)
1055 if (fep->separate_dvdl[i])
1057 fprintf(outfile, "%7s", print_names[i]);
1059 else if ((i == efptTEMPERATURE) && bSimTemp)
1061 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1064 fprintf(outfile, " Count ");
1065 if (expand->elamstats == elamstatsMINVAR)
1067 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1071 fprintf(outfile, "G(in kT) dG(in kT)\n");
1073 for (ifep = 0; ifep < nlim; ifep++)
1084 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1085 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1086 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1087 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1090 fprintf(outfile, "%3d", (ifep+1));
1091 for (i = 0; i < efptNR; i++)
1093 if (fep->separate_dvdl[i])
1095 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1097 else if (i == efptTEMPERATURE && bSimTemp)
1099 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1102 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1104 if (expand->elamstats == elamstatsWL)
1106 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1110 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1113 else /* we have equilibrated weights */
1115 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1117 if (expand->elamstats == elamstatsMINVAR)
1119 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1123 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1125 if (ifep == fep_state)
1127 fprintf(outfile, " <<\n");
1131 fprintf(outfile, " \n");
1134 fprintf(outfile, "\n");
1136 if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1138 fprintf(outfile, " Transition Matrix\n");
1139 for (ifep = 0; ifep < nlim; ifep++)
1141 fprintf(outfile, "%12d", (ifep+1));
1143 fprintf(outfile, "\n");
1144 for (ifep = 0; ifep < nlim; ifep++)
1146 for (jfep = 0; jfep < nlim; jfep++)
1148 if (dfhist->n_at_lam[ifep] > 0)
1150 if (expand->bSymmetrizedTMatrix)
1152 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1156 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1163 fprintf(outfile, "%12.8f", Tprint);
1165 fprintf(outfile, "%3d\n", (ifep+1));
1168 fprintf(outfile, " Empirical Transition Matrix\n");
1169 for (ifep = 0; ifep < nlim; ifep++)
1171 fprintf(outfile, "%12d", (ifep+1));
1173 fprintf(outfile, "\n");
1174 for (ifep = 0; ifep < nlim; ifep++)
1176 for (jfep = 0; jfep < nlim; jfep++)
1178 if (dfhist->n_at_lam[ifep] > 0)
1180 if (expand->bSymmetrizedTMatrix)
1182 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1186 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1193 fprintf(outfile, "%12.8f", Tprint);
1195 fprintf(outfile, "%3d\n", (ifep+1));
1201 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1202 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1204 rvec *v, t_mdatoms *mdatoms)
1205 /* Note that the state variable is only needed for simulated tempering, not
1206 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1208 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1210 int i, nlim, lamnew, totalsamples;
1211 real oneovert, maxscaled = 0, maxweighted = 0;
1214 double *temperature_lambdas;
1215 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1217 expand = ir->expandedvals;
1218 simtemp = ir->simtempvals;
1219 nlim = ir->fepvals->n_lambda;
1221 snew(scaled_lamee, nlim);
1222 snew(weighted_lamee, nlim);
1223 snew(pfep_lamee, nlim);
1226 /* update the count at the current lambda*/
1227 dfhist->n_at_lam[fep_state]++;
1229 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1230 pressure controlled.*/
1233 where does this PV term go?
1234 for (i=0;i<nlim;i++)
1236 fep_lamee[i] += pVTerm;
1240 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1241 /* we don't need to include the pressure term, since the volume is the same between the two.
1242 is there some term we are neglecting, however? */
1244 if (ir->efep != efepNO)
1246 for (i = 0; i < nlim; i++)
1250 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1251 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1252 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1256 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1257 /* mc_temp is currently set to the system reft unless otherwise defined */
1260 /* save these energies for printing, so they don't get overwritten by the next step */
1261 /* they aren't overwritten in the non-free energy case, but we always print with these
1269 for (i = 0; i < nlim; i++)
1271 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1276 for (i = 0; i < nlim; i++)
1278 pfep_lamee[i] = scaled_lamee[i];
1280 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1283 maxscaled = scaled_lamee[i];
1284 maxweighted = weighted_lamee[i];
1288 if (scaled_lamee[i] > maxscaled)
1290 maxscaled = scaled_lamee[i];
1292 if (weighted_lamee[i] > maxweighted)
1294 maxweighted = weighted_lamee[i];
1299 for (i = 0; i < nlim; i++)
1301 scaled_lamee[i] -= maxscaled;
1302 weighted_lamee[i] -= maxweighted;
1305 /* update weights - we decide whether or not to actually do this inside */
1307 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1308 if (bDoneEquilibrating)
1312 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1316 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1317 ir->expandedvals->lmc_seed, step);
1318 /* if using simulated tempering, we need to adjust the temperatures */
1319 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1324 int nstart, nend, gt;
1326 snew(buf_ngtc, ir->opts.ngtc);
1328 for (i = 0; i < ir->opts.ngtc; i++)
1330 if (ir->opts.ref_t[i] > 0)
1332 told = ir->opts.ref_t[i];
1333 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1334 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1338 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1341 nend = mdatoms->homenr;
1342 for (n = nstart; n < nend; n++)
1347 gt = mdatoms->cTC[n];
1349 for (d = 0; d < DIM; d++)
1351 v[n][d] *= buf_ngtc[gt];
1355 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1357 /* we need to recalculate the masses if the temperature has changed */
1358 init_npt_masses(ir, state, MassQ, FALSE);
1359 for (i = 0; i < state->nnhpres; i++)
1361 for (j = 0; j < ir->opts.nhchainlength; j++)
1363 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1366 for (i = 0; i < ir->opts.ngtc; i++)
1368 for (j = 0; j < ir->opts.nhchainlength; j++)
1370 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1377 /* now check on the Wang-Landau updating critera */
1379 if (EWL(expand->elamstats))
1381 bSwitchtoOneOverT = FALSE;
1382 if (expand->bWLoneovert)
1385 for (i = 0; i < nlim; i++)
1387 totalsamples += dfhist->n_at_lam[i];
1389 oneovert = (1.0*nlim)/totalsamples;
1390 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1391 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1392 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1393 (dfhist->wl_delta < expand->init_wl_delta))
1395 bSwitchtoOneOverT = TRUE;
1398 if (bSwitchtoOneOverT)
1400 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1404 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1407 for (i = 0; i < nlim; i++)
1409 dfhist->wl_histo[i] = 0;
1411 dfhist->wl_delta *= expand->wl_scale;
1414 fprintf(log, "\nStep %d: weights are now:", (int)step);
1415 for (i = 0; i < nlim; i++)
1417 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1425 sfree(scaled_lamee);
1426 sfree(weighted_lamee);