2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012-2018, The GROMACS development team.
5 * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/math/units.h"
44 #include "gromacs/math/utilities.h"
45 #include "gromacs/mdtypes/enerdata.h"
46 #include "gromacs/mdtypes/forcerec.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/mdtypes/state.h"
50 #include "gromacs/random/threefry.h"
51 #include "gromacs/random/uniformrealdistribution.h"
52 #include "gromacs/utility/enumerationhelpers.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 #include "expanded_internal.h"
58 static void init_df_history_weights(df_history_t* dfhist, const t_expanded* expand, int nlim)
61 dfhist->wl_delta = expand->init_wl_delta;
62 for (i = 0; i < nlim; i++)
64 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
65 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
69 /* Eventually should contain all the functions needed to initialize expanded ensemble
70 before the md loop starts */
71 void init_expanded_ensemble(gmx_bool bStateFromCP, const t_inputrec* ir, df_history_t* dfhist)
75 init_df_history_weights(dfhist, ir->expandedvals.get(), ir->fepvals->n_lambda);
79 static void GenerateGibbsProbabilities(const real* ene, double* p_k, double* pks, int minfep, int maxfep)
87 /* find the maximum value */
88 for (i = minfep; i <= maxfep; i++)
95 /* find the denominator */
96 for (i = minfep; i <= maxfep; i++)
98 *pks += std::exp(ene[i] - maxene);
101 for (i = minfep; i <= maxfep; i++)
103 p_k[i] = std::exp(ene[i] - maxene) / *pks;
108 GenerateWeightedGibbsProbabilities(const real* ene, double* p_k, double* pks, int nlim, real* nvals, real delta)
117 for (i = 0; i < nlim; i++)
121 /* add the delta, since we need to make sure it's greater than zero, and
122 we need a non-arbitrary number? */
123 nene[i] = ene[i] + std::log(nvals[i] + delta);
127 nene[i] = ene[i] + std::log(nvals[i]);
131 /* find the maximum value */
133 for (i = 0; i < nlim; i++)
135 if (nene[i] > maxene)
141 /* subtract off the maximum, avoiding overflow */
142 for (i = 0; i < nlim; i++)
147 /* find the denominator */
148 for (i = 0; i < nlim; i++)
150 *pks += std::exp(nene[i]);
154 for (i = 0; i < nlim; i++)
156 p_k[i] = std::exp(nene[i]) / *pks;
161 static int FindMinimum(const real* min_metric, int N)
168 min_val = min_metric[0];
170 for (nval = 0; nval < N; nval++)
172 if (min_metric[nval] < min_val)
174 min_val = min_metric[nval];
181 static gmx_bool CheckHistogramRatios(int nhisto, const real* histo, real ratio)
189 for (i = 0; i < nhisto; i++)
196 /* no samples! is bad!*/
200 nmean /= static_cast<real>(nhisto);
203 for (i = 0; i < nhisto; i++)
205 /* make sure that all points are in the ratio < x < 1/ratio range */
206 if (!((histo[i] / nmean < 1.0 / ratio) && (histo[i] / nmean > ratio)))
215 static gmx_bool CheckIfDoneEquilibrating(int nlim, const t_expanded* expand, const df_history_t* dfhist, int64_t step)
219 gmx_bool bDoneEquilibrating = TRUE;
222 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
223 if (expand->lmc_forced_nstart > 0)
225 for (i = 0; i < nlim; i++)
227 if (dfhist->n_at_lam[i]
228 < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're
229 definitely not done equilibrating*/
231 bDoneEquilibrating = FALSE;
238 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
239 bDoneEquilibrating = TRUE;
241 /* calculate the total number of samples */
242 switch (expand->elmceq)
244 case LambdaWeightWillReachEquilibrium::No:
245 /* We have not equilibrated, and won't, ever. */
246 bDoneEquilibrating = FALSE;
248 case LambdaWeightWillReachEquilibrium::Yes:
249 /* we have equilibrated -- we're done */
250 bDoneEquilibrating = TRUE;
252 case LambdaWeightWillReachEquilibrium::Steps:
253 /* first, check if we are equilibrating by steps, if we're still under */
254 if (step < expand->equil_steps)
256 bDoneEquilibrating = FALSE;
259 case LambdaWeightWillReachEquilibrium::Samples:
261 for (i = 0; i < nlim; i++)
263 totalsamples += dfhist->n_at_lam[i];
265 if (totalsamples < expand->equil_samples)
267 bDoneEquilibrating = FALSE;
270 case LambdaWeightWillReachEquilibrium::NumAtLambda:
271 for (i = 0; i < nlim; i++)
273 if (dfhist->n_at_lam[i]
274 < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're
275 definitely not done equilibrating*/
277 bDoneEquilibrating = FALSE;
282 case LambdaWeightWillReachEquilibrium::WLDelta:
283 if (EWL(expand->elamstats)) /* This check is in readir as well, but
286 if (dfhist->wl_delta > expand->equil_wl_delta)
288 bDoneEquilibrating = FALSE;
292 case LambdaWeightWillReachEquilibrium::Ratio:
293 /* we can use the flatness as a judge of good weights, as long as
294 we're not doing minvar, or Wang-Landau.
295 But turn off for now until we figure out exactly how we do this.
298 if (!(EWL(expand->elamstats) || expand->elamstats == LambdaWeightCalculation::Minvar))
300 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need
301 to convert to floats for this histogram function. */
304 snew(modhisto, nlim);
305 for (i = 0; i < nlim; i++)
307 modhisto[i] = 1.0 * (dfhist->n_at_lam[i] - expand->lmc_forced_nstart);
309 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
313 bDoneEquilibrating = FALSE;
317 default: bDoneEquilibrating = TRUE; break;
320 return bDoneEquilibrating;
323 static gmx_bool UpdateWeights(int nlim,
325 df_history_t* dfhist,
327 const real* scaled_lamee,
328 const real* weighted_lamee,
331 gmx_bool bSufficientSamples;
332 real acceptanceWeight;
334 int min_nvalm, min_nvalp, maxc;
335 real omega_m1_0, omega_p1_0;
336 real zero_sum_weights;
337 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array,
338 *dwp_array, *dwm_array;
339 real clam_varm, clam_varp, clam_osum, clam_weightsm, clam_weightsp, clam_minvar;
340 real * lam_variance, *lam_dg;
344 /* Future potential todos for this function (see #3848):
345 * - Update the names in the dhist structure to be clearer. Not done for now since this
346 * a bugfix update and we are mininizing other code changes.
347 * - Modularize the code some more.
348 * - potentially merge with accelerated weight histogram functionality, since it's very similar.
350 /* if we have equilibrated the expanded ensemble weights, we are not updating them, so 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
374 == LambdaWeightCalculation::WL) /* Using standard Wang-Landau for weight updates */
376 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
377 dfhist->wl_histo[fep_state] += 1.0;
379 else if (expand->elamstats == LambdaWeightCalculation::WWL)
380 /* Using weighted Wang-Landau for weight updates.
381 * Very closly equivalent to accelerated weight histogram approach
382 * applied to expanded ensemble. */
386 /* first increment count */
387 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim - 1);
388 for (i = 0; i < nlim; i++)
390 dfhist->wl_histo[i] += static_cast<real>(p_k[i]);
393 /* then increment weights (uses count) */
395 GenerateWeightedGibbsProbabilities(
396 weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
398 for (i = 0; i < nlim; i++)
400 dfhist->sum_weights[i] -= dfhist->wl_delta * static_cast<real>(p_k[i]);
402 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
407 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
408 dfhist->sum_weights[i] -= log(di);
414 zero_sum_weights = dfhist->sum_weights[0];
415 for (i = 0; i < nlim; i++)
417 dfhist->sum_weights[i] -= zero_sum_weights;
421 if (expand->elamstats == LambdaWeightCalculation::Barker
422 || expand->elamstats == LambdaWeightCalculation::Metropolis
423 || expand->elamstats == LambdaWeightCalculation::Minvar)
425 maxc = 2 * expand->c_range + 1;
428 snew(lam_variance, nlim);
430 snew(omegap_array, maxc);
431 snew(weightsp_array, maxc);
432 snew(varp_array, maxc);
433 snew(dwp_array, maxc);
435 snew(omegam_array, maxc);
436 snew(weightsm_array, maxc);
437 snew(varm_array, maxc);
438 snew(dwm_array, maxc);
440 /* unpack the values of the free energy differences and the
441 * variance in their estimates between nearby lambdas. We will
442 * only actually update 2 of these, the state we are currently
443 * at and the one we end up moving to
446 for (i = 0; i < nlim - 1; i++)
447 { /* only through the second to last */
448 lam_dg[i] = dfhist->sum_dg[i + 1] - dfhist->sum_dg[i];
450 gmx::square(dfhist->sum_variance[i + 1]) - gmx::square(dfhist->sum_variance[i]);
453 /* accumulate running averages of thermodynamic averages for Bennett Acceptance Ratio-based
454 * estimates of the free energy .
455 * Rather than peforming self-consistent estimation of the free energies at each step,
456 * we keep track of an array of possible different free energies (cnvals),
457 * and we self-consistently choose the best one. The one that leads to a free energy estimate
458 * that is closest to itself is the best estimate of the free energy. It is essentially a
459 * parallellized version of self-consistent iteration. maxc is the number of these constants. */
461 for (int nval = 0; nval < maxc; nval++)
463 const real cnval = static_cast<real>(nval - expand->c_range);
465 /* Compute acceptance criterion weight to the state below this one for use in averages.
466 * Note we do not have to have just moved from that state to use this free energy
467 * estimate; these are essentially "virtual" moves. */
471 const auto lambdaEnergyDifference =
472 cnval - (scaled_lamee[fep_state] - scaled_lamee[fep_state - 1]);
474 gmx::calculateAcceptanceWeight(expand->elamstats, lambdaEnergyDifference);
475 dfhist->accum_m[fep_state][nval] += acceptanceWeight;
476 dfhist->accum_m2[fep_state][nval] += acceptanceWeight * acceptanceWeight;
479 // Compute acceptance criterion weight to transition to the next state
480 if (fep_state < nlim - 1)
482 const auto lambdaEnergyDifference =
483 -cnval + (scaled_lamee[fep_state + 1] - scaled_lamee[fep_state]);
485 gmx::calculateAcceptanceWeight(expand->elamstats, lambdaEnergyDifference);
486 dfhist->accum_p[fep_state][nval] += acceptanceWeight;
487 dfhist->accum_p2[fep_state][nval] += acceptanceWeight * acceptanceWeight;
490 /* Determination of Metropolis transition and Barker transition weights */
492 int numObservationsCurrentState = dfhist->n_at_lam[fep_state];
493 /* determine the number of observations above and below the current state */
494 int numObservationsLowerState = 0;
497 numObservationsLowerState = dfhist->n_at_lam[fep_state - 1];
499 int numObservationsHigherState = 0;
500 if (fep_state < nlim - 1)
502 numObservationsHigherState = dfhist->n_at_lam[fep_state + 1];
505 /* Calculate the biases for each expanded ensemble state that minimize the total
506 * variance, as implemented in Martinez-Veracoechea and Escobedo,
507 * J. Phys. Chem. B 2008, 112, 8120-8128
509 * The variance associated with the free energy estimate between two states i and j
511 * Var(i,j) = {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} / numObservations(i->j)
512 * + {avg[xi(j->i)^2] / avg[xi(j->i)]^2 - 1} / numObservations(j->i)
513 * where xi(i->j) is the acceptance factor / weight associated with moving from state i to j
514 * As we are calculating the acceptance factor to the neighbors every time we're visiting
515 * a state, numObservations(i->j) == numObservations(i) and numObservations(j->i) == numObservations(j)
518 /* Accumulation of acceptance weight averages between the current state and the
519 * states +1 (p1) and -1 (m1), averaged at current state (0)
521 real avgAcceptanceCurrentToLower = 0;
522 real avgAcceptanceCurrentToHigher = 0;
523 /* Accumulation of acceptance weight averages quantities between states 0
524 * and states +1 and -1, squared
526 real avgAcceptanceCurrentToLowerSquared = 0;
527 real avgAcceptanceCurrentToHigherSquared = 0;
528 /* Accumulation of free energy quantities from lower state (m1) to current state (0) and squared */
529 real avgAcceptanceLowerToCurrent = 0;
530 real avgAcceptanceLowerToCurrentSquared = 0;
531 /* Accumulation of free energy quantities from upper state (p1) to current state (0) and squared */
532 real avgAcceptanceHigherToCurrent = 0;
533 real avgAcceptanceHigherToCurrentSquared = 0;
535 if (numObservationsCurrentState > 0)
537 avgAcceptanceCurrentToLower = dfhist->accum_m[fep_state][nval] / numObservationsCurrentState;
538 avgAcceptanceCurrentToHigher =
539 dfhist->accum_p[fep_state][nval] / numObservationsCurrentState;
540 avgAcceptanceCurrentToLowerSquared =
541 dfhist->accum_m2[fep_state][nval] / numObservationsCurrentState;
542 avgAcceptanceCurrentToHigherSquared =
543 dfhist->accum_p2[fep_state][nval] / numObservationsCurrentState;
546 if ((fep_state > 0) && (numObservationsLowerState > 0))
548 avgAcceptanceLowerToCurrent =
549 dfhist->accum_p[fep_state - 1][nval] / numObservationsLowerState;
550 avgAcceptanceLowerToCurrentSquared =
551 dfhist->accum_p2[fep_state - 1][nval] / numObservationsLowerState;
554 if ((fep_state < nlim - 1) && (numObservationsHigherState > 0))
556 avgAcceptanceHigherToCurrent =
557 dfhist->accum_m[fep_state + 1][nval] / numObservationsHigherState;
558 avgAcceptanceHigherToCurrentSquared =
559 dfhist->accum_m2[fep_state + 1][nval] / numObservationsHigherState;
561 /* These are accumulation of positive values (see definition of acceptance functions
562 * above), or of squares of positive values.
563 * We're taking this for granted in the following calculation, so make sure
564 * here that nothing weird happened. Although technically all values should be positive,
565 * because of floating point precisions, they might be numerically zero. */
567 avgAcceptanceCurrentToLower >= 0 && avgAcceptanceCurrentToLowerSquared >= 0
568 && avgAcceptanceCurrentToHigher >= 0
569 && avgAcceptanceCurrentToHigherSquared >= 0 && avgAcceptanceLowerToCurrent >= 0
570 && avgAcceptanceLowerToCurrentSquared >= 0 && avgAcceptanceHigherToCurrent >= 0
571 && avgAcceptanceHigherToCurrentSquared >= 0,
572 "By definition, the acceptance factors should all be nonnegative.");
574 real varianceCurrentToLower = 0;
575 real varianceCurrentToHigher = 0;
576 real weightDifferenceToLower = 0;
577 real weightDifferenceToHigher = 0;
578 real varianceToLower = 0;
579 real varianceToHigher = 0;
583 if (numObservationsCurrentState > 0)
585 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
587 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
588 * acceptances are all positive!), and hence
589 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
590 * We're catching that case explicitly to avoid numerical
591 * problems dividing by zero when the overlap between states is small (#3304)
593 if (avgAcceptanceCurrentToLower > 0)
595 varianceCurrentToLower =
596 avgAcceptanceCurrentToLowerSquared
597 / (avgAcceptanceCurrentToLower * avgAcceptanceCurrentToLower)
600 if (numObservationsLowerState > 0)
602 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
604 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
605 * acceptances are all positive!), and hence
606 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
607 * We're catching that case explicitly to avoid numerical
608 * problems dividing by zero when the overlap between states is small (#3304)
610 real varianceLowerToCurrent = 0;
611 if (avgAcceptanceLowerToCurrent > 0)
613 varianceLowerToCurrent =
614 avgAcceptanceLowerToCurrentSquared
615 / (avgAcceptanceLowerToCurrent * avgAcceptanceLowerToCurrent)
618 /* Free energy difference to the state one state lower */
619 /* if these either of these quantities are zero, the energies are */
620 /* way too large for the dynamic range. We need an alternate guesstimate */
621 if ((avgAcceptanceCurrentToLower == 0) || (avgAcceptanceLowerToCurrent == 0))
623 weightDifferenceToLower =
624 (scaled_lamee[fep_state] - scaled_lamee[fep_state - 1]);
628 weightDifferenceToLower = (std::log(avgAcceptanceCurrentToLower)
629 - std::log(avgAcceptanceLowerToCurrent))
632 /* Variance of the free energy difference to the one state lower */
634 (1.0 / numObservationsCurrentState) * (varianceCurrentToLower)
635 + (1.0 / numObservationsLowerState) * (varianceLowerToCurrent);
640 if (fep_state < nlim - 1)
642 if (numObservationsCurrentState > 0)
644 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
646 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
647 * acceptances are all positive!), and hence
648 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
649 * We're catching that case explicitly to avoid numerical
650 * problems dividing by zero when the overlap between states is small (#3304)
653 if (avgAcceptanceCurrentToHigher < 0)
655 varianceCurrentToHigher =
656 avgAcceptanceCurrentToHigherSquared
657 / (avgAcceptanceCurrentToHigher * avgAcceptanceCurrentToHigher)
660 if (numObservationsHigherState > 0)
662 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
664 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
665 * acceptances are all positive!), and hence
666 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
667 * We're catching that case explicitly to avoid numerical
668 * problems dividing by zero when the overlap between states is small (#3304)
670 real varianceHigherToCurrent = 0;
671 if (avgAcceptanceHigherToCurrent > 0)
673 varianceHigherToCurrent =
674 avgAcceptanceHigherToCurrentSquared
675 / (avgAcceptanceHigherToCurrent * avgAcceptanceHigherToCurrent)
678 /* Free energy difference to the state one state higher */
679 /* if these either of these quantities are zero, the energies are */
680 /* way too large for the dynamic range. We need an alternate guesstimate */
681 if ((avgAcceptanceHigherToCurrent == 0) || (avgAcceptanceCurrentToHigher == 0))
683 weightDifferenceToHigher =
684 (scaled_lamee[fep_state + 1] - scaled_lamee[fep_state]);
688 weightDifferenceToHigher = (std::log(avgAcceptanceHigherToCurrent)
689 - std::log(avgAcceptanceCurrentToHigher))
692 /* Variance of the free energy difference to the one state higher */
694 (1.0 / numObservationsHigherState) * (varianceHigherToCurrent)
695 + (1.0 / numObservationsCurrentState) * (varianceCurrentToHigher);
700 if (numObservationsCurrentState > 0)
702 omegam_array[nval] = varianceCurrentToLower;
706 omegam_array[nval] = 0;
708 weightsm_array[nval] = weightDifferenceToLower;
709 varm_array[nval] = varianceToLower;
710 if (numObservationsLowerState > 0)
713 fabs((cnval + std::log((1.0 * numObservationsCurrentState) / numObservationsLowerState))
714 - lam_dg[fep_state - 1]);
718 dwm_array[nval] = std::fabs(cnval - lam_dg[fep_state - 1]);
721 if (numObservationsCurrentState > 0)
723 omegap_array[nval] = varianceCurrentToHigher;
727 omegap_array[nval] = 0;
729 weightsp_array[nval] = weightDifferenceToHigher;
730 varp_array[nval] = varianceToHigher;
731 if ((numObservationsHigherState > 0) && (numObservationsCurrentState > 0))
734 fabs((cnval + std::log((1.0 * numObservationsHigherState) / numObservationsCurrentState))
735 - lam_dg[fep_state]);
739 dwp_array[nval] = std::fabs(cnval - lam_dg[fep_state]);
743 /* find the free energy estimate closest to the guessed weight's value */
745 min_nvalm = FindMinimum(dwm_array, maxc);
746 omega_m1_0 = omegam_array[min_nvalm];
747 clam_weightsm = weightsm_array[min_nvalm];
748 clam_varm = varm_array[min_nvalm];
750 min_nvalp = FindMinimum(dwp_array, maxc);
751 omega_p1_0 = omegap_array[min_nvalp];
752 clam_weightsp = weightsp_array[min_nvalp];
753 clam_varp = varp_array[min_nvalp];
755 clam_osum = omega_m1_0 + omega_p1_0;
759 clam_minvar = 0.5 * std::log(clam_osum);
764 lam_dg[fep_state - 1] = clam_weightsm;
765 lam_variance[fep_state - 1] = clam_varm;
768 if (fep_state < nlim - 1)
770 lam_dg[fep_state] = clam_weightsp;
771 lam_variance[fep_state] = clam_varp;
774 if (expand->elamstats == LambdaWeightCalculation::Minvar)
776 bSufficientSamples = TRUE;
777 /* make sure the number of samples in each state are all
778 * past a user-specified threshold
780 for (i = 0; i < nlim; i++)
782 if (dfhist->n_at_lam[i] < expand->minvarmin)
784 bSufficientSamples = FALSE;
787 if (bSufficientSamples)
789 dfhist->sum_minvar[fep_state] = clam_minvar;
792 for (i = 0; i < nlim; i++)
794 dfhist->sum_minvar[i] += (expand->minvar_const - clam_minvar);
796 expand->minvar_const = clam_minvar;
797 dfhist->sum_minvar[fep_state] = 0.0;
801 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
806 /* we need to rezero minvar now, since it could change at fep_state = 0 */
807 dfhist->sum_dg[0] = 0.0;
808 dfhist->sum_variance[0] = 0.0;
809 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
811 for (i = 1; i < nlim; i++)
813 dfhist->sum_dg[i] = lam_dg[i - 1] + dfhist->sum_dg[i - 1];
814 dfhist->sum_variance[i] =
815 std::sqrt(lam_variance[i - 1] + gmx::square(dfhist->sum_variance[i - 1]));
816 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
823 sfree(weightsm_array);
828 sfree(weightsp_array);
835 static int ChooseNewLambda(int nlim,
836 const t_expanded* expand,
837 df_history_t* dfhist,
839 const real* weighted_lamee,
844 /* Choose new lambda value, and update transition matrix */
846 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
847 real r1, r2, de, trialprob, tprob = 0;
848 double * propose, *accept, *remainder;
851 gmx::ThreeFry2x64<0> rng(
852 seed, gmx::RandomDomain::ExpandedEnsemble); // We only draw once, so zero bits internal counter is fine
853 gmx::UniformRealDistribution<real> dist;
855 starting_fep_state = fep_state;
856 lamnew = fep_state; /* so that there is a default setting -- stays the same */
858 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
860 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim - 1] <= expand->lmc_forced_nstart))
862 /* Use a marching method to run through the lambdas and get preliminary free energy data,
863 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
865 /* if we have enough at this lambda, move on to the next one */
867 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
869 lamnew = fep_state + 1;
870 if (lamnew == nlim) /* whoops, stepped too far! */
885 snew(remainder, nlim);
887 for (i = 0; i < expand->lmc_repeats; i++)
889 rng.restart(step, i);
892 for (ifep = 0; ifep < nlim; ifep++)
898 if ((expand->elmcmove == LambdaMoveCalculation::Gibbs)
899 || (expand->elmcmove == LambdaMoveCalculation::MetropolisGibbs))
901 /* use the Gibbs sampler, with restricted range */
902 if (expand->gibbsdeltalam < 0)
909 minfep = fep_state - expand->gibbsdeltalam;
910 maxfep = fep_state + expand->gibbsdeltalam;
915 if (maxfep > nlim - 1)
921 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
923 if (expand->elmcmove == LambdaMoveCalculation::Gibbs)
925 for (ifep = minfep; ifep <= maxfep; ifep++)
927 propose[ifep] = p_k[ifep];
932 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
934 if (r1 <= p_k[lamnew])
941 else if (expand->elmcmove == LambdaMoveCalculation::MetropolisGibbs)
944 /* Metropolized Gibbs sampling */
945 for (ifep = minfep; ifep <= maxfep; ifep++)
947 remainder[ifep] = 1 - p_k[ifep];
950 /* find the proposal probabilities */
952 if (remainder[fep_state] == 0)
954 /* only the current state has any probability */
955 /* we have to stay at the current state */
960 for (ifep = minfep; ifep <= maxfep; ifep++)
962 if (ifep != fep_state)
964 propose[ifep] = p_k[ifep] / remainder[fep_state];
973 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
975 pnorm = p_k[lamtrial] / remainder[fep_state];
976 if (lamtrial != fep_state)
986 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
988 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
989 trialprob = (remainder[fep_state]) / (remainder[lamtrial]);
990 if (trialprob < tprob)
1005 /* now figure out the acceptance probability for each */
1006 for (ifep = minfep; ifep <= maxfep; ifep++)
1009 if (remainder[ifep] != 0)
1011 trialprob = (remainder[fep_state]) / (remainder[ifep]);
1015 trialprob = 1.0; /* this state is the only choice! */
1017 if (trialprob < tprob)
1021 /* probability for fep_state=0, but that's fine, it's never proposed! */
1022 accept[ifep] = tprob;
1026 if (lamnew > maxfep)
1028 /* it's possible some rounding is failing */
1029 if (gmx_within_tol(remainder[fep_state], 0, 50 * GMX_DOUBLE_EPS))
1031 /* numerical rounding error -- no state other than the original has weight */
1036 /* probably not a numerical issue */
1038 int nerror = 200 + (maxfep - minfep + 1) * 60;
1040 snew(errorstr, nerror);
1041 /* if its greater than maxfep, then something went wrong -- probably underflow
1042 in the calculation of sum weights. Generated detailed info for failure */
1045 "Something wrong in choosing new lambda state with a Gibbs move -- "
1046 "probably underflow in weight determination.\nDenominator is: "
1047 "%3d%17.10e\n i dE numerator weights\n",
1050 for (ifep = minfep; ifep <= maxfep; ifep++)
1052 loc += sprintf(&errorstr[loc],
1053 "%3d %17.10e%17.10e%17.10e\n",
1055 weighted_lamee[ifep],
1057 dfhist->sum_weights[ifep]);
1059 gmx_fatal(FARGS, "%s", errorstr);
1063 else if ((expand->elmcmove == LambdaMoveCalculation::Metropolis)
1064 || (expand->elmcmove == LambdaMoveCalculation::Barker))
1066 /* use the metropolis sampler with trial +/- 1 */
1072 lamtrial = fep_state;
1076 lamtrial = fep_state - 1;
1081 if (fep_state == nlim - 1)
1083 lamtrial = fep_state;
1087 lamtrial = fep_state + 1;
1091 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
1092 if (expand->elmcmove == LambdaMoveCalculation::Metropolis)
1097 tprob = std::exp(de);
1099 propose[fep_state] = 0;
1100 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
1102 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
1103 accept[lamtrial] = tprob;
1105 else if (expand->elmcmove == LambdaMoveCalculation::Barker)
1107 if (de > 0) /* Numerically stable version */
1109 tprob = 1.0 / (1.0 + std::exp(-de));
1113 tprob = std::exp(de) / (std::exp(de) + 1.0);
1115 propose[fep_state] = (1 - tprob);
1116 propose[lamtrial] +=
1117 tprob; /* we add, to account for the fact that at the end, they might be the same point */
1118 accept[fep_state] = 1.0;
1119 accept[lamtrial] = 1.0;
1133 for (ifep = 0; ifep < nlim; ifep++)
1135 dfhist->Tij[fep_state][ifep] += propose[ifep] * accept[ifep];
1136 dfhist->Tij[fep_state][fep_state] += propose[ifep] * (1.0 - accept[ifep]);
1141 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1150 /* print out the weights to the log, along with current state */
1151 void PrintFreeEnergyInfoToFile(FILE* outfile,
1152 const t_lambda* fep,
1153 const t_expanded* expand,
1154 const t_simtemp* simtemp,
1155 const df_history_t* dfhist,
1160 int nlim, ifep, jfep;
1161 real dw, dg, dv, Tprint;
1162 gmx_bool bSimTemp = FALSE;
1164 nlim = fep->n_lambda;
1165 if (simtemp != nullptr)
1170 if (step % frequency == 0)
1172 fprintf(outfile, " MC-lambda information\n");
1173 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1175 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1177 fprintf(outfile, " N");
1178 for (auto i : keysOf(fep->separate_dvdl))
1180 if (fep->separate_dvdl[i])
1182 fprintf(outfile, "%7s", enumValueToString(i));
1184 else if ((i == FreeEnergyPerturbationCouplingType::Temperature) && bSimTemp)
1186 fprintf(outfile, "%10s", enumValueToString(i)); /* more space for temperature formats */
1189 fprintf(outfile, " Count ");
1190 if (expand->elamstats == LambdaWeightCalculation::Minvar)
1192 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1196 fprintf(outfile, "G(in kT) dG(in kT)\n");
1198 for (ifep = 0; ifep < nlim; ifep++)
1200 if (ifep == nlim - 1)
1208 dw = dfhist->sum_weights[ifep + 1] - dfhist->sum_weights[ifep];
1209 dg = dfhist->sum_dg[ifep + 1] - dfhist->sum_dg[ifep];
1210 dv = std::sqrt(gmx::square(dfhist->sum_variance[ifep + 1])
1211 - gmx::square(dfhist->sum_variance[ifep]));
1213 fprintf(outfile, "%3d", (ifep + 1));
1214 for (auto i : keysOf(fep->separate_dvdl))
1216 if (fep->separate_dvdl[i])
1218 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1220 else if (i == FreeEnergyPerturbationCouplingType::Temperature && bSimTemp)
1222 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1225 if (EWL(expand->elamstats)
1226 && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1228 if (expand->elamstats == LambdaWeightCalculation::WL)
1230 fprintf(outfile, " %8d", static_cast<int>(dfhist->wl_histo[ifep]));
1234 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1237 else /* we have equilibrated weights */
1239 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1241 if (expand->elamstats == LambdaWeightCalculation::Minvar)
1244 " %10.5f %10.5f %10.5f %10.5f",
1245 dfhist->sum_weights[ifep],
1246 dfhist->sum_dg[ifep],
1252 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1254 if (ifep == fep_state)
1256 fprintf(outfile, " <<\n");
1260 fprintf(outfile, " \n");
1263 fprintf(outfile, "\n");
1265 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1267 fprintf(outfile, " Transition Matrix\n");
1268 for (ifep = 0; ifep < nlim; ifep++)
1270 fprintf(outfile, "%12d", (ifep + 1));
1272 fprintf(outfile, "\n");
1273 for (ifep = 0; ifep < nlim; ifep++)
1275 for (jfep = 0; jfep < nlim; jfep++)
1277 if (dfhist->n_at_lam[ifep] > 0)
1279 if (expand->bSymmetrizedTMatrix)
1281 Tprint = (dfhist->Tij[ifep][jfep] + dfhist->Tij[jfep][ifep])
1282 / (dfhist->n_at_lam[ifep] + dfhist->n_at_lam[jfep]);
1286 Tprint = (dfhist->Tij[ifep][jfep]) / (dfhist->n_at_lam[ifep]);
1293 fprintf(outfile, "%12.8f", Tprint);
1295 fprintf(outfile, "%3d\n", (ifep + 1));
1298 fprintf(outfile, " Empirical Transition Matrix\n");
1299 for (ifep = 0; ifep < nlim; ifep++)
1301 fprintf(outfile, "%12d", (ifep + 1));
1303 fprintf(outfile, "\n");
1304 for (ifep = 0; ifep < nlim; ifep++)
1306 for (jfep = 0; jfep < nlim; jfep++)
1308 if (dfhist->n_at_lam[ifep] > 0)
1310 if (expand->bSymmetrizedTMatrix)
1312 Tprint = (dfhist->Tij_empirical[ifep][jfep] + dfhist->Tij_empirical[jfep][ifep])
1313 / (dfhist->n_at_lam[ifep] + dfhist->n_at_lam[jfep]);
1317 Tprint = dfhist->Tij_empirical[ifep][jfep] / (dfhist->n_at_lam[ifep]);
1324 fprintf(outfile, "%12.8f", Tprint);
1326 fprintf(outfile, "%3d\n", (ifep + 1));
1332 int expandedEnsembleUpdateLambdaState(FILE* log,
1333 const t_inputrec* ir,
1334 const gmx_enerdata_t* enerd,
1336 df_history_t* dfhist,
1339 real * pfep_lamee, *scaled_lamee, *weighted_lamee;
1341 int i, nlim, lamnew, totalsamples;
1342 real oneovert, maxscaled = 0, maxweighted = 0;
1345 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1347 expand = ir->expandedvals.get();
1348 simtemp = ir->simtempvals.get();
1349 nlim = ir->fepvals->n_lambda;
1351 snew(scaled_lamee, nlim);
1352 snew(weighted_lamee, nlim);
1353 snew(pfep_lamee, nlim);
1356 /* update the count at the current lambda*/
1357 dfhist->n_at_lam[fep_state]++;
1359 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda
1360 state that's pressure controlled.*/
1363 where does this PV term go?
1364 for (i=0;i<nlim;i++)
1366 fep_lamee[i] += pVTerm;
1370 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1371 /* we don't need to include the pressure term, since the volume is the same between the two.
1372 is there some term we are neglecting, however? */
1374 if (ir->efep != FreeEnergyPerturbationType::No)
1376 for (i = 0; i < nlim; i++)
1380 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1382 enerd->foreignLambdaTerms.deltaH(i) / (simtemp->temperatures[i] * gmx::c_boltz)
1383 + enerd->term[F_EPOT]
1384 * (1.0 / (simtemp->temperatures[i])
1385 - 1.0 / (simtemp->temperatures[fep_state]))
1390 scaled_lamee[i] = enerd->foreignLambdaTerms.deltaH(i) / (expand->mc_temp * gmx::c_boltz);
1391 /* mc_temp is currently set to the system reft unless otherwise defined */
1394 /* save these energies for printing, so they don't get overwritten by the next step */
1395 /* they aren't overwritten in the non-free energy case, but we always print with these
1403 for (i = 0; i < nlim; i++)
1407 * (1.0 / simtemp->temperatures[i] - 1.0 / simtemp->temperatures[fep_state])
1413 for (i = 0; i < nlim; i++)
1415 pfep_lamee[i] = scaled_lamee[i];
1417 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1420 maxscaled = scaled_lamee[i];
1421 maxweighted = weighted_lamee[i];
1425 if (scaled_lamee[i] > maxscaled)
1427 maxscaled = scaled_lamee[i];
1429 if (weighted_lamee[i] > maxweighted)
1431 maxweighted = weighted_lamee[i];
1436 for (i = 0; i < nlim; i++)
1438 scaled_lamee[i] -= maxscaled;
1439 weighted_lamee[i] -= maxweighted;
1442 /* update weights - we decide whether or not to actually do this inside */
1444 bDoneEquilibrating =
1445 UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1446 if (bDoneEquilibrating)
1451 "\nStep %" PRId64 ": Weights have equilibrated, using criteria: %s\n",
1453 enumValueToString(expand->elmceq));
1457 lamnew = ChooseNewLambda(
1458 nlim, expand, dfhist, fep_state, weighted_lamee, p_k, ir->expandedvals->lmc_seed, step);
1460 /* now check on the Wang-Landau updating critera */
1462 if (EWL(expand->elamstats))
1464 bSwitchtoOneOverT = FALSE;
1465 if (expand->bWLoneovert)
1468 for (i = 0; i < nlim; i++)
1470 totalsamples += dfhist->n_at_lam[i];
1472 oneovert = (1.0 * nlim) / totalsamples;
1473 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1474 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1475 if ((dfhist->wl_delta <= ((totalsamples) / (totalsamples - 1.00001)) * oneovert)
1476 && (dfhist->wl_delta < expand->init_wl_delta))
1478 bSwitchtoOneOverT = TRUE;
1481 if (bSwitchtoOneOverT)
1484 oneovert; /* now we reduce by this each time, instead of only at flatness */
1488 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1491 for (i = 0; i < nlim; i++)
1493 dfhist->wl_histo[i] = 0;
1495 dfhist->wl_delta *= expand->wl_scale;
1498 fprintf(log, "\nStep %d: weights are now:", static_cast<int>(step));
1499 for (i = 0; i < nlim; i++)
1501 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1509 sfree(scaled_lamee);
1510 sfree(weighted_lamee);
1516 //! Update reference temperature for simulated tempering state change
1517 static void simulatedTemperingUpdateTemperature(t_inputrec* ir,
1522 gmx::ArrayRef<const unsigned short> cTC,
1525 const t_simtemp* simtemp = ir->simtempvals.get();
1526 std::vector<real> buf_ngtc(ir->opts.ngtc);
1528 for (int i = 0; i < ir->opts.ngtc; i++)
1530 if (ir->opts.ref_t[i] > 0)
1532 real told = ir->opts.ref_t[i];
1533 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1534 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i] / told); /* using the buffer as temperature scaling */
1538 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1540 for (int n = 0; n < homenr; n++)
1542 const int gt = cTC.empty() ? 0 : cTC[n];
1543 for (int d = 0; d < DIM; d++)
1545 v[n][d] *= buf_ngtc[gt];
1549 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1551 /* we need to recalculate the masses if the temperature has changed */
1552 init_npt_masses(ir, state, MassQ, FALSE);
1553 for (int i = 0; i < state->nnhpres; i++)
1555 for (int j = 0; j < ir->opts.nhchainlength; j++)
1557 state->nhpres_vxi[i + j] *= buf_ngtc[i];
1560 for (int i = 0; i < ir->opts.ngtc; i++)
1562 for (int j = 0; j < ir->opts.nhchainlength; j++)
1564 state->nosehoover_vxi[i + j] *= buf_ngtc[i];
1570 int ExpandedEnsembleDynamics(FILE* log,
1572 const gmx_enerdata_t* enerd,
1576 df_history_t* dfhist,
1580 gmx::ArrayRef<const unsigned short> cTC)
1581 /* Note that the state variable is only needed for simulated tempering, not
1582 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1584 const int newLambda = expandedEnsembleUpdateLambdaState(log, ir, enerd, fep_state, dfhist, step);
1585 // if using simulated tempering, we need to adjust the temperatures
1586 // only need to change the temperatures if we change the state
1587 if (ir->bSimTemp && (newLambda != fep_state))
1589 simulatedTemperingUpdateTemperature(ir, state, MassQ, v, homenr, cTC, newLambda);