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.
45 #include "gromacs/domdec/domdec.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/xtcio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/listed_forces/disre.h"
52 #include "gromacs/listed_forces/orires.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/calcmu.h"
58 #include "gromacs/mdlib/constr.h"
59 #include "gromacs/mdlib/force.h"
60 #include "gromacs/mdlib/update.h"
61 #include "gromacs/mdtypes/enerdata.h"
62 #include "gromacs/mdtypes/forcerec.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/mdatom.h"
66 #include "gromacs/mdtypes/state.h"
67 #include "gromacs/random/threefry.h"
68 #include "gromacs/random/uniformrealdistribution.h"
69 #include "gromacs/timing/wallcycle.h"
70 #include "gromacs/utility/enumerationhelpers.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxmpi.h"
73 #include "gromacs/utility/smalloc.h"
75 #include "expanded_internal.h"
77 static void init_df_history_weights(df_history_t* dfhist, const t_expanded* expand, int nlim)
80 dfhist->wl_delta = expand->init_wl_delta;
81 for (i = 0; i < nlim; i++)
83 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
84 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
88 /* Eventually should contain all the functions needed to initialize expanded ensemble
89 before the md loop starts */
90 void init_expanded_ensemble(gmx_bool bStateFromCP, const t_inputrec* ir, df_history_t* dfhist)
94 init_df_history_weights(dfhist, ir->expandedvals.get(), ir->fepvals->n_lambda);
98 static void GenerateGibbsProbabilities(const real* ene, double* p_k, double* pks, int minfep, int maxfep)
105 maxene = ene[minfep];
106 /* find the maximum value */
107 for (i = minfep; i <= maxfep; i++)
114 /* find the denominator */
115 for (i = minfep; i <= maxfep; i++)
117 *pks += std::exp(ene[i] - maxene);
120 for (i = minfep; i <= maxfep; i++)
122 p_k[i] = std::exp(ene[i] - maxene) / *pks;
127 GenerateWeightedGibbsProbabilities(const real* ene, double* p_k, double* pks, int nlim, real* nvals, real delta)
136 for (i = 0; i < nlim; i++)
140 /* add the delta, since we need to make sure it's greater than zero, and
141 we need a non-arbitrary number? */
142 nene[i] = ene[i] + std::log(nvals[i] + delta);
146 nene[i] = ene[i] + std::log(nvals[i]);
150 /* find the maximum value */
152 for (i = 0; i < nlim; i++)
154 if (nene[i] > maxene)
160 /* subtract off the maximum, avoiding overflow */
161 for (i = 0; i < nlim; i++)
166 /* find the denominator */
167 for (i = 0; i < nlim; i++)
169 *pks += std::exp(nene[i]);
173 for (i = 0; i < nlim; i++)
175 p_k[i] = std::exp(nene[i]) / *pks;
180 static int FindMinimum(const real* min_metric, int N)
187 min_val = min_metric[0];
189 for (nval = 0; nval < N; nval++)
191 if (min_metric[nval] < min_val)
193 min_val = min_metric[nval];
200 static gmx_bool CheckHistogramRatios(int nhisto, const real* histo, real ratio)
208 for (i = 0; i < nhisto; i++)
215 /* no samples! is bad!*/
219 nmean /= static_cast<real>(nhisto);
222 for (i = 0; i < nhisto; i++)
224 /* make sure that all points are in the ratio < x < 1/ratio range */
225 if (!((histo[i] / nmean < 1.0 / ratio) && (histo[i] / nmean > ratio)))
234 static gmx_bool CheckIfDoneEquilibrating(int nlim, const t_expanded* expand, const df_history_t* dfhist, int64_t step)
238 gmx_bool bDoneEquilibrating = TRUE;
241 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
242 if (expand->lmc_forced_nstart > 0)
244 for (i = 0; i < nlim; i++)
246 if (dfhist->n_at_lam[i]
247 < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're
248 definitely not done equilibrating*/
250 bDoneEquilibrating = FALSE;
257 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
258 bDoneEquilibrating = TRUE;
260 /* calculate the total number of samples */
261 switch (expand->elmceq)
263 case LambdaWeightWillReachEquilibrium::No:
264 /* We have not equilibrated, and won't, ever. */
265 bDoneEquilibrating = FALSE;
267 case LambdaWeightWillReachEquilibrium::Yes:
268 /* we have equilibrated -- we're done */
269 bDoneEquilibrating = TRUE;
271 case LambdaWeightWillReachEquilibrium::Steps:
272 /* first, check if we are equilibrating by steps, if we're still under */
273 if (step < expand->equil_steps)
275 bDoneEquilibrating = FALSE;
278 case LambdaWeightWillReachEquilibrium::Samples:
280 for (i = 0; i < nlim; i++)
282 totalsamples += dfhist->n_at_lam[i];
284 if (totalsamples < expand->equil_samples)
286 bDoneEquilibrating = FALSE;
289 case LambdaWeightWillReachEquilibrium::NumAtLambda:
290 for (i = 0; i < nlim; i++)
292 if (dfhist->n_at_lam[i]
293 < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're
294 definitely not done equilibrating*/
296 bDoneEquilibrating = FALSE;
301 case LambdaWeightWillReachEquilibrium::WLDelta:
302 if (EWL(expand->elamstats)) /* This check is in readir as well, but
305 if (dfhist->wl_delta > expand->equil_wl_delta)
307 bDoneEquilibrating = FALSE;
311 case LambdaWeightWillReachEquilibrium::Ratio:
312 /* we can use the flatness as a judge of good weights, as long as
313 we're not doing minvar, or Wang-Landau.
314 But turn off for now until we figure out exactly how we do this.
317 if (!(EWL(expand->elamstats) || expand->elamstats == LambdaWeightCalculation::Minvar))
319 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need
320 to convert to floats for this histogram function. */
323 snew(modhisto, nlim);
324 for (i = 0; i < nlim; i++)
326 modhisto[i] = 1.0 * (dfhist->n_at_lam[i] - expand->lmc_forced_nstart);
328 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
332 bDoneEquilibrating = FALSE;
336 default: bDoneEquilibrating = TRUE; break;
339 return bDoneEquilibrating;
342 static gmx_bool UpdateWeights(int nlim,
344 df_history_t* dfhist,
346 const real* scaled_lamee,
347 const real* weighted_lamee,
350 gmx_bool bSufficientSamples;
351 real acceptanceWeight;
353 int min_nvalm, min_nvalp, maxc;
354 real omega_m1_0, omega_p1_0;
355 real zero_sum_weights;
356 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array,
357 *dwp_array, *dwm_array;
358 real clam_varm, clam_varp, clam_osum, clam_weightsm, clam_weightsp, clam_minvar;
359 real * lam_variance, *lam_dg;
363 /* Future potential todos for this function (see #3848):
364 * - Update the names in the dhist structure to be clearer. Not done for now since this
365 * a bugfix update and we are mininizing other code changes.
366 * - Modularize the code some more.
367 * - potentially merge with accelerated weight histogram functionality, since it's very similar.
369 /* if we have equilibrated the expanded ensemble weights, we are not updating them, so exit now */
375 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
377 dfhist->bEquil = TRUE;
378 /* zero out the visited states so we know how many equilibrated states we have
380 for (i = 0; i < nlim; i++)
382 dfhist->n_at_lam[i] = 0;
387 /* If we reached this far, we have not equilibrated yet, keep on
388 going resetting the weights */
390 if (EWL(expand->elamstats))
392 if (expand->elamstats
393 == LambdaWeightCalculation::WL) /* Using standard Wang-Landau for weight updates */
395 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
396 dfhist->wl_histo[fep_state] += 1.0;
398 else if (expand->elamstats == LambdaWeightCalculation::WWL)
399 /* Using weighted Wang-Landau for weight updates.
400 * Very closly equivalent to accelerated weight histogram approach
401 * applied to expanded ensemble. */
405 /* first increment count */
406 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim - 1);
407 for (i = 0; i < nlim; i++)
409 dfhist->wl_histo[i] += static_cast<real>(p_k[i]);
412 /* then increment weights (uses count) */
414 GenerateWeightedGibbsProbabilities(
415 weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
417 for (i = 0; i < nlim; i++)
419 dfhist->sum_weights[i] -= dfhist->wl_delta * static_cast<real>(p_k[i]);
421 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
426 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
427 dfhist->sum_weights[i] -= log(di);
433 zero_sum_weights = dfhist->sum_weights[0];
434 for (i = 0; i < nlim; i++)
436 dfhist->sum_weights[i] -= zero_sum_weights;
440 if (expand->elamstats == LambdaWeightCalculation::Barker
441 || expand->elamstats == LambdaWeightCalculation::Metropolis
442 || expand->elamstats == LambdaWeightCalculation::Minvar)
444 maxc = 2 * expand->c_range + 1;
447 snew(lam_variance, nlim);
449 snew(omegap_array, maxc);
450 snew(weightsp_array, maxc);
451 snew(varp_array, maxc);
452 snew(dwp_array, maxc);
454 snew(omegam_array, maxc);
455 snew(weightsm_array, maxc);
456 snew(varm_array, maxc);
457 snew(dwm_array, maxc);
459 /* unpack the values of the free energy differences and the
460 * variance in their estimates between nearby lambdas. We will
461 * only actually update 2 of these, the state we are currently
462 * at and the one we end up moving to
465 for (i = 0; i < nlim - 1; i++)
466 { /* only through the second to last */
467 lam_dg[i] = dfhist->sum_dg[i + 1] - dfhist->sum_dg[i];
469 gmx::square(dfhist->sum_variance[i + 1]) - gmx::square(dfhist->sum_variance[i]);
472 /* accumulate running averages of thermodynamic averages for Bennett Acceptance Ratio-based
473 * estimates of the free energy .
474 * Rather than peforming self-consistent estimation of the free energies at each step,
475 * we keep track of an array of possible different free energies (cnvals),
476 * and we self-consistently choose the best one. The one that leads to a free energy estimate
477 * that is closest to itself is the best estimate of the free energy. It is essentially a
478 * parallellized version of self-consistent iteration. maxc is the number of these constants. */
480 for (int nval = 0; nval < maxc; nval++)
482 const real cnval = static_cast<real>(nval - expand->c_range);
484 /* Compute acceptance criterion weight to the state below this one for use in averages.
485 * Note we do not have to have just moved from that state to use this free energy
486 * estimate; these are essentially "virtual" moves. */
490 const auto lambdaEnergyDifference =
491 cnval - (scaled_lamee[fep_state] - scaled_lamee[fep_state - 1]);
493 gmx::calculateAcceptanceWeight(expand->elamstats, lambdaEnergyDifference);
494 dfhist->accum_m[fep_state][nval] += acceptanceWeight;
495 dfhist->accum_m2[fep_state][nval] += acceptanceWeight * acceptanceWeight;
498 // Compute acceptance criterion weight to transition to the next state
499 if (fep_state < nlim - 1)
501 const auto lambdaEnergyDifference =
502 -cnval + (scaled_lamee[fep_state + 1] - scaled_lamee[fep_state]);
504 gmx::calculateAcceptanceWeight(expand->elamstats, lambdaEnergyDifference);
505 dfhist->accum_p[fep_state][nval] += acceptanceWeight;
506 dfhist->accum_p2[fep_state][nval] += acceptanceWeight * acceptanceWeight;
509 /* Determination of Metropolis transition and Barker transition weights */
511 int numObservationsCurrentState = dfhist->n_at_lam[fep_state];
512 /* determine the number of observations above and below the current state */
513 int numObservationsLowerState = 0;
516 numObservationsLowerState = dfhist->n_at_lam[fep_state - 1];
518 int numObservationsHigherState = 0;
519 if (fep_state < nlim - 1)
521 numObservationsHigherState = dfhist->n_at_lam[fep_state + 1];
524 /* Calculate the biases for each expanded ensemble state that minimize the total
525 * variance, as implemented in Martinez-Veracoechea and Escobedo,
526 * J. Phys. Chem. B 2008, 112, 8120-8128
528 * The variance associated with the free energy estimate between two states i and j
530 * Var(i,j) = {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} / numObservations(i->j)
531 * + {avg[xi(j->i)^2] / avg[xi(j->i)]^2 - 1} / numObservations(j->i)
532 * where xi(i->j) is the acceptance factor / weight associated with moving from state i to j
533 * As we are calculating the acceptance factor to the neighbors every time we're visiting
534 * a state, numObservations(i->j) == numObservations(i) and numObservations(j->i) == numObservations(j)
537 /* Accumulation of acceptance weight averages between the current state and the
538 * states +1 (p1) and -1 (m1), averaged at current state (0)
540 real avgAcceptanceCurrentToLower = 0;
541 real avgAcceptanceCurrentToHigher = 0;
542 /* Accumulation of acceptance weight averages quantities between states 0
543 * and states +1 and -1, squared
545 real avgAcceptanceCurrentToLowerSquared = 0;
546 real avgAcceptanceCurrentToHigherSquared = 0;
547 /* Accumulation of free energy quantities from lower state (m1) to current state (0) and squared */
548 real avgAcceptanceLowerToCurrent = 0;
549 real avgAcceptanceLowerToCurrentSquared = 0;
550 /* Accumulation of free energy quantities from upper state (p1) to current state (0) and squared */
551 real avgAcceptanceHigherToCurrent = 0;
552 real avgAcceptanceHigherToCurrentSquared = 0;
554 if (numObservationsCurrentState > 0)
556 avgAcceptanceCurrentToLower = dfhist->accum_m[fep_state][nval] / numObservationsCurrentState;
557 avgAcceptanceCurrentToHigher =
558 dfhist->accum_p[fep_state][nval] / numObservationsCurrentState;
559 avgAcceptanceCurrentToLowerSquared =
560 dfhist->accum_m2[fep_state][nval] / numObservationsCurrentState;
561 avgAcceptanceCurrentToHigherSquared =
562 dfhist->accum_p2[fep_state][nval] / numObservationsCurrentState;
565 if ((fep_state > 0) && (numObservationsLowerState > 0))
567 avgAcceptanceLowerToCurrent =
568 dfhist->accum_p[fep_state - 1][nval] / numObservationsLowerState;
569 avgAcceptanceLowerToCurrentSquared =
570 dfhist->accum_p2[fep_state - 1][nval] / numObservationsLowerState;
573 if ((fep_state < nlim - 1) && (numObservationsHigherState > 0))
575 avgAcceptanceHigherToCurrent =
576 dfhist->accum_m[fep_state + 1][nval] / numObservationsHigherState;
577 avgAcceptanceHigherToCurrentSquared =
578 dfhist->accum_m2[fep_state + 1][nval] / numObservationsHigherState;
580 /* These are accumulation of positive values (see definition of acceptance functions
581 * above), or of squares of positive values.
582 * We're taking this for granted in the following calculation, so make sure
583 * here that nothing weird happened. Although technically all values should be positive,
584 * because of floating point precisions, they might be numerically zero. */
586 avgAcceptanceCurrentToLower >= 0 && avgAcceptanceCurrentToLowerSquared >= 0
587 && avgAcceptanceCurrentToHigher >= 0
588 && avgAcceptanceCurrentToHigherSquared >= 0 && avgAcceptanceLowerToCurrent >= 0
589 && avgAcceptanceLowerToCurrentSquared >= 0 && avgAcceptanceHigherToCurrent >= 0
590 && avgAcceptanceHigherToCurrentSquared >= 0,
591 "By definition, the acceptance factors should all be nonnegative.");
593 real varianceCurrentToLower = 0;
594 real varianceCurrentToHigher = 0;
595 real weightDifferenceToLower = 0;
596 real weightDifferenceToHigher = 0;
597 real varianceToLower = 0;
598 real varianceToHigher = 0;
602 if (numObservationsCurrentState > 0)
604 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
606 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
607 * acceptances are all positive!), and hence
608 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
609 * We're catching that case explicitly to avoid numerical
610 * problems dividing by zero when the overlap between states is small (#3304)
612 if (avgAcceptanceCurrentToLower > 0)
614 varianceCurrentToLower =
615 avgAcceptanceCurrentToLowerSquared
616 / (avgAcceptanceCurrentToLower * avgAcceptanceCurrentToLower)
619 if (numObservationsLowerState > 0)
621 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
623 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
624 * acceptances are all positive!), and hence
625 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
626 * We're catching that case explicitly to avoid numerical
627 * problems dividing by zero when the overlap between states is small (#3304)
629 real varianceLowerToCurrent = 0;
630 if (avgAcceptanceLowerToCurrent > 0)
632 varianceLowerToCurrent =
633 avgAcceptanceLowerToCurrentSquared
634 / (avgAcceptanceLowerToCurrent * avgAcceptanceLowerToCurrent)
637 /* Free energy difference to the state one state lower */
638 /* if these either of these quantities are zero, the energies are */
639 /* way too large for the dynamic range. We need an alternate guesstimate */
640 if ((avgAcceptanceCurrentToLower == 0) || (avgAcceptanceLowerToCurrent == 0))
642 weightDifferenceToLower =
643 (scaled_lamee[fep_state] - scaled_lamee[fep_state - 1]);
647 weightDifferenceToLower = (std::log(avgAcceptanceCurrentToLower)
648 - std::log(avgAcceptanceLowerToCurrent))
651 /* Variance of the free energy difference to the one state lower */
653 (1.0 / numObservationsCurrentState) * (varianceCurrentToLower)
654 + (1.0 / numObservationsLowerState) * (varianceLowerToCurrent);
659 if (fep_state < nlim - 1)
661 if (numObservationsCurrentState > 0)
663 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
665 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
666 * acceptances are all positive!), and hence
667 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
668 * We're catching that case explicitly to avoid numerical
669 * problems dividing by zero when the overlap between states is small (#3304)
672 if (avgAcceptanceCurrentToHigher < 0)
674 varianceCurrentToHigher =
675 avgAcceptanceCurrentToHigherSquared
676 / (avgAcceptanceCurrentToHigher * avgAcceptanceCurrentToHigher)
679 if (numObservationsHigherState > 0)
681 /* Calculate {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1}
683 * Note that if avg[xi(i->j)] == 0, also avg[xi(i->j)^2] == 0 (since the
684 * acceptances are all positive!), and hence
685 * {avg[xi(i->j)^2] / avg[xi(i->j)]^2 - 1} -> 0 for avg[xi(i->j)] -> 0
686 * We're catching that case explicitly to avoid numerical
687 * problems dividing by zero when the overlap between states is small (#3304)
689 real varianceHigherToCurrent = 0;
690 if (avgAcceptanceHigherToCurrent > 0)
692 varianceHigherToCurrent =
693 avgAcceptanceHigherToCurrentSquared
694 / (avgAcceptanceHigherToCurrent * avgAcceptanceHigherToCurrent)
697 /* Free energy difference to the state one state higher */
698 /* if these either of these quantities are zero, the energies are */
699 /* way too large for the dynamic range. We need an alternate guesstimate */
700 if ((avgAcceptanceHigherToCurrent == 0) || (avgAcceptanceCurrentToHigher == 0))
702 weightDifferenceToHigher =
703 (scaled_lamee[fep_state + 1] - scaled_lamee[fep_state]);
707 weightDifferenceToHigher = (std::log(avgAcceptanceHigherToCurrent)
708 - std::log(avgAcceptanceCurrentToHigher))
711 /* Variance of the free energy difference to the one state higher */
713 (1.0 / numObservationsHigherState) * (varianceHigherToCurrent)
714 + (1.0 / numObservationsCurrentState) * (varianceCurrentToHigher);
719 if (numObservationsCurrentState > 0)
721 omegam_array[nval] = varianceCurrentToLower;
725 omegam_array[nval] = 0;
727 weightsm_array[nval] = weightDifferenceToLower;
728 varm_array[nval] = varianceToLower;
729 if (numObservationsLowerState > 0)
732 fabs((cnval + std::log((1.0 * numObservationsCurrentState) / numObservationsLowerState))
733 - lam_dg[fep_state - 1]);
737 dwm_array[nval] = std::fabs(cnval - lam_dg[fep_state - 1]);
740 if (numObservationsCurrentState > 0)
742 omegap_array[nval] = varianceCurrentToHigher;
746 omegap_array[nval] = 0;
748 weightsp_array[nval] = weightDifferenceToHigher;
749 varp_array[nval] = varianceToHigher;
750 if ((numObservationsHigherState > 0) && (numObservationsCurrentState > 0))
753 fabs((cnval + std::log((1.0 * numObservationsHigherState) / numObservationsCurrentState))
754 - lam_dg[fep_state]);
758 dwp_array[nval] = std::fabs(cnval - lam_dg[fep_state]);
762 /* find the free energy estimate closest to the guessed weight's value */
764 min_nvalm = FindMinimum(dwm_array, maxc);
765 omega_m1_0 = omegam_array[min_nvalm];
766 clam_weightsm = weightsm_array[min_nvalm];
767 clam_varm = varm_array[min_nvalm];
769 min_nvalp = FindMinimum(dwp_array, maxc);
770 omega_p1_0 = omegap_array[min_nvalp];
771 clam_weightsp = weightsp_array[min_nvalp];
772 clam_varp = varp_array[min_nvalp];
774 clam_osum = omega_m1_0 + omega_p1_0;
778 clam_minvar = 0.5 * std::log(clam_osum);
783 lam_dg[fep_state - 1] = clam_weightsm;
784 lam_variance[fep_state - 1] = clam_varm;
787 if (fep_state < nlim - 1)
789 lam_dg[fep_state] = clam_weightsp;
790 lam_variance[fep_state] = clam_varp;
793 if (expand->elamstats == LambdaWeightCalculation::Minvar)
795 bSufficientSamples = TRUE;
796 /* make sure the number of samples in each state are all
797 * past a user-specified threshold
799 for (i = 0; i < nlim; i++)
801 if (dfhist->n_at_lam[i] < expand->minvarmin)
803 bSufficientSamples = FALSE;
806 if (bSufficientSamples)
808 dfhist->sum_minvar[fep_state] = clam_minvar;
811 for (i = 0; i < nlim; i++)
813 dfhist->sum_minvar[i] += (expand->minvar_const - clam_minvar);
815 expand->minvar_const = clam_minvar;
816 dfhist->sum_minvar[fep_state] = 0.0;
820 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
825 /* we need to rezero minvar now, since it could change at fep_state = 0 */
826 dfhist->sum_dg[0] = 0.0;
827 dfhist->sum_variance[0] = 0.0;
828 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
830 for (i = 1; i < nlim; i++)
832 dfhist->sum_dg[i] = lam_dg[i - 1] + dfhist->sum_dg[i - 1];
833 dfhist->sum_variance[i] =
834 std::sqrt(lam_variance[i - 1] + gmx::square(dfhist->sum_variance[i - 1]));
835 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
842 sfree(weightsm_array);
847 sfree(weightsp_array);
854 static int ChooseNewLambda(int nlim,
855 const t_expanded* expand,
856 df_history_t* dfhist,
858 const real* weighted_lamee,
863 /* Choose new lambda value, and update transition matrix */
865 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
866 real r1, r2, de, trialprob, tprob = 0;
867 double * propose, *accept, *remainder;
870 gmx::ThreeFry2x64<0> rng(
871 seed, gmx::RandomDomain::ExpandedEnsemble); // We only draw once, so zero bits internal counter is fine
872 gmx::UniformRealDistribution<real> dist;
874 starting_fep_state = fep_state;
875 lamnew = fep_state; /* so that there is a default setting -- stays the same */
877 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
879 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim - 1] <= expand->lmc_forced_nstart))
881 /* Use a marching method to run through the lambdas and get preliminary free energy data,
882 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
884 /* if we have enough at this lambda, move on to the next one */
886 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
888 lamnew = fep_state + 1;
889 if (lamnew == nlim) /* whoops, stepped too far! */
904 snew(remainder, nlim);
906 for (i = 0; i < expand->lmc_repeats; i++)
908 rng.restart(step, i);
911 for (ifep = 0; ifep < nlim; ifep++)
917 if ((expand->elmcmove == LambdaMoveCalculation::Gibbs)
918 || (expand->elmcmove == LambdaMoveCalculation::MetropolisGibbs))
920 /* use the Gibbs sampler, with restricted range */
921 if (expand->gibbsdeltalam < 0)
928 minfep = fep_state - expand->gibbsdeltalam;
929 maxfep = fep_state + expand->gibbsdeltalam;
934 if (maxfep > nlim - 1)
940 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
942 if (expand->elmcmove == LambdaMoveCalculation::Gibbs)
944 for (ifep = minfep; ifep <= maxfep; ifep++)
946 propose[ifep] = p_k[ifep];
951 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
953 if (r1 <= p_k[lamnew])
960 else if (expand->elmcmove == LambdaMoveCalculation::MetropolisGibbs)
963 /* Metropolized Gibbs sampling */
964 for (ifep = minfep; ifep <= maxfep; ifep++)
966 remainder[ifep] = 1 - p_k[ifep];
969 /* find the proposal probabilities */
971 if (remainder[fep_state] == 0)
973 /* only the current state has any probability */
974 /* we have to stay at the current state */
979 for (ifep = minfep; ifep <= maxfep; ifep++)
981 if (ifep != fep_state)
983 propose[ifep] = p_k[ifep] / remainder[fep_state];
992 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
994 pnorm = p_k[lamtrial] / remainder[fep_state];
995 if (lamtrial != fep_state)
1005 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
1007 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
1008 trialprob = (remainder[fep_state]) / (remainder[lamtrial]);
1009 if (trialprob < tprob)
1024 /* now figure out the acceptance probability for each */
1025 for (ifep = minfep; ifep <= maxfep; ifep++)
1028 if (remainder[ifep] != 0)
1030 trialprob = (remainder[fep_state]) / (remainder[ifep]);
1034 trialprob = 1.0; /* this state is the only choice! */
1036 if (trialprob < tprob)
1040 /* probability for fep_state=0, but that's fine, it's never proposed! */
1041 accept[ifep] = tprob;
1045 if (lamnew > maxfep)
1047 /* it's possible some rounding is failing */
1048 if (gmx_within_tol(remainder[fep_state], 0, 50 * GMX_DOUBLE_EPS))
1050 /* numerical rounding error -- no state other than the original has weight */
1055 /* probably not a numerical issue */
1057 int nerror = 200 + (maxfep - minfep + 1) * 60;
1059 snew(errorstr, nerror);
1060 /* if its greater than maxfep, then something went wrong -- probably underflow
1061 in the calculation of sum weights. Generated detailed info for failure */
1064 "Something wrong in choosing new lambda state with a Gibbs move -- "
1065 "probably underflow in weight determination.\nDenominator is: "
1066 "%3d%17.10e\n i dE numerator weights\n",
1069 for (ifep = minfep; ifep <= maxfep; ifep++)
1071 loc += sprintf(&errorstr[loc],
1072 "%3d %17.10e%17.10e%17.10e\n",
1074 weighted_lamee[ifep],
1076 dfhist->sum_weights[ifep]);
1078 gmx_fatal(FARGS, "%s", errorstr);
1082 else if ((expand->elmcmove == LambdaMoveCalculation::Metropolis)
1083 || (expand->elmcmove == LambdaMoveCalculation::Barker))
1085 /* use the metropolis sampler with trial +/- 1 */
1091 lamtrial = fep_state;
1095 lamtrial = fep_state - 1;
1100 if (fep_state == nlim - 1)
1102 lamtrial = fep_state;
1106 lamtrial = fep_state + 1;
1110 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
1111 if (expand->elmcmove == LambdaMoveCalculation::Metropolis)
1116 tprob = std::exp(de);
1118 propose[fep_state] = 0;
1119 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
1121 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
1122 accept[lamtrial] = tprob;
1124 else if (expand->elmcmove == LambdaMoveCalculation::Barker)
1126 if (de > 0) /* Numerically stable version */
1128 tprob = 1.0 / (1.0 + std::exp(-de));
1132 tprob = std::exp(de) / (std::exp(de) + 1.0);
1134 propose[fep_state] = (1 - tprob);
1135 propose[lamtrial] +=
1136 tprob; /* we add, to account for the fact that at the end, they might be the same point */
1137 accept[fep_state] = 1.0;
1138 accept[lamtrial] = 1.0;
1152 for (ifep = 0; ifep < nlim; ifep++)
1154 dfhist->Tij[fep_state][ifep] += propose[ifep] * accept[ifep];
1155 dfhist->Tij[fep_state][fep_state] += propose[ifep] * (1.0 - accept[ifep]);
1160 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1169 /* print out the weights to the log, along with current state */
1170 void PrintFreeEnergyInfoToFile(FILE* outfile,
1171 const t_lambda* fep,
1172 const t_expanded* expand,
1173 const t_simtemp* simtemp,
1174 const df_history_t* dfhist,
1179 int nlim, ifep, jfep;
1180 real dw, dg, dv, Tprint;
1181 gmx_bool bSimTemp = FALSE;
1183 nlim = fep->n_lambda;
1184 if (simtemp != nullptr)
1189 if (step % frequency == 0)
1191 fprintf(outfile, " MC-lambda information\n");
1192 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1194 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1196 fprintf(outfile, " N");
1197 for (auto i : keysOf(fep->separate_dvdl))
1199 if (fep->separate_dvdl[i])
1201 fprintf(outfile, "%7s", enumValueToString(i));
1203 else if ((i == FreeEnergyPerturbationCouplingType::Temperature) && bSimTemp)
1205 fprintf(outfile, "%10s", enumValueToString(i)); /* more space for temperature formats */
1208 fprintf(outfile, " Count ");
1209 if (expand->elamstats == LambdaWeightCalculation::Minvar)
1211 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1215 fprintf(outfile, "G(in kT) dG(in kT)\n");
1217 for (ifep = 0; ifep < nlim; ifep++)
1219 if (ifep == nlim - 1)
1227 dw = dfhist->sum_weights[ifep + 1] - dfhist->sum_weights[ifep];
1228 dg = dfhist->sum_dg[ifep + 1] - dfhist->sum_dg[ifep];
1229 dv = std::sqrt(gmx::square(dfhist->sum_variance[ifep + 1])
1230 - gmx::square(dfhist->sum_variance[ifep]));
1232 fprintf(outfile, "%3d", (ifep + 1));
1233 for (auto i : keysOf(fep->separate_dvdl))
1235 if (fep->separate_dvdl[i])
1237 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1239 else if (i == FreeEnergyPerturbationCouplingType::Temperature && bSimTemp)
1241 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1244 if (EWL(expand->elamstats)
1245 && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1247 if (expand->elamstats == LambdaWeightCalculation::WL)
1249 fprintf(outfile, " %8d", static_cast<int>(dfhist->wl_histo[ifep]));
1253 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1256 else /* we have equilibrated weights */
1258 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1260 if (expand->elamstats == LambdaWeightCalculation::Minvar)
1263 " %10.5f %10.5f %10.5f %10.5f",
1264 dfhist->sum_weights[ifep],
1265 dfhist->sum_dg[ifep],
1271 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1273 if (ifep == fep_state)
1275 fprintf(outfile, " <<\n");
1279 fprintf(outfile, " \n");
1282 fprintf(outfile, "\n");
1284 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1286 fprintf(outfile, " Transition Matrix\n");
1287 for (ifep = 0; ifep < nlim; ifep++)
1289 fprintf(outfile, "%12d", (ifep + 1));
1291 fprintf(outfile, "\n");
1292 for (ifep = 0; ifep < nlim; ifep++)
1294 for (jfep = 0; jfep < nlim; jfep++)
1296 if (dfhist->n_at_lam[ifep] > 0)
1298 if (expand->bSymmetrizedTMatrix)
1300 Tprint = (dfhist->Tij[ifep][jfep] + dfhist->Tij[jfep][ifep])
1301 / (dfhist->n_at_lam[ifep] + dfhist->n_at_lam[jfep]);
1305 Tprint = (dfhist->Tij[ifep][jfep]) / (dfhist->n_at_lam[ifep]);
1312 fprintf(outfile, "%12.8f", Tprint);
1314 fprintf(outfile, "%3d\n", (ifep + 1));
1317 fprintf(outfile, " Empirical Transition Matrix\n");
1318 for (ifep = 0; ifep < nlim; ifep++)
1320 fprintf(outfile, "%12d", (ifep + 1));
1322 fprintf(outfile, "\n");
1323 for (ifep = 0; ifep < nlim; ifep++)
1325 for (jfep = 0; jfep < nlim; jfep++)
1327 if (dfhist->n_at_lam[ifep] > 0)
1329 if (expand->bSymmetrizedTMatrix)
1331 Tprint = (dfhist->Tij_empirical[ifep][jfep] + dfhist->Tij_empirical[jfep][ifep])
1332 / (dfhist->n_at_lam[ifep] + dfhist->n_at_lam[jfep]);
1336 Tprint = dfhist->Tij_empirical[ifep][jfep] / (dfhist->n_at_lam[ifep]);
1343 fprintf(outfile, "%12.8f", Tprint);
1345 fprintf(outfile, "%3d\n", (ifep + 1));
1351 int ExpandedEnsembleDynamics(FILE* log,
1353 const gmx_enerdata_t* enerd,
1357 df_history_t* dfhist,
1360 const t_mdatoms* mdatoms)
1361 /* Note that the state variable is only needed for simulated tempering, not
1362 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1364 real * pfep_lamee, *scaled_lamee, *weighted_lamee;
1366 int i, nlim, lamnew, totalsamples;
1367 real oneovert, maxscaled = 0, maxweighted = 0;
1370 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1372 expand = ir->expandedvals.get();
1373 simtemp = ir->simtempvals.get();
1374 nlim = ir->fepvals->n_lambda;
1376 snew(scaled_lamee, nlim);
1377 snew(weighted_lamee, nlim);
1378 snew(pfep_lamee, nlim);
1381 /* update the count at the current lambda*/
1382 dfhist->n_at_lam[fep_state]++;
1384 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda
1385 state that's pressure controlled.*/
1388 where does this PV term go?
1389 for (i=0;i<nlim;i++)
1391 fep_lamee[i] += pVTerm;
1395 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1396 /* we don't need to include the pressure term, since the volume is the same between the two.
1397 is there some term we are neglecting, however? */
1399 if (ir->efep != FreeEnergyPerturbationType::No)
1401 for (i = 0; i < nlim; i++)
1405 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1407 enerd->foreignLambdaTerms.deltaH(i) / (simtemp->temperatures[i] * gmx::c_boltz)
1408 + enerd->term[F_EPOT]
1409 * (1.0 / (simtemp->temperatures[i])
1410 - 1.0 / (simtemp->temperatures[fep_state]))
1415 scaled_lamee[i] = enerd->foreignLambdaTerms.deltaH(i) / (expand->mc_temp * gmx::c_boltz);
1416 /* mc_temp is currently set to the system reft unless otherwise defined */
1419 /* save these energies for printing, so they don't get overwritten by the next step */
1420 /* they aren't overwritten in the non-free energy case, but we always print with these
1428 for (i = 0; i < nlim; i++)
1432 * (1.0 / simtemp->temperatures[i] - 1.0 / simtemp->temperatures[fep_state])
1438 for (i = 0; i < nlim; i++)
1440 pfep_lamee[i] = scaled_lamee[i];
1442 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1445 maxscaled = scaled_lamee[i];
1446 maxweighted = weighted_lamee[i];
1450 if (scaled_lamee[i] > maxscaled)
1452 maxscaled = scaled_lamee[i];
1454 if (weighted_lamee[i] > maxweighted)
1456 maxweighted = weighted_lamee[i];
1461 for (i = 0; i < nlim; i++)
1463 scaled_lamee[i] -= maxscaled;
1464 weighted_lamee[i] -= maxweighted;
1467 /* update weights - we decide whether or not to actually do this inside */
1469 bDoneEquilibrating =
1470 UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1471 if (bDoneEquilibrating)
1476 "\nStep %" PRId64 ": Weights have equilibrated, using criteria: %s\n",
1478 enumValueToString(expand->elmceq));
1482 lamnew = ChooseNewLambda(
1483 nlim, expand, dfhist, fep_state, weighted_lamee, p_k, ir->expandedvals->lmc_seed, step);
1484 /* if using simulated tempering, we need to adjust the temperatures */
1485 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1490 int nstart, nend, gt;
1492 snew(buf_ngtc, ir->opts.ngtc);
1494 for (i = 0; i < ir->opts.ngtc; i++)
1496 if (ir->opts.ref_t[i] > 0)
1498 told = ir->opts.ref_t[i];
1499 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1500 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i] / told); /* using the buffer as temperature scaling */
1504 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1507 nend = mdatoms->homenr;
1508 for (n = nstart; n < nend; n++)
1513 gt = mdatoms->cTC[n];
1515 for (d = 0; d < DIM; d++)
1517 v[n][d] *= buf_ngtc[gt];
1521 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1523 /* we need to recalculate the masses if the temperature has changed */
1524 init_npt_masses(ir, state, MassQ, FALSE);
1525 for (i = 0; i < state->nnhpres; i++)
1527 for (j = 0; j < ir->opts.nhchainlength; j++)
1529 state->nhpres_vxi[i + j] *= buf_ngtc[i];
1532 for (i = 0; i < ir->opts.ngtc; i++)
1534 for (j = 0; j < ir->opts.nhchainlength; j++)
1536 state->nosehoover_vxi[i + j] *= buf_ngtc[i];
1543 /* now check on the Wang-Landau updating critera */
1545 if (EWL(expand->elamstats))
1547 bSwitchtoOneOverT = FALSE;
1548 if (expand->bWLoneovert)
1551 for (i = 0; i < nlim; i++)
1553 totalsamples += dfhist->n_at_lam[i];
1555 oneovert = (1.0 * nlim) / totalsamples;
1556 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1557 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1558 if ((dfhist->wl_delta <= ((totalsamples) / (totalsamples - 1.00001)) * oneovert)
1559 && (dfhist->wl_delta < expand->init_wl_delta))
1561 bSwitchtoOneOverT = TRUE;
1564 if (bSwitchtoOneOverT)
1567 oneovert; /* now we reduce by this each time, instead of only at flatness */
1571 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1574 for (i = 0; i < nlim; i++)
1576 dfhist->wl_histo[i] = 0;
1578 dfhist->wl_delta *= expand->wl_scale;
1581 fprintf(log, "\nStep %d: weights are now:", static_cast<int>(step));
1582 for (i = 0; i < nlim; i++)
1584 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1592 sfree(scaled_lamee);
1593 sfree(weighted_lamee);