2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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/mdtypes/state.h"
63 #include "gromacs/random/threefry.h"
64 #include "gromacs/random/uniformrealdistribution.h"
65 #include "gromacs/timing/wallcycle.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxmpi.h"
68 #include "gromacs/utility/smalloc.h"
70 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
73 dfhist->wl_delta = expand->init_wl_delta;
74 for (i = 0; i < nlim; i++)
76 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
77 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
81 /* Eventually should contain all the functions needed to initialize expanded ensemble
82 before the md loop starts */
83 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
87 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
91 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
99 /* find the maximum value */
100 for (i = minfep; i <= maxfep; i++)
107 /* find the denominator */
108 for (i = minfep; i <= maxfep; i++)
110 *pks += std::exp(ene[i]-maxene);
113 for (i = minfep; i <= maxfep; i++)
115 p_k[i] = std::exp(ene[i]-maxene) / *pks;
119 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
128 for (i = 0; i < nlim; i++)
132 /* add the delta, since we need to make sure it's greater than zero, and
133 we need a non-arbitrary number? */
134 nene[i] = ene[i] + std::log(nvals[i]+delta);
138 nene[i] = ene[i] + std::log(nvals[i]);
142 /* find the maximum value */
144 for (i = 0; i < nlim; i++)
146 if (nene[i] > maxene)
152 /* subtract off the maximum, avoiding overflow */
153 for (i = 0; i < nlim; i++)
158 /* find the denominator */
159 for (i = 0; i < nlim; i++)
161 *pks += std::exp(nene[i]);
165 for (i = 0; i < nlim; i++)
167 p_k[i] = std::exp(nene[i]) / *pks;
172 static int FindMinimum(real *min_metric, int N)
179 min_val = min_metric[0];
181 for (nval = 0; nval < N; nval++)
183 if (min_metric[nval] < min_val)
185 min_val = min_metric[nval];
192 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
200 for (i = 0; i < nhisto; i++)
207 /* no samples! is bad!*/
211 nmean /= (real)nhisto;
214 for (i = 0; i < nhisto; i++)
216 /* make sure that all points are in the ratio < x < 1/ratio range */
217 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
226 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
230 gmx_bool bDoneEquilibrating = TRUE;
233 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
234 if (expand->lmc_forced_nstart > 0)
236 for (i = 0; i < nlim; i++)
238 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
241 bDoneEquilibrating = FALSE;
248 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
249 bDoneEquilibrating = TRUE;
251 /* calculate the total number of samples */
252 switch (expand->elmceq)
255 /* We have not equilibrated, and won't, ever. */
256 bDoneEquilibrating = FALSE;
259 /* we have equilibrated -- we're done */
260 bDoneEquilibrating = TRUE;
263 /* first, check if we are equilibrating by steps, if we're still under */
264 if (step < expand->equil_steps)
266 bDoneEquilibrating = FALSE;
271 for (i = 0; i < nlim; i++)
273 totalsamples += dfhist->n_at_lam[i];
275 if (totalsamples < expand->equil_samples)
277 bDoneEquilibrating = FALSE;
281 for (i = 0; i < nlim; i++)
283 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
286 bDoneEquilibrating = FALSE;
292 if (EWL(expand->elamstats)) /* This check is in readir as well, but
295 if (dfhist->wl_delta > expand->equil_wl_delta)
297 bDoneEquilibrating = FALSE;
302 /* we can use the flatness as a judge of good weights, as long as
303 we're not doing minvar, or Wang-Landau.
304 But turn off for now until we figure out exactly how we do this.
307 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
309 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
310 floats for this histogram function. */
313 snew(modhisto, nlim);
314 for (i = 0; i < nlim; i++)
316 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
318 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
322 bDoneEquilibrating = FALSE;
327 bDoneEquilibrating = TRUE;
331 return bDoneEquilibrating;
334 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
335 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
337 gmx_bool bSufficientSamples;
339 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
340 real omega_m1_0, omega_p1_0, clam_osum;
341 real de, de_function;
342 real cnval, zero_sum_weights;
343 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
344 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
345 real *lam_variance, *lam_dg;
348 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;
350 /* if we have equilibrated the weights, exit now */
356 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
358 dfhist->bEquil = TRUE;
359 /* zero out the visited states so we know how many equilibrated states we have
361 for (i = 0; i < nlim; i++)
363 dfhist->n_at_lam[i] = 0;
368 /* If we reached this far, we have not equilibrated yet, keep on
369 going resetting the weights */
371 if (EWL(expand->elamstats))
373 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
375 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
376 dfhist->wl_histo[fep_state] += 1.0;
378 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
382 /* first increment count */
383 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
384 for (i = 0; i < nlim; i++)
386 dfhist->wl_histo[i] += (real)p_k[i];
389 /* then increment weights (uses count) */
391 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
393 for (i = 0; i < nlim; i++)
395 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
397 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
402 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
403 dfhist->sum_weights[i] -= log(di);
409 zero_sum_weights = dfhist->sum_weights[0];
410 for (i = 0; i < nlim; i++)
412 dfhist->sum_weights[i] -= zero_sum_weights;
416 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
419 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
420 maxc = 2*expand->c_range+1;
423 snew(lam_variance, nlim);
425 snew(omegap_array, maxc);
426 snew(weightsp_array, maxc);
427 snew(varp_array, maxc);
428 snew(dwp_array, maxc);
430 snew(omegam_array, maxc);
431 snew(weightsm_array, maxc);
432 snew(varm_array, maxc);
433 snew(dwm_array, maxc);
435 /* unpack the current lambdas -- we will only update 2 of these */
437 for (i = 0; i < nlim-1; i++)
438 { /* only through the second to last */
439 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
440 lam_variance[i] = gmx::square(dfhist->sum_variance[i+1]) - gmx::square(dfhist->sum_variance[i]);
443 /* accumulate running averages */
444 for (nval = 0; nval < maxc; nval++)
446 /* constants for later use */
447 cnval = (real)(nval-expand->c_range);
448 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
451 de = std::exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
452 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
454 de_function = 1.0/(1.0+de);
456 else if (expand->elamstats == elamstatsMETROPOLIS)
464 de_function = 1.0/de;
467 dfhist->accum_m[fep_state][nval] += de_function;
468 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
471 if (fep_state < nlim-1)
473 de = std::exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
474 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
476 de_function = 1.0/(1.0+de);
478 else if (expand->elamstats == elamstatsMETROPOLIS)
486 de_function = 1.0/de;
489 dfhist->accum_p[fep_state][nval] += de_function;
490 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
493 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
495 n0 = dfhist->n_at_lam[fep_state];
498 nm1 = dfhist->n_at_lam[fep_state-1];
504 if (fep_state < nlim-1)
506 np1 = dfhist->n_at_lam[fep_state+1];
513 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
514 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;
518 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
519 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
520 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
521 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
524 if ((fep_state > 0 ) && (nm1 > 0))
526 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
527 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
530 if ((fep_state < nlim-1) && (np1 > 0))
532 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
533 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
547 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
550 real omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
551 clam_weightsm = (std::log(chi_m1_0) - std::log(chi_p1_m1)) + cnval;
552 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
557 if (fep_state < nlim-1)
561 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
564 real omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
565 clam_weightsp = (std::log(chi_m1_p1) - std::log(chi_p1_0)) + cnval;
566 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
573 omegam_array[nval] = omega_m1_0;
577 omegam_array[nval] = 0;
579 weightsm_array[nval] = clam_weightsm;
580 varm_array[nval] = clam_varm;
583 dwm_array[nval] = fabs( (cnval + std::log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
587 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
592 omegap_array[nval] = omega_p1_0;
596 omegap_array[nval] = 0;
598 weightsp_array[nval] = clam_weightsp;
599 varp_array[nval] = clam_varp;
600 if ((np1 > 0) && (n0 > 0))
602 dwp_array[nval] = fabs( (cnval + std::log((1.0*np1)/n0)) - lam_dg[fep_state] );
606 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
611 /* find the C's closest to the old weights value */
613 min_nvalm = FindMinimum(dwm_array, maxc);
614 omega_m1_0 = omegam_array[min_nvalm];
615 clam_weightsm = weightsm_array[min_nvalm];
616 clam_varm = varm_array[min_nvalm];
618 min_nvalp = FindMinimum(dwp_array, maxc);
619 omega_p1_0 = omegap_array[min_nvalp];
620 clam_weightsp = weightsp_array[min_nvalp];
621 clam_varp = varp_array[min_nvalp];
623 clam_osum = omega_m1_0 + omega_p1_0;
627 clam_minvar = 0.5*std::log(clam_osum);
632 lam_dg[fep_state-1] = clam_weightsm;
633 lam_variance[fep_state-1] = clam_varm;
636 if (fep_state < nlim-1)
638 lam_dg[fep_state] = clam_weightsp;
639 lam_variance[fep_state] = clam_varp;
642 if (expand->elamstats == elamstatsMINVAR)
644 bSufficientSamples = TRUE;
645 /* make sure they are all past a threshold */
646 for (i = 0; i < nlim; i++)
648 if (dfhist->n_at_lam[i] < expand->minvarmin)
650 bSufficientSamples = FALSE;
653 if (bSufficientSamples)
655 dfhist->sum_minvar[fep_state] = clam_minvar;
658 for (i = 0; i < nlim; i++)
660 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
662 expand->minvar_const = clam_minvar;
663 dfhist->sum_minvar[fep_state] = 0.0;
667 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
672 /* we need to rezero minvar now, since it could change at fep_state = 0 */
673 dfhist->sum_dg[0] = 0.0;
674 dfhist->sum_variance[0] = 0.0;
675 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
677 for (i = 1; i < nlim; i++)
679 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
680 dfhist->sum_variance[i] = std::sqrt(lam_variance[i-1] + gmx::square(dfhist->sum_variance[i-1]));
681 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
688 sfree(weightsm_array);
693 sfree(weightsp_array);
700 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
701 gmx_int64_t seed, gmx_int64_t step)
703 /* Choose new lambda value, and update transition matrix */
705 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
706 real r1, r2, de, trialprob, tprob = 0;
707 double *propose, *accept, *remainder;
710 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::ExpandedEnsemble); // We only draw once, so zero bits internal counter is fine
711 gmx::UniformRealDistribution<real> dist;
713 starting_fep_state = fep_state;
714 lamnew = fep_state; /* so that there is a default setting -- stays the same */
716 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
718 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
720 /* Use a marching method to run through the lambdas and get preliminary free energy data,
721 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
723 /* if we have enough at this lambda, move on to the next one */
725 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
727 lamnew = fep_state+1;
728 if (lamnew == nlim) /* whoops, stepped too far! */
743 snew(remainder, nlim);
745 for (i = 0; i < expand->lmc_repeats; i++)
747 rng.restart(step, i);
750 for (ifep = 0; ifep < nlim; ifep++)
756 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
758 /* use the Gibbs sampler, with restricted range */
759 if (expand->gibbsdeltalam < 0)
766 minfep = fep_state - expand->gibbsdeltalam;
767 maxfep = fep_state + expand->gibbsdeltalam;
778 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
780 if (expand->elmcmove == elmcmoveGIBBS)
782 for (ifep = minfep; ifep <= maxfep; ifep++)
784 propose[ifep] = p_k[ifep];
789 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
791 if (r1 <= p_k[lamnew])
798 else if (expand->elmcmove == elmcmoveMETGIBBS)
801 /* Metropolized Gibbs sampling */
802 for (ifep = minfep; ifep <= maxfep; ifep++)
804 remainder[ifep] = 1 - p_k[ifep];
807 /* find the proposal probabilities */
809 if (remainder[fep_state] == 0)
811 /* only the current state has any probability */
812 /* we have to stay at the current state */
817 for (ifep = minfep; ifep <= maxfep; ifep++)
819 if (ifep != fep_state)
821 propose[ifep] = p_k[ifep]/remainder[fep_state];
830 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
832 pnorm = p_k[lamtrial]/remainder[fep_state];
833 if (lamtrial != fep_state)
843 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
845 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
846 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
847 if (trialprob < tprob)
862 /* now figure out the acceptance probability for each */
863 for (ifep = minfep; ifep <= maxfep; ifep++)
866 if (remainder[ifep] != 0)
868 trialprob = (remainder[fep_state])/(remainder[ifep]);
872 trialprob = 1.0; /* this state is the only choice! */
874 if (trialprob < tprob)
878 /* probability for fep_state=0, but that's fine, it's never proposed! */
879 accept[ifep] = tprob;
885 /* it's possible some rounding is failing */
886 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
888 /* numerical rounding error -- no state other than the original has weight */
893 /* probably not a numerical issue */
895 int nerror = 200+(maxfep-minfep+1)*60;
897 snew(errorstr, nerror);
898 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
899 of sum weights. Generated detailed info for failure */
900 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);
901 for (ifep = minfep; ifep <= maxfep; ifep++)
903 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
905 gmx_fatal(FARGS, errorstr);
909 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
911 /* use the metropolis sampler with trial +/- 1 */
917 lamtrial = fep_state;
921 lamtrial = fep_state-1;
926 if (fep_state == nlim-1)
928 lamtrial = fep_state;
932 lamtrial = fep_state+1;
936 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
937 if (expand->elmcmove == elmcmoveMETROPOLIS)
940 trialprob = std::exp(de);
941 if (trialprob < tprob)
945 propose[fep_state] = 0;
946 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
947 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
948 accept[lamtrial] = tprob;
951 else if (expand->elmcmove == elmcmoveBARKER)
953 tprob = 1.0/(1.0+std::exp(-de));
955 propose[fep_state] = (1-tprob);
956 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
957 accept[fep_state] = 1.0;
958 accept[lamtrial] = 1.0;
972 for (ifep = 0; ifep < nlim; ifep++)
974 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
975 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
980 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
989 /* print out the weights to the log, along with current state */
990 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
991 int fep_state, int frequency, gmx_int64_t step)
993 int nlim, i, ifep, jfep;
994 real dw, dg, dv, Tprint;
995 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
996 gmx_bool bSimTemp = FALSE;
998 nlim = fep->n_lambda;
999 if (simtemp != nullptr)
1004 if (step % frequency == 0)
1006 fprintf(outfile, " MC-lambda information\n");
1007 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1009 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1011 fprintf(outfile, " N");
1012 for (i = 0; i < efptNR; i++)
1014 if (fep->separate_dvdl[i])
1016 fprintf(outfile, "%7s", print_names[i]);
1018 else if ((i == efptTEMPERATURE) && bSimTemp)
1020 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1023 fprintf(outfile, " Count ");
1024 if (expand->elamstats == elamstatsMINVAR)
1026 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1030 fprintf(outfile, "G(in kT) dG(in kT)\n");
1032 for (ifep = 0; ifep < nlim; ifep++)
1042 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1043 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1044 dv = std::sqrt(gmx::square(dfhist->sum_variance[ifep+1]) - gmx::square(dfhist->sum_variance[ifep]));
1046 fprintf(outfile, "%3d", (ifep+1));
1047 for (i = 0; i < efptNR; i++)
1049 if (fep->separate_dvdl[i])
1051 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1053 else if (i == efptTEMPERATURE && bSimTemp)
1055 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1058 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1060 if (expand->elamstats == elamstatsWL)
1062 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1066 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1069 else /* we have equilibrated weights */
1071 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1073 if (expand->elamstats == elamstatsMINVAR)
1075 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1079 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1081 if (ifep == fep_state)
1083 fprintf(outfile, " <<\n");
1087 fprintf(outfile, " \n");
1090 fprintf(outfile, "\n");
1092 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1094 fprintf(outfile, " Transition Matrix\n");
1095 for (ifep = 0; ifep < nlim; ifep++)
1097 fprintf(outfile, "%12d", (ifep+1));
1099 fprintf(outfile, "\n");
1100 for (ifep = 0; ifep < nlim; ifep++)
1102 for (jfep = 0; jfep < nlim; jfep++)
1104 if (dfhist->n_at_lam[ifep] > 0)
1106 if (expand->bSymmetrizedTMatrix)
1108 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1112 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1119 fprintf(outfile, "%12.8f", Tprint);
1121 fprintf(outfile, "%3d\n", (ifep+1));
1124 fprintf(outfile, " Empirical Transition Matrix\n");
1125 for (ifep = 0; ifep < nlim; ifep++)
1127 fprintf(outfile, "%12d", (ifep+1));
1129 fprintf(outfile, "\n");
1130 for (ifep = 0; ifep < nlim; ifep++)
1132 for (jfep = 0; jfep < nlim; jfep++)
1134 if (dfhist->n_at_lam[ifep] > 0)
1136 if (expand->bSymmetrizedTMatrix)
1138 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1142 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1149 fprintf(outfile, "%12.8f", Tprint);
1151 fprintf(outfile, "%3d\n", (ifep+1));
1157 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1158 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1160 rvec *v, t_mdatoms *mdatoms)
1161 /* Note that the state variable is only needed for simulated tempering, not
1162 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1164 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1166 int i, nlim, lamnew, totalsamples;
1167 real oneovert, maxscaled = 0, maxweighted = 0;
1170 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1172 expand = ir->expandedvals;
1173 simtemp = ir->simtempvals;
1174 nlim = ir->fepvals->n_lambda;
1176 snew(scaled_lamee, nlim);
1177 snew(weighted_lamee, nlim);
1178 snew(pfep_lamee, nlim);
1181 /* update the count at the current lambda*/
1182 dfhist->n_at_lam[fep_state]++;
1184 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1185 pressure controlled.*/
1188 where does this PV term go?
1189 for (i=0;i<nlim;i++)
1191 fep_lamee[i] += pVTerm;
1195 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1196 /* we don't need to include the pressure term, since the volume is the same between the two.
1197 is there some term we are neglecting, however? */
1199 if (ir->efep != efepNO)
1201 for (i = 0; i < nlim; i++)
1205 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1206 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1207 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1211 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1212 /* mc_temp is currently set to the system reft unless otherwise defined */
1215 /* save these energies for printing, so they don't get overwritten by the next step */
1216 /* they aren't overwritten in the non-free energy case, but we always print with these
1224 for (i = 0; i < nlim; i++)
1226 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1231 for (i = 0; i < nlim; i++)
1233 pfep_lamee[i] = scaled_lamee[i];
1235 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1238 maxscaled = scaled_lamee[i];
1239 maxweighted = weighted_lamee[i];
1243 if (scaled_lamee[i] > maxscaled)
1245 maxscaled = scaled_lamee[i];
1247 if (weighted_lamee[i] > maxweighted)
1249 maxweighted = weighted_lamee[i];
1254 for (i = 0; i < nlim; i++)
1256 scaled_lamee[i] -= maxscaled;
1257 weighted_lamee[i] -= maxweighted;
1260 /* update weights - we decide whether or not to actually do this inside */
1262 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1263 if (bDoneEquilibrating)
1267 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1271 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1272 ir->expandedvals->lmc_seed, step);
1273 /* if using simulated tempering, we need to adjust the temperatures */
1274 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1279 int nstart, nend, gt;
1281 snew(buf_ngtc, ir->opts.ngtc);
1283 for (i = 0; i < ir->opts.ngtc; i++)
1285 if (ir->opts.ref_t[i] > 0)
1287 told = ir->opts.ref_t[i];
1288 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1289 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1293 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1296 nend = mdatoms->homenr;
1297 for (n = nstart; n < nend; n++)
1302 gt = mdatoms->cTC[n];
1304 for (d = 0; d < DIM; d++)
1306 v[n][d] *= buf_ngtc[gt];
1310 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1312 /* we need to recalculate the masses if the temperature has changed */
1313 init_npt_masses(ir, state, MassQ, FALSE);
1314 for (i = 0; i < state->nnhpres; i++)
1316 for (j = 0; j < ir->opts.nhchainlength; j++)
1318 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1321 for (i = 0; i < ir->opts.ngtc; i++)
1323 for (j = 0; j < ir->opts.nhchainlength; j++)
1325 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1332 /* now check on the Wang-Landau updating critera */
1334 if (EWL(expand->elamstats))
1336 bSwitchtoOneOverT = FALSE;
1337 if (expand->bWLoneovert)
1340 for (i = 0; i < nlim; i++)
1342 totalsamples += dfhist->n_at_lam[i];
1344 oneovert = (1.0*nlim)/totalsamples;
1345 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1346 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1347 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1348 (dfhist->wl_delta < expand->init_wl_delta))
1350 bSwitchtoOneOverT = TRUE;
1353 if (bSwitchtoOneOverT)
1355 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1359 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1362 for (i = 0; i < nlim; i++)
1364 dfhist->wl_histo[i] = 0;
1366 dfhist->wl_delta *= expand->wl_scale;
1369 fprintf(log, "\nStep %d: weights are now:", (int)step);
1370 for (i = 0; i < nlim; i++)
1372 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1380 sfree(scaled_lamee);
1381 sfree(weighted_lamee);