Split txtdump.*, part 1
[gromacs.git] / src / gromacs / mdlib / expanded.cpp
blobf79f9f375599b6dd039a6cc6c6f0a0201d92e30a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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.
35 #include "gmxpre.h"
37 #include <stdio.h>
39 #include <cmath>
41 #include <algorithm>
43 #include "gromacs/domdec/domdec.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/fileio/xtcio.h"
47 #include "gromacs/gmxlib/chargegroup.h"
48 #include "gromacs/gmxlib/disre.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/gmxlib/orires.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdlib/calcmu.h"
55 #include "gromacs/mdlib/constr.h"
56 #include "gromacs/mdlib/force.h"
57 #include "gromacs/mdlib/mdrun.h"
58 #include "gromacs/mdlib/update.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/random/random.h"
62 #include "gromacs/timing/wallcycle.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxmpi.h"
65 #include "gromacs/utility/smalloc.h"
67 static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
69 int i;
70 dfhist->wl_delta = expand->init_wl_delta;
71 for (i = 0; i < nlim; i++)
73 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
74 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
78 /* Eventually should contain all the functions needed to initialize expanded ensemble
79 before the md loop starts */
80 extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
82 if (!bStateFromCP)
84 init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
88 static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
91 int i;
92 real maxene;
94 *pks = 0.0;
95 maxene = ene[minfep];
96 /* find the maximum value */
97 for (i = minfep; i <= maxfep; i++)
99 if (ene[i] > maxene)
101 maxene = ene[i];
104 /* find the denominator */
105 for (i = minfep; i <= maxfep; i++)
107 *pks += std::exp(ene[i]-maxene);
109 /*numerators*/
110 for (i = minfep; i <= maxfep; i++)
112 p_k[i] = std::exp(ene[i]-maxene) / *pks;
116 static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
119 int i;
120 real maxene;
121 real *nene;
122 *pks = 0.0;
124 snew(nene, nlim);
125 for (i = 0; i < nlim; i++)
127 if (nvals[i] == 0)
129 /* add the delta, since we need to make sure it's greater than zero, and
130 we need a non-arbitrary number? */
131 nene[i] = ene[i] + std::log(nvals[i]+delta);
133 else
135 nene[i] = ene[i] + std::log(nvals[i]);
139 /* find the maximum value */
140 maxene = nene[0];
141 for (i = 0; i < nlim; i++)
143 if (nene[i] > maxene)
145 maxene = nene[i];
149 /* subtract off the maximum, avoiding overflow */
150 for (i = 0; i < nlim; i++)
152 nene[i] -= maxene;
155 /* find the denominator */
156 for (i = 0; i < nlim; i++)
158 *pks += std::exp(nene[i]);
161 /*numerators*/
162 for (i = 0; i < nlim; i++)
164 p_k[i] = std::exp(nene[i]) / *pks;
166 sfree(nene);
169 real do_logsum(int N, real *a_n)
172 /* RETURN VALUE */
173 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
174 real maxarg;
175 real sum;
176 int i;
177 real logsum;
178 /* compute maximum argument to exp(.) */
180 maxarg = a_n[0];
181 for (i = 1; i < N; i++)
183 maxarg = std::max(maxarg, a_n[i]);
186 /* compute sum of exp(a_n - maxarg) */
187 sum = 0.0;
188 for (i = 0; i < N; i++)
190 sum = sum + std::exp(a_n[i] - maxarg);
193 /* compute log sum */
194 logsum = std::log(sum) + maxarg;
195 return logsum;
198 int FindMinimum(real *min_metric, int N)
201 real min_val;
202 int min_nval, nval;
204 min_nval = 0;
205 min_val = min_metric[0];
207 for (nval = 0; nval < N; nval++)
209 if (min_metric[nval] < min_val)
211 min_val = min_metric[nval];
212 min_nval = nval;
215 return min_nval;
218 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
221 int i;
222 real nmean;
223 gmx_bool bIfFlat;
225 nmean = 0;
226 for (i = 0; i < nhisto; i++)
228 nmean += histo[i];
231 if (nmean == 0)
233 /* no samples! is bad!*/
234 bIfFlat = FALSE;
235 return bIfFlat;
237 nmean /= (real)nhisto;
239 bIfFlat = TRUE;
240 for (i = 0; i < nhisto; i++)
242 /* make sure that all points are in the ratio < x < 1/ratio range */
243 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
245 bIfFlat = FALSE;
246 break;
249 return bIfFlat;
252 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
255 int i, totalsamples;
256 gmx_bool bDoneEquilibrating = TRUE;
257 gmx_bool bIfFlat;
259 /* If we are doing slow growth to get initial values, we haven't finished equilibrating */
260 if (expand->lmc_forced_nstart > 0)
262 for (i = 0; i < nlim; i++)
264 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
265 done equilibrating*/
267 bDoneEquilibrating = FALSE;
268 break;
272 else
274 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
275 bDoneEquilibrating = TRUE;
277 /* calculate the total number of samples */
278 switch (expand->elmceq)
280 case elmceqNO:
281 /* We have not equilibrated, and won't, ever. */
282 bDoneEquilibrating = FALSE;
283 break;
284 case elmceqYES:
285 /* we have equilibrated -- we're done */
286 bDoneEquilibrating = TRUE;
287 break;
288 case elmceqSTEPS:
289 /* first, check if we are equilibrating by steps, if we're still under */
290 if (step < expand->equil_steps)
292 bDoneEquilibrating = FALSE;
294 break;
295 case elmceqSAMPLES:
296 totalsamples = 0;
297 for (i = 0; i < nlim; i++)
299 totalsamples += dfhist->n_at_lam[i];
301 if (totalsamples < expand->equil_samples)
303 bDoneEquilibrating = FALSE;
305 break;
306 case elmceqNUMATLAM:
307 for (i = 0; i < nlim; i++)
309 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
310 done equilibrating*/
312 bDoneEquilibrating = FALSE;
313 break;
316 break;
317 case elmceqWLDELTA:
318 if (EWL(expand->elamstats)) /* This check is in readir as well, but
319 just to be sure */
321 if (dfhist->wl_delta > expand->equil_wl_delta)
323 bDoneEquilibrating = FALSE;
326 break;
327 case elmceqRATIO:
328 /* we can use the flatness as a judge of good weights, as long as
329 we're not doing minvar, or Wang-Landau.
330 But turn off for now until we figure out exactly how we do this.
333 if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
335 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
336 floats for this histogram function. */
338 real *modhisto;
339 snew(modhisto, nlim);
340 for (i = 0; i < nlim; i++)
342 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
344 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
345 sfree(modhisto);
346 if (!bIfFlat)
348 bDoneEquilibrating = FALSE;
351 break;
352 default:
353 bDoneEquilibrating = TRUE;
354 break;
357 return bDoneEquilibrating;
360 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
361 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
363 gmx_bool bSufficientSamples;
364 int i;
365 int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
366 real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
367 real de, de_function;
368 real cnval, zero_sum_weights;
369 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
370 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
371 real *lam_variance, *lam_dg;
372 double *p_k;
373 double pks = 0;
374 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;
376 /* if we have equilibrated the weights, exit now */
377 if (dfhist->bEquil)
379 return FALSE;
382 if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
384 dfhist->bEquil = TRUE;
385 /* zero out the visited states so we know how many equilibrated states we have
386 from here on out.*/
387 for (i = 0; i < nlim; i++)
389 dfhist->n_at_lam[i] = 0;
391 return TRUE;
394 /* If we reached this far, we have not equilibrated yet, keep on
395 going resetting the weights */
397 if (EWL(expand->elamstats))
399 if (expand->elamstats == elamstatsWL) /* Standard Wang-Landau */
401 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
402 dfhist->wl_histo[fep_state] += 1.0;
404 else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
406 snew(p_k, nlim);
408 /* first increment count */
409 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
410 for (i = 0; i < nlim; i++)
412 dfhist->wl_histo[i] += (real)p_k[i];
415 /* then increment weights (uses count) */
416 pks = 0.0;
417 GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
419 for (i = 0; i < nlim; i++)
421 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
423 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
425 real di;
426 for (i=0;i<nlim;i++)
428 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
429 dfhist->sum_weights[i] -= log(di);
432 sfree(p_k);
435 zero_sum_weights = dfhist->sum_weights[0];
436 for (i = 0; i < nlim; i++)
438 dfhist->sum_weights[i] -= zero_sum_weights;
442 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
445 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
446 maxc = 2*expand->c_range+1;
448 snew(lam_dg, nlim);
449 snew(lam_variance, nlim);
451 snew(omegap_array, maxc);
452 snew(weightsp_array, maxc);
453 snew(varp_array, maxc);
454 snew(dwp_array, maxc);
456 snew(omegam_array, maxc);
457 snew(weightsm_array, maxc);
458 snew(varm_array, maxc);
459 snew(dwm_array, maxc);
461 /* unpack the current lambdas -- we will only update 2 of these */
463 for (i = 0; i < nlim-1; i++)
464 { /* only through the second to last */
465 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
466 lam_variance[i] = sqr(dfhist->sum_variance[i+1]) - sqr(dfhist->sum_variance[i]);
469 /* accumulate running averages */
470 for (nval = 0; nval < maxc; nval++)
472 /* constants for later use */
473 cnval = (real)(nval-expand->c_range);
474 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
475 if (fep_state > 0)
477 de = std::exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
478 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
480 de_function = 1.0/(1.0+de);
482 else if (expand->elamstats == elamstatsMETROPOLIS)
484 if (de < 1.0)
486 de_function = 1.0;
488 else
490 de_function = 1.0/de;
493 dfhist->accum_m[fep_state][nval] += de_function;
494 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
497 if (fep_state < nlim-1)
499 de = std::exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
500 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
502 de_function = 1.0/(1.0+de);
504 else if (expand->elamstats == elamstatsMETROPOLIS)
506 if (de < 1.0)
508 de_function = 1.0;
510 else
512 de_function = 1.0/de;
515 dfhist->accum_p[fep_state][nval] += de_function;
516 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
519 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
521 n0 = dfhist->n_at_lam[fep_state];
522 if (fep_state > 0)
524 nm1 = dfhist->n_at_lam[fep_state-1];
526 else
528 nm1 = 0;
530 if (fep_state < nlim-1)
532 np1 = dfhist->n_at_lam[fep_state+1];
534 else
536 np1 = 0;
539 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
540 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;
542 if (n0 > 0)
544 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
545 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
546 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
547 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
550 if ((fep_state > 0 ) && (nm1 > 0))
552 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
553 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
556 if ((fep_state < nlim-1) && (np1 > 0))
558 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
559 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
562 omega_m1_0 = 0;
563 omega_p1_0 = 0;
564 clam_weightsm = 0;
565 clam_weightsp = 0;
566 clam_varm = 0;
567 clam_varp = 0;
569 if (fep_state > 0)
571 if (n0 > 0)
573 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
575 if (nm1 > 0)
577 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
579 if ((n0 > 0) && (nm1 > 0))
581 clam_weightsm = (std::log(chi_m1_0) - std::log(chi_p1_m1)) + cnval;
582 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
586 if (fep_state < nlim-1)
588 if (n0 > 0)
590 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
592 if (np1 > 0)
594 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
596 if ((n0 > 0) && (np1 > 0))
598 clam_weightsp = (std::log(chi_m1_p1) - std::log(chi_p1_0)) + cnval;
599 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
603 if (n0 > 0)
605 omegam_array[nval] = omega_m1_0;
607 else
609 omegam_array[nval] = 0;
611 weightsm_array[nval] = clam_weightsm;
612 varm_array[nval] = clam_varm;
613 if (nm1 > 0)
615 dwm_array[nval] = fabs( (cnval + std::log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
617 else
619 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
622 if (n0 > 0)
624 omegap_array[nval] = omega_p1_0;
626 else
628 omegap_array[nval] = 0;
630 weightsp_array[nval] = clam_weightsp;
631 varp_array[nval] = clam_varp;
632 if ((np1 > 0) && (n0 > 0))
634 dwp_array[nval] = fabs( (cnval + std::log((1.0*np1)/n0)) - lam_dg[fep_state] );
636 else
638 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
643 /* find the C's closest to the old weights value */
645 min_nvalm = FindMinimum(dwm_array, maxc);
646 omega_m1_0 = omegam_array[min_nvalm];
647 clam_weightsm = weightsm_array[min_nvalm];
648 clam_varm = varm_array[min_nvalm];
650 min_nvalp = FindMinimum(dwp_array, maxc);
651 omega_p1_0 = omegap_array[min_nvalp];
652 clam_weightsp = weightsp_array[min_nvalp];
653 clam_varp = varp_array[min_nvalp];
655 clam_osum = omega_m1_0 + omega_p1_0;
656 clam_minvar = 0;
657 if (clam_osum > 0)
659 clam_minvar = 0.5*std::log(clam_osum);
662 if (fep_state > 0)
664 lam_dg[fep_state-1] = clam_weightsm;
665 lam_variance[fep_state-1] = clam_varm;
668 if (fep_state < nlim-1)
670 lam_dg[fep_state] = clam_weightsp;
671 lam_variance[fep_state] = clam_varp;
674 if (expand->elamstats == elamstatsMINVAR)
676 bSufficientSamples = TRUE;
677 /* make sure they are all past a threshold */
678 for (i = 0; i < nlim; i++)
680 if (dfhist->n_at_lam[i] < expand->minvarmin)
682 bSufficientSamples = FALSE;
685 if (bSufficientSamples)
687 dfhist->sum_minvar[fep_state] = clam_minvar;
688 if (fep_state == 0)
690 for (i = 0; i < nlim; i++)
692 dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
694 expand->minvar_const = clam_minvar;
695 dfhist->sum_minvar[fep_state] = 0.0;
697 else
699 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
704 /* we need to rezero minvar now, since it could change at fep_state = 0 */
705 dfhist->sum_dg[0] = 0.0;
706 dfhist->sum_variance[0] = 0.0;
707 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
709 for (i = 1; i < nlim; i++)
711 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
712 dfhist->sum_variance[i] = std::sqrt(lam_variance[i-1] + sqr(dfhist->sum_variance[i-1]));
713 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
716 sfree(lam_dg);
717 sfree(lam_variance);
719 sfree(omegam_array);
720 sfree(weightsm_array);
721 sfree(varm_array);
722 sfree(dwm_array);
724 sfree(omegap_array);
725 sfree(weightsp_array);
726 sfree(varp_array);
727 sfree(dwp_array);
729 return FALSE;
732 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
733 gmx_int64_t seed, gmx_int64_t step)
735 /* Choose new lambda value, and update transition matrix */
737 int i, ifep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
738 real r1, r2, de, trialprob, tprob = 0;
739 double *propose, *accept, *remainder;
740 double pks;
741 real pnorm;
743 starting_fep_state = fep_state;
744 lamnew = fep_state; /* so that there is a default setting -- stays the same */
746 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
748 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
750 /* Use a marching method to run through the lambdas and get preliminary free energy data,
751 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
753 /* if we have enough at this lambda, move on to the next one */
755 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
757 lamnew = fep_state+1;
758 if (lamnew == nlim) /* whoops, stepped too far! */
760 lamnew -= 1;
763 else
765 lamnew = fep_state;
767 return lamnew;
771 snew(propose, nlim);
772 snew(accept, nlim);
773 snew(remainder, nlim);
775 for (i = 0; i < expand->lmc_repeats; i++)
777 double rnd[2];
779 gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
781 for (ifep = 0; ifep < nlim; ifep++)
783 propose[ifep] = 0;
784 accept[ifep] = 0;
787 if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
789 /* use the Gibbs sampler, with restricted range */
790 if (expand->gibbsdeltalam < 0)
792 minfep = 0;
793 maxfep = nlim-1;
795 else
797 minfep = fep_state - expand->gibbsdeltalam;
798 maxfep = fep_state + expand->gibbsdeltalam;
799 if (minfep < 0)
801 minfep = 0;
803 if (maxfep > nlim-1)
805 maxfep = nlim-1;
809 GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
811 if (expand->elmcmove == elmcmoveGIBBS)
813 for (ifep = minfep; ifep <= maxfep; ifep++)
815 propose[ifep] = p_k[ifep];
816 accept[ifep] = 1.0;
818 /* Gibbs sampling */
819 r1 = rnd[0];
820 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
822 if (r1 <= p_k[lamnew])
824 break;
826 r1 -= p_k[lamnew];
829 else if (expand->elmcmove == elmcmoveMETGIBBS)
832 /* Metropolized Gibbs sampling */
833 for (ifep = minfep; ifep <= maxfep; ifep++)
835 remainder[ifep] = 1 - p_k[ifep];
838 /* find the proposal probabilities */
840 if (remainder[fep_state] == 0)
842 /* only the current state has any probability */
843 /* we have to stay at the current state */
844 lamnew = fep_state;
846 else
848 for (ifep = minfep; ifep <= maxfep; ifep++)
850 if (ifep != fep_state)
852 propose[ifep] = p_k[ifep]/remainder[fep_state];
854 else
856 propose[ifep] = 0;
860 r1 = rnd[0];
861 for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
863 pnorm = p_k[lamtrial]/remainder[fep_state];
864 if (lamtrial != fep_state)
866 if (r1 <= pnorm)
868 break;
870 r1 -= pnorm;
874 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
875 tprob = 1.0;
876 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
877 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
878 if (trialprob < tprob)
880 tprob = trialprob;
882 r2 = rnd[1];
883 if (r2 < tprob)
885 lamnew = lamtrial;
887 else
889 lamnew = fep_state;
893 /* now figure out the acceptance probability for each */
894 for (ifep = minfep; ifep <= maxfep; ifep++)
896 tprob = 1.0;
897 if (remainder[ifep] != 0)
899 trialprob = (remainder[fep_state])/(remainder[ifep]);
901 else
903 trialprob = 1.0; /* this state is the only choice! */
905 if (trialprob < tprob)
907 tprob = trialprob;
909 /* probability for fep_state=0, but that's fine, it's never proposed! */
910 accept[ifep] = tprob;
914 if (lamnew > maxfep)
916 /* it's possible some rounding is failing */
917 if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
919 /* numerical rounding error -- no state other than the original has weight */
920 lamnew = fep_state;
922 else
924 /* probably not a numerical issue */
925 int loc = 0;
926 int nerror = 200+(maxfep-minfep+1)*60;
927 char *errorstr;
928 snew(errorstr, nerror);
929 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
930 of sum weights. Generated detailed info for failure */
931 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);
932 for (ifep = minfep; ifep <= maxfep; ifep++)
934 loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
936 gmx_fatal(FARGS, errorstr);
940 else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
942 /* use the metropolis sampler with trial +/- 1 */
943 r1 = rnd[0];
944 if (r1 < 0.5)
946 if (fep_state == 0)
948 lamtrial = fep_state;
950 else
952 lamtrial = fep_state-1;
955 else
957 if (fep_state == nlim-1)
959 lamtrial = fep_state;
961 else
963 lamtrial = fep_state+1;
967 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
968 if (expand->elmcmove == elmcmoveMETROPOLIS)
970 tprob = 1.0;
971 trialprob = std::exp(de);
972 if (trialprob < tprob)
974 tprob = trialprob;
976 propose[fep_state] = 0;
977 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
978 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
979 accept[lamtrial] = tprob;
982 else if (expand->elmcmove == elmcmoveBARKER)
984 tprob = 1.0/(1.0+std::exp(-de));
986 propose[fep_state] = (1-tprob);
987 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
988 accept[fep_state] = 1.0;
989 accept[lamtrial] = 1.0;
992 r2 = rnd[1];
993 if (r2 < tprob)
995 lamnew = lamtrial;
997 else
999 lamnew = fep_state;
1003 for (ifep = 0; ifep < nlim; ifep++)
1005 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
1006 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1008 fep_state = lamnew;
1011 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1013 sfree(propose);
1014 sfree(accept);
1015 sfree(remainder);
1017 return lamnew;
1020 /* print out the weights to the log, along with current state */
1021 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1022 int fep_state, int frequency, gmx_int64_t step)
1024 int nlim, i, ifep, jfep;
1025 real dw, dg, dv, Tprint;
1026 const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1027 gmx_bool bSimTemp = FALSE;
1029 nlim = fep->n_lambda;
1030 if (simtemp != NULL)
1032 bSimTemp = TRUE;
1035 if (step % frequency == 0)
1037 fprintf(outfile, " MC-lambda information\n");
1038 if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1040 fprintf(outfile, " Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1042 fprintf(outfile, " N");
1043 for (i = 0; i < efptNR; i++)
1045 if (fep->separate_dvdl[i])
1047 fprintf(outfile, "%7s", print_names[i]);
1049 else if ((i == efptTEMPERATURE) && bSimTemp)
1051 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1054 fprintf(outfile, " Count ");
1055 if (expand->elamstats == elamstatsMINVAR)
1057 fprintf(outfile, "W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1059 else
1061 fprintf(outfile, "G(in kT) dG(in kT)\n");
1063 for (ifep = 0; ifep < nlim; ifep++)
1065 if (ifep == nlim-1)
1067 dw = 0.0;
1068 dg = 0.0;
1069 dv = 0.0;
1071 else
1073 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1074 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1075 dv = std::sqrt(sqr(dfhist->sum_variance[ifep+1]) - sqr(dfhist->sum_variance[ifep]));
1077 fprintf(outfile, "%3d", (ifep+1));
1078 for (i = 0; i < efptNR; i++)
1080 if (fep->separate_dvdl[i])
1082 fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1084 else if (i == efptTEMPERATURE && bSimTemp)
1086 fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1089 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1091 if (expand->elamstats == elamstatsWL)
1093 fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1095 else
1097 fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1100 else /* we have equilibrated weights */
1102 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1104 if (expand->elamstats == elamstatsMINVAR)
1106 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1108 else
1110 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1112 if (ifep == fep_state)
1114 fprintf(outfile, " <<\n");
1116 else
1118 fprintf(outfile, " \n");
1121 fprintf(outfile, "\n");
1123 if ((step % expand->nstTij == 0) && (expand->nstTij > 0) && (step > 0))
1125 fprintf(outfile, " Transition Matrix\n");
1126 for (ifep = 0; ifep < nlim; ifep++)
1128 fprintf(outfile, "%12d", (ifep+1));
1130 fprintf(outfile, "\n");
1131 for (ifep = 0; ifep < nlim; ifep++)
1133 for (jfep = 0; jfep < nlim; jfep++)
1135 if (dfhist->n_at_lam[ifep] > 0)
1137 if (expand->bSymmetrizedTMatrix)
1139 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1141 else
1143 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1146 else
1148 Tprint = 0.0;
1150 fprintf(outfile, "%12.8f", Tprint);
1152 fprintf(outfile, "%3d\n", (ifep+1));
1155 fprintf(outfile, " Empirical Transition Matrix\n");
1156 for (ifep = 0; ifep < nlim; ifep++)
1158 fprintf(outfile, "%12d", (ifep+1));
1160 fprintf(outfile, "\n");
1161 for (ifep = 0; ifep < nlim; ifep++)
1163 for (jfep = 0; jfep < nlim; jfep++)
1165 if (dfhist->n_at_lam[ifep] > 0)
1167 if (expand->bSymmetrizedTMatrix)
1169 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1171 else
1173 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1176 else
1178 Tprint = 0.0;
1180 fprintf(outfile, "%12.8f", Tprint);
1182 fprintf(outfile, "%3d\n", (ifep+1));
1188 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1189 t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
1190 gmx_int64_t step,
1191 rvec *v, t_mdatoms *mdatoms)
1192 /* Note that the state variable is only needed for simulated tempering, not
1193 Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
1195 real *pfep_lamee, *scaled_lamee, *weighted_lamee;
1196 double *p_k;
1197 int i, nlim, lamnew, totalsamples;
1198 real oneovert, maxscaled = 0, maxweighted = 0;
1199 t_expanded *expand;
1200 t_simtemp *simtemp;
1201 gmx_bool bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1203 expand = ir->expandedvals;
1204 simtemp = ir->simtempvals;
1205 nlim = ir->fepvals->n_lambda;
1207 snew(scaled_lamee, nlim);
1208 snew(weighted_lamee, nlim);
1209 snew(pfep_lamee, nlim);
1210 snew(p_k, nlim);
1212 /* update the count at the current lambda*/
1213 dfhist->n_at_lam[fep_state]++;
1215 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1216 pressure controlled.*/
1218 pVTerm = 0;
1219 where does this PV term go?
1220 for (i=0;i<nlim;i++)
1222 fep_lamee[i] += pVTerm;
1226 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1227 /* we don't need to include the pressure term, since the volume is the same between the two.
1228 is there some term we are neglecting, however? */
1230 if (ir->efep != efepNO)
1232 for (i = 0; i < nlim; i++)
1234 if (ir->bSimTemp)
1236 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1237 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1238 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
1240 else
1242 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1243 /* mc_temp is currently set to the system reft unless otherwise defined */
1246 /* save these energies for printing, so they don't get overwritten by the next step */
1247 /* they aren't overwritten in the non-free energy case, but we always print with these
1248 for simplicity */
1251 else
1253 if (ir->bSimTemp)
1255 for (i = 0; i < nlim; i++)
1257 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
1262 for (i = 0; i < nlim; i++)
1264 pfep_lamee[i] = scaled_lamee[i];
1266 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1267 if (i == 0)
1269 maxscaled = scaled_lamee[i];
1270 maxweighted = weighted_lamee[i];
1272 else
1274 if (scaled_lamee[i] > maxscaled)
1276 maxscaled = scaled_lamee[i];
1278 if (weighted_lamee[i] > maxweighted)
1280 maxweighted = weighted_lamee[i];
1285 for (i = 0; i < nlim; i++)
1287 scaled_lamee[i] -= maxscaled;
1288 weighted_lamee[i] -= maxweighted;
1291 /* update weights - we decide whether or not to actually do this inside */
1293 bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
1294 if (bDoneEquilibrating)
1296 if (log)
1298 fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1302 lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
1303 ir->expandedvals->lmc_seed, step);
1304 /* if using simulated tempering, we need to adjust the temperatures */
1305 if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
1307 int i, j, n, d;
1308 real *buf_ngtc;
1309 real told;
1310 int nstart, nend, gt;
1312 snew(buf_ngtc, ir->opts.ngtc);
1314 for (i = 0; i < ir->opts.ngtc; i++)
1316 if (ir->opts.ref_t[i] > 0)
1318 told = ir->opts.ref_t[i];
1319 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1320 buf_ngtc[i] = std::sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1324 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1326 nstart = 0;
1327 nend = mdatoms->homenr;
1328 for (n = nstart; n < nend; n++)
1330 gt = 0;
1331 if (mdatoms->cTC)
1333 gt = mdatoms->cTC[n];
1335 for (d = 0; d < DIM; d++)
1337 v[n][d] *= buf_ngtc[gt];
1341 if (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir))
1343 /* we need to recalculate the masses if the temperature has changed */
1344 init_npt_masses(ir, state, MassQ, FALSE);
1345 for (i = 0; i < state->nnhpres; i++)
1347 for (j = 0; j < ir->opts.nhchainlength; j++)
1349 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1352 for (i = 0; i < ir->opts.ngtc; i++)
1354 for (j = 0; j < ir->opts.nhchainlength; j++)
1356 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1360 sfree(buf_ngtc);
1363 /* now check on the Wang-Landau updating critera */
1365 if (EWL(expand->elamstats))
1367 bSwitchtoOneOverT = FALSE;
1368 if (expand->bWLoneovert)
1370 totalsamples = 0;
1371 for (i = 0; i < nlim; i++)
1373 totalsamples += dfhist->n_at_lam[i];
1375 oneovert = (1.0*nlim)/totalsamples;
1376 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1377 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1378 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1379 (dfhist->wl_delta < expand->init_wl_delta))
1381 bSwitchtoOneOverT = TRUE;
1384 if (bSwitchtoOneOverT)
1386 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1388 else
1390 bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1391 if (bIfReset)
1393 for (i = 0; i < nlim; i++)
1395 dfhist->wl_histo[i] = 0;
1397 dfhist->wl_delta *= expand->wl_scale;
1398 if (log)
1400 fprintf(log, "\nStep %d: weights are now:", (int)step);
1401 for (i = 0; i < nlim; i++)
1403 fprintf(log, " %.5f", dfhist->sum_weights[i]);
1405 fprintf(log, "\n");
1410 sfree(pfep_lamee);
1411 sfree(scaled_lamee);
1412 sfree(weighted_lamee);
1413 sfree(p_k);
1415 return lamnew;