2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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.
44 #include "gromacs/domdec/domdec.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/xtcio.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/update.h"
59 #include "gromacs/mdtypes/enerdata.h"
60 #include "gromacs/mdtypes/forcerec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/mdatom.h"
64 #include "gromacs/mdtypes/state.h"
65 #include "gromacs/random/threefry.h"
66 #include "gromacs/random/uniformrealdistribution.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxmpi.h"
70 #include "gromacs/utility/smalloc.h"
72 static void init_df_history_weights(df_history_t *dfhist, const 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 void init_expanded_ensemble(gmx_bool bStateFromCP, const t_inputrec *ir, df_history_t *dfhist)
89 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
93 static void GenerateGibbsProbabilities(const 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 += std::exp(ene[i]-maxene);
115 for (i = minfep; i <= maxfep; i++)
117 p_k[i] = std::exp(ene[i]-maxene) / *pks;
121 static void GenerateWeightedGibbsProbabilities(const 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] + std::log(nvals[i]+delta);
140 nene[i] = ene[i] + std::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 += std::exp(nene[i]);
167 for (i = 0; i < nlim; i++)
169 p_k[i] = std::exp(nene[i]) / *pks;
174 static int FindMinimum(const real *min_metric, int N)
181 min_val = min_metric[0];
183 for (nval = 0; nval < N; nval++)
185 if (min_metric[nval] < min_val)
187 min_val = min_metric[nval];
194 static gmx_bool CheckHistogramRatios(int nhisto, const real *histo, real ratio)
202 for (i = 0; i < nhisto; i++)
209 /* no samples! is bad!*/
213 nmean /= static_cast<real>(nhisto);
216 for (i = 0; i < nhisto; i++)
218 /* make sure that all points are in the ratio < x < 1/ratio range */
219 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
228 static gmx_bool CheckIfDoneEquilibrating(int nlim, const t_expanded *expand, const df_history_t *dfhist, int64_t step)
232 gmx_bool bDoneEquilibrating = TRUE;
235 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
236 if (expand->lmc_forced_nstart > 0)
238 for (i = 0; i < nlim; i++)
240 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
243 bDoneEquilibrating = FALSE;
250 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
251 bDoneEquilibrating = TRUE;
253 /* calculate the total number of samples */
254 switch (expand->elmceq)
257 /* We have not equilibrated, and won't, ever. */
258 bDoneEquilibrating = FALSE;
261 /* we have equilibrated -- we're done */
262 bDoneEquilibrating = TRUE;
265 /* first, check if we are equilibrating by steps, if we're still under */
266 if (step < expand->equil_steps)
268 bDoneEquilibrating = FALSE;
273 for (i = 0; i < nlim; i++)
275 totalsamples += dfhist->n_at_lam[i];
277 if (totalsamples < expand->equil_samples)
279 bDoneEquilibrating = FALSE;
283 for (i = 0; i < nlim; i++)
285 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
288 bDoneEquilibrating = FALSE;
294 if (EWL(expand->elamstats)) /* This check is in readir as well, but
297 if (dfhist->wl_delta > expand->equil_wl_delta)
299 bDoneEquilibrating = FALSE;
304 /* we can use the flatness as a judge of good weights, as long as
305 we're not doing minvar, or Wang-Landau.
306 But turn off for now until we figure out exactly how we do this.
309 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
311 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
312 floats for this histogram function. */
315 snew(modhisto, nlim);
316 for (i = 0; i < nlim; i++)
318 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
320 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
324 bDoneEquilibrating = FALSE;
329 bDoneEquilibrating = TRUE;
333 return bDoneEquilibrating;
336 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
337 int fep_state, const real *scaled_lamee, const real *weighted_lamee, int64_t step)
339 gmx_bool bSufficientSamples;
341 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
342 real omega_m1_0, omega_p1_0, clam_osum;
343 real de, de_function;
344 real cnval, zero_sum_weights;
345 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
346 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
347 real *lam_variance, *lam_dg;
350 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;
352 /* if we have equilibrated the weights, exit now */
358 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
360 dfhist->bEquil = TRUE;
361 /* zero out the visited states so we know how many equilibrated states we have
363 for (i = 0; i < nlim; i++)
365 dfhist->n_at_lam[i] = 0;
370 /* If we reached this far, we have not equilibrated yet, keep on
371 going resetting the weights */
373 if (EWL(expand->elamstats))
375 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
377 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
378 dfhist->wl_histo[fep_state] += 1.0;
380 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
384 /* first increment count */
385 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
386 for (i = 0; i < nlim; i++)
388 dfhist->wl_histo[i] += static_cast<real>(p_k[i]);
391 /* then increment weights (uses count) */
393 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
395 for (i = 0; i < nlim; i++)
397 dfhist->sum_weights[i] -= dfhist->wl_delta*static_cast<real>(p_k[i]);
399 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
404 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
405 dfhist->sum_weights[i] -= log(di);
411 zero_sum_weights = dfhist->sum_weights[0];
412 for (i = 0; i < nlim; i++)
414 dfhist->sum_weights[i] -= zero_sum_weights;
418 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
421 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
422 maxc = 2*expand->c_range+1;
425 snew(lam_variance, nlim);
427 snew(omegap_array, maxc);
428 snew(weightsp_array, maxc);
429 snew(varp_array, maxc);
430 snew(dwp_array, maxc);
432 snew(omegam_array, maxc);
433 snew(weightsm_array, maxc);
434 snew(varm_array, maxc);
435 snew(dwm_array, maxc);
437 /* unpack the current lambdas -- we will only update 2 of these */
439 for (i = 0; i < nlim-1; i++)
440 { /* only through the second to last */
441 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
442 lam_variance[i] = gmx::square(dfhist->sum_variance[i+1]) - gmx::square(dfhist->sum_variance[i]);
445 /* accumulate running averages */
446 for (nval = 0; nval < maxc; nval++)
448 /* constants for later use */
449 cnval = static_cast<real>(nval-expand->c_range);
450 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
453 de = std::exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
454 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
456 de_function = 1.0/(1.0+de);
458 else if (expand->elamstats == elamstatsMETROPOLIS)
466 de_function = 1.0/de;
469 dfhist->accum_m[fep_state][nval] += de_function;
470 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
473 if (fep_state < nlim-1)
475 de = std::exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
476 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
478 de_function = 1.0/(1.0+de);
480 else if (expand->elamstats == elamstatsMETROPOLIS)
488 de_function = 1.0/de;
491 dfhist->accum_p[fep_state][nval] += de_function;
492 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
495 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
497 n0 = dfhist->n_at_lam[fep_state];
500 nm1 = dfhist->n_at_lam[fep_state-1];
506 if (fep_state < nlim-1)
508 np1 = dfhist->n_at_lam[fep_state+1];
515 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
516 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;
520 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
521 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
522 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
523 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
526 if ((fep_state > 0 ) && (nm1 > 0))
528 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
529 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
532 if ((fep_state < nlim-1) && (np1 > 0))
534 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
535 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
549 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
552 real omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
553 clam_weightsm = (std::log(chi_m1_0) - std::log(chi_p1_m1)) + cnval;
554 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
559 if (fep_state < nlim-1)
563 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
566 real omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
567 clam_weightsp = (std::log(chi_m1_p1) - std::log(chi_p1_0)) + cnval;
568 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
575 omegam_array[nval] = omega_m1_0;
579 omegam_array[nval] = 0;
581 weightsm_array[nval] = clam_weightsm;
582 varm_array[nval] = clam_varm;
585 dwm_array[nval] = fabs( (cnval + std::log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
589 dwm_array[nval] = std::fabs( cnval - lam_dg[fep_state-1] );
594 omegap_array[nval] = omega_p1_0;
598 omegap_array[nval] = 0;
600 weightsp_array[nval] = clam_weightsp;
601 varp_array[nval] = clam_varp;
602 if ((np1 > 0) && (n0 > 0))
604 dwp_array[nval] = fabs( (cnval + std::log((1.0*np1)/n0)) - lam_dg[fep_state] );
608 dwp_array[nval] = std::fabs( cnval - lam_dg[fep_state] );
613 /* find the C's closest to the old weights value */
615 min_nvalm = FindMinimum(dwm_array, maxc);
616 omega_m1_0 = omegam_array[min_nvalm];
617 clam_weightsm = weightsm_array[min_nvalm];
618 clam_varm = varm_array[min_nvalm];
620 min_nvalp = FindMinimum(dwp_array, maxc);
621 omega_p1_0 = omegap_array[min_nvalp];
622 clam_weightsp = weightsp_array[min_nvalp];
623 clam_varp = varp_array[min_nvalp];
625 clam_osum = omega_m1_0 + omega_p1_0;
629 clam_minvar = 0.5*std::log(clam_osum);
634 lam_dg[fep_state-1] = clam_weightsm;
635 lam_variance[fep_state-1] = clam_varm;
638 if (fep_state < nlim-1)
640 lam_dg[fep_state] = clam_weightsp;
641 lam_variance[fep_state] = clam_varp;
644 if (expand->elamstats == elamstatsMINVAR)
646 bSufficientSamples = TRUE;
647 /* make sure they are all past a threshold */
648 for (i = 0; i < nlim; i++)
650 if (dfhist->n_at_lam[i] < expand->minvarmin)
652 bSufficientSamples = FALSE;
655 if (bSufficientSamples)
657 dfhist->sum_minvar[fep_state] = clam_minvar;
660 for (i = 0; i < nlim; i++)
662 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
664 expand->minvar_const = clam_minvar;
665 dfhist->sum_minvar[fep_state] = 0.0;
669 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
674 /* we need to rezero minvar now, since it could change at fep_state = 0 */
675 dfhist->sum_dg[0] = 0.0;
676 dfhist->sum_variance[0] = 0.0;
677 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
679 for (i = 1; i < nlim; i++)
681 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
682 dfhist->sum_variance[i] = std::sqrt(lam_variance[i-1] + gmx::square(dfhist->sum_variance[i-1]));
683 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
690 sfree(weightsm_array);
695 sfree(weightsp_array);
702 static int ChooseNewLambda(int nlim, const t_expanded *expand, df_history_t *dfhist, int fep_state,
703 const real *weighted_lamee, double *p_k,
704 int64_t seed, int64_t step)
706 /* Choose new lambda value, and update transition matrix */
708 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
709 real r1, r2, de, trialprob, tprob = 0;
710 double *propose, *accept, *remainder;
713 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::ExpandedEnsemble); // We only draw once, so zero bits internal counter is fine
714 gmx::UniformRealDistribution<real> dist;
716 starting_fep_state = fep_state;
717 lamnew = fep_state; /* so that there is a default setting -- stays the same */
719 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
721 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
723 /* Use a marching method to run through the lambdas and get preliminary free energy data,
724 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
726 /* if we have enough at this lambda, move on to the next one */
728 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
730 lamnew = fep_state+1;
731 if (lamnew == nlim) /* whoops, stepped too far! */
746 snew(remainder, nlim);
748 for (i = 0; i < expand->lmc_repeats; i++)
750 rng.restart(step, i);
753 for (ifep = 0; ifep < nlim; ifep++)
759 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
761 /* use the Gibbs sampler, with restricted range */
762 if (expand->gibbsdeltalam < 0)
769 minfep = fep_state - expand->gibbsdeltalam;
770 maxfep = fep_state + expand->gibbsdeltalam;
781 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
783 if (expand->elmcmove == elmcmoveGIBBS)
785 for (ifep = minfep; ifep <= maxfep; ifep++)
787 propose[ifep] = p_k[ifep];
792 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
794 if (r1 <= p_k[lamnew])
801 else if (expand->elmcmove == elmcmoveMETGIBBS)
804 /* Metropolized Gibbs sampling */
805 for (ifep = minfep; ifep <= maxfep; ifep++)
807 remainder[ifep] = 1 - p_k[ifep];
810 /* find the proposal probabilities */
812 if (remainder[fep_state] == 0)
814 /* only the current state has any probability */
815 /* we have to stay at the current state */
820 for (ifep = minfep; ifep <= maxfep; ifep++)
822 if (ifep != fep_state)
824 propose[ifep] = p_k[ifep]/remainder[fep_state];
833 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
835 pnorm = p_k[lamtrial]/remainder[fep_state];
836 if (lamtrial != fep_state)
846 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
848 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
849 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
850 if (trialprob < tprob)
865 /* now figure out the acceptance probability for each */
866 for (ifep = minfep; ifep <= maxfep; ifep++)
869 if (remainder[ifep] != 0)
871 trialprob = (remainder[fep_state])/(remainder[ifep]);
875 trialprob = 1.0; /* this state is the only choice! */
877 if (trialprob < tprob)
881 /* probability for fep_state=0, but that's fine, it's never proposed! */
882 accept[ifep] = tprob;
888 /* it's possible some rounding is failing */
889 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
891 /* numerical rounding error -- no state other than the original has weight */
896 /* probably not a numerical issue */
898 int nerror = 200+(maxfep-minfep+1)*60;
900 snew(errorstr, nerror);
901 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
902 of sum weights. Generated detailed info for failure */
903 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);
904 for (ifep = minfep; ifep <= maxfep; ifep++)
906 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
908 gmx_fatal(FARGS, "%s", errorstr);
912 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
914 /* use the metropolis sampler with trial +/- 1 */
920 lamtrial = fep_state;
924 lamtrial = fep_state-1;
929 if (fep_state == nlim-1)
931 lamtrial = fep_state;
935 lamtrial = fep_state+1;
939 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
940 if (expand->elmcmove == elmcmoveMETROPOLIS)
943 trialprob = std::exp(de);
944 if (trialprob < tprob)
948 propose[fep_state] = 0;
949 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
950 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
951 accept[lamtrial] = tprob;
954 else if (expand->elmcmove == elmcmoveBARKER)
956 tprob = 1.0/(1.0+std::exp(-de));
958 propose[fep_state] = (1-tprob);
959 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
960 accept[fep_state] = 1.0;
961 accept[lamtrial] = 1.0;
975 for (ifep = 0; ifep < nlim; ifep++)
977 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
978 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
983 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
992 /* print out the weights to the log, along with current state */
993 void PrintFreeEnergyInfoToFile(FILE *outfile, const t_lambda *fep, const t_expanded *expand,
994 const t_simtemp *simtemp, const df_history_t *dfhist,
995 int fep_state, int frequency, int64_t step)
997 int nlim, i, ifep, jfep;
998 real dw, dg, dv, Tprint;
999 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1000 gmx_bool bSimTemp = FALSE;
1002 nlim = fep->n_lambda;
1003 if (simtemp != nullptr)
1008 if (step % frequency == 0)
1010 fprintf(outfile, " MC-lambda information\n");
1011 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1013 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1015 fprintf(outfile, " N");
1016 for (i = 0; i < efptNR; i++)
1018 if (fep->separate_dvdl[i])
1020 fprintf(outfile, "%7s", print_names[i]);
1022 else if ((i == efptTEMPERATURE) && bSimTemp)
1024 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1027 fprintf(outfile, " Count ");
1028 if (expand->elamstats == elamstatsMINVAR)
1030 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1034 fprintf(outfile, "G(in kT) dG(in kT)\n");
1036 for (ifep = 0; ifep < nlim; ifep++)
1046 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1047 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1048 dv = std::sqrt(gmx::square(dfhist->sum_variance[ifep+1]) - gmx::square(dfhist->sum_variance[ifep]));
1050 fprintf(outfile, "%3d", (ifep+1));
1051 for (i = 0; i < efptNR; i++)
1053 if (fep->separate_dvdl[i])
1055 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1057 else if (i == efptTEMPERATURE && bSimTemp)
1059 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1062 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1064 if (expand->elamstats == elamstatsWL)
1066 fprintf(outfile, " %8d", static_cast<int>(dfhist->wl_histo[ifep]));
1070 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1073 else /* we have equilibrated weights */
1075 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1077 if (expand->elamstats == elamstatsMINVAR)
1079 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1083 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1085 if (ifep == fep_state)
1087 fprintf(outfile, " <<\n");
1091 fprintf(outfile, " \n");
1094 fprintf(outfile, "\n");
1096 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1098 fprintf(outfile, " Transition Matrix\n");
1099 for (ifep = 0; ifep < nlim; ifep++)
1101 fprintf(outfile, "%12d", (ifep+1));
1103 fprintf(outfile, "\n");
1104 for (ifep = 0; ifep < nlim; ifep++)
1106 for (jfep = 0; jfep < nlim; jfep++)
1108 if (dfhist->n_at_lam[ifep] > 0)
1110 if (expand->bSymmetrizedTMatrix)
1112 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1116 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1123 fprintf(outfile, "%12.8f", Tprint);
1125 fprintf(outfile, "%3d\n", (ifep+1));
1128 fprintf(outfile, " Empirical Transition Matrix\n");
1129 for (ifep = 0; ifep < nlim; ifep++)
1131 fprintf(outfile, "%12d", (ifep+1));
1133 fprintf(outfile, "\n");
1134 for (ifep = 0; ifep < nlim; ifep++)
1136 for (jfep = 0; jfep < nlim; jfep++)
1138 if (dfhist->n_at_lam[ifep] > 0)
1140 if (expand->bSymmetrizedTMatrix)
1142 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1146 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1153 fprintf(outfile, "%12.8f", Tprint);
1155 fprintf(outfile, "%3d\n", (ifep+1));
1161 int ExpandedEnsembleDynamics(FILE *log, const t_inputrec *ir, const gmx_enerdata_t *enerd,
1162 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1164 rvec *v, const t_mdatoms *mdatoms)
1165 /* Note that the state variable is only needed for simulated tempering, not
1166 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1168 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1170 int i, nlim, lamnew, totalsamples;
1171 real oneovert, maxscaled = 0, maxweighted = 0;
1174 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1176 expand = ir->expandedvals;
1177 simtemp = ir->simtempvals;
1178 nlim = ir->fepvals->n_lambda;
1180 snew(scaled_lamee, nlim);
1181 snew(weighted_lamee, nlim);
1182 snew(pfep_lamee, nlim);
1185 /* update the count at the current lambda*/
1186 dfhist->n_at_lam[fep_state]++;
1188 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1189 pressure controlled.*/
1192 where does this PV term go?
1193 for (i=0;i<nlim;i++)
1195 fep_lamee[i] += pVTerm;
1199 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1200 /* we don't need to include the pressure term, since the volume is the same between the two.
1201 is there some term we are neglecting, however? */
1203 if (ir->efep != efepNO)
1205 for (i = 0; i < nlim; i++)
1209 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1210 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1211 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1215 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1216 /* mc_temp is currently set to the system reft unless otherwise defined */
1219 /* save these energies for printing, so they don't get overwritten by the next step */
1220 /* they aren't overwritten in the non-free energy case, but we always print with these
1228 for (i = 0; i < nlim; i++)
1230 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1235 for (i = 0; i < nlim; i++)
1237 pfep_lamee[i] = scaled_lamee[i];
1239 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1242 maxscaled = scaled_lamee[i];
1243 maxweighted = weighted_lamee[i];
1247 if (scaled_lamee[i] > maxscaled)
1249 maxscaled = scaled_lamee[i];
1251 if (weighted_lamee[i] > maxweighted)
1253 maxweighted = weighted_lamee[i];
1258 for (i = 0; i < nlim; i++)
1260 scaled_lamee[i] -= maxscaled;
1261 weighted_lamee[i] -= maxweighted;
1264 /* update weights - we decide whether or not to actually do this inside */
1266 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1267 if (bDoneEquilibrating)
1271 fprintf(log, "\nStep %" PRId64 ": Weights have equilibrated, using criteria: %s\n", step, elmceq_names[expand->elmceq]);
1275 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1276 ir->expandedvals->lmc_seed, step);
1277 /* if using simulated tempering, we need to adjust the temperatures */
1278 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1283 int nstart, nend, gt;
1285 snew(buf_ngtc, ir->opts.ngtc);
1287 for (i = 0; i < ir->opts.ngtc; i++)
1289 if (ir->opts.ref_t[i] > 0)
1291 told = ir->opts.ref_t[i];
1292 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1293 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1297 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1300 nend = mdatoms->homenr;
1301 for (n = nstart; n < nend; n++)
1306 gt = mdatoms->cTC[n];
1308 for (d = 0; d < DIM; d++)
1310 v[n][d] *= buf_ngtc[gt];
1314 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1316 /* we need to recalculate the masses if the temperature has changed */
1317 init_npt_masses(ir, state, MassQ, FALSE);
1318 for (i = 0; i < state->nnhpres; i++)
1320 for (j = 0; j < ir->opts.nhchainlength; j++)
1322 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1325 for (i = 0; i < ir->opts.ngtc; i++)
1327 for (j = 0; j < ir->opts.nhchainlength; j++)
1329 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1336 /* now check on the Wang-Landau updating critera */
1338 if (EWL(expand->elamstats))
1340 bSwitchtoOneOverT = FALSE;
1341 if (expand->bWLoneovert)
1344 for (i = 0; i < nlim; i++)
1346 totalsamples += dfhist->n_at_lam[i];
1348 oneovert = (1.0*nlim)/totalsamples;
1349 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1350 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1351 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1352 (dfhist->wl_delta < expand->init_wl_delta))
1354 bSwitchtoOneOverT = TRUE;
1357 if (bSwitchtoOneOverT)
1359 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1363 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1366 for (i = 0; i < nlim; i++)
1368 dfhist->wl_histo[i] = 0;
1370 dfhist->wl_delta *= expand->wl_scale;
1373 fprintf(log, "\nStep %d: weights are now:", static_cast<int>(step));
1374 for (i = 0; i < nlim; i++)
1376 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1384 sfree(scaled_lamee);
1385 sfree(weighted_lamee);