1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROwing Monsters And Cloning Shrimps
40 #include<catamount/dclock.h>
46 #ifdef HAVE_SYS_TIME_H
59 #include "chargegroup.h"
80 #include "gmx_random.h"
83 #include "gmx_wallcycle.h"
93 void GenerateGibbsProbabilities(real *ene, real *p_k, real *pks, int minfep, int maxfep) {
100 /* find the maximum value */
101 for (i=minfep;i<=maxfep;i++)
108 /* find the denominator */
109 for (i=minfep;i<=maxfep;i++)
111 *pks += exp(ene[i]-maxene);
114 for (i=minfep;i<=maxfep;i++)
116 p_k[i] = exp(ene[i]-maxene) / *pks;
120 void GenerateWeightedGibbsProbabilities(real *ene, real *p_k, real *pks, int nlim, real *nvals,real delta) {
128 for (i=0;i<nlim;i++) {
130 /* add the delta, since we need to make sure it's greater than zero, and
131 we need a non-arbitrary number? */
132 nene[i] = ene[i] + log(nvals[i]+delta);
134 nene[i] = ene[i] + log(nvals[i]);
138 /* find the maximum value */
142 if (nene[i] > maxene) {
147 /* subtract off the maximum, avoiding overflow */
153 /* find the denominator */
156 *pks += exp(nene[i]);
162 p_k[i] = exp(nene[i]) / *pks;
167 real do_logsum(int N, real *a_n) {
170 /* log(\sum_{i=0}^(N-1) exp[a_n]) */
175 /* compute maximum argument to exp(.) */
180 maxarg = max(maxarg,a_n[i]);
183 /* compute sum of exp(a_n - maxarg) */
187 sum = sum + exp(a_n[i] - maxarg);
190 /* compute log sum */
191 logsum = log(sum) + maxarg;
195 int FindMinimum(real *min_metric, int N) {
201 min_val = min_metric[0];
203 for (nval=0; nval<N; nval++)
205 if (min_metric[nval] < min_val)
207 min_val = min_metric[nval];
214 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
222 for (i=0;i<nhisto;i++)
229 /* no samples! is bad!*/
233 nmean /= (real)nhisto;
236 for (i=0;i<nhisto;i++)
238 /* make sure that all points are in the ratio < x < 1/ratio range */
239 if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
248 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_large_int_t step)
252 gmx_bool bDoneEquilibrating=TRUE;
255 /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
257 /* calculate the total number of samples */
258 switch (expand->elmceq)
261 /* We have not equilibrated, and won't, ever. */
264 /* we have equilibrated -- we're done */
267 /* first, check if we are equilibrating by steps, if we're still under */
268 if (step < expand->equil_steps)
270 bDoneEquilibrating = FALSE;
277 totalsamples += dfhist->n_at_lam[i];
279 if (totalsamples < expand->equil_samples)
281 bDoneEquilibrating = FALSE;
287 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
290 bDoneEquilibrating = FALSE;
296 if (EWL(expand->elamstats)) /* This check is in readir as well, but
299 if (dfhist->wl_delta > expand->equil_wl_delta)
301 bDoneEquilibrating = FALSE;
306 /* we can use the flatness as a judge of good weights, as long as
307 we're not doing minvar, or Wang-Landau.
308 But turn off for now until we figure out exactly how we do this.
311 if (!(EWL(expand->elamstats) || expand->elamstats==elamstatsMINVAR))
313 /* we want to use flatness -avoiding- the forced-through samples. Plus, we need to convert to
314 floats for this histogram function. */
320 modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
322 bIfFlat = CheckHistogramRatios(nlim,modhisto,expand->equil_ratio);
326 bDoneEquilibrating = FALSE;
330 bDoneEquilibrating = TRUE;
332 /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
334 if (expand->lmc_forced_nstart > 0)
338 if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
341 bDoneEquilibrating = FALSE;
346 return bDoneEquilibrating;
349 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
350 int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_large_int_t step)
352 real maxdiff = 0.000000001;
353 gmx_bool bSufficientSamples;
354 int i, k, n, nz, indexi, indexk, min_n, max_n, nlam, totali;
355 int n0,np1,nm1,nval,min_nvalm,min_nvalp,maxc;
356 real omega_m1_0,omega_p1_m1,omega_m1_p1,omega_p1_0,clam_osum;
357 real de,de_function,dr,denom,maxdr,pks=0;
358 real min_val,cnval,zero_sum_weights;
359 real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
360 real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
361 real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg, *p_k;
362 real *numweighted_lamee, *logfrac;
364 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;
366 /* if we have equilibrated the weights, exit now */
372 if (CheckIfDoneEquilibrating(nlim,expand,dfhist,step))
374 dfhist->bEquil = TRUE;
375 /* zero out the visited states so we know how many equilibrated states we have
379 dfhist->n_at_lam[i] = 0;
384 /* If we reached this far, we have not equilibrated yet, keep on
385 going resetting the weights */
387 if (EWL(expand->elamstats))
389 if (expand->elamstats==elamstatsWL) /* Standard Wang-Landau */
391 dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
392 dfhist->wl_histo[fep_state] += 1.0;
394 else if (expand->elamstats==elamstatsWWL) /* Weighted Wang-Landau */
398 /* first increment count */
399 GenerateGibbsProbabilities(weighted_lamee,p_k,&pks,0,nlim-1);
400 for (i=0;i<nlim;i++) {
401 dfhist->wl_histo[i] += p_k[i];
404 /* then increment weights (uses count) */
406 GenerateWeightedGibbsProbabilities(weighted_lamee,p_k,&pks,nlim,dfhist->wl_histo,dfhist->wl_delta);
410 dfhist->sum_weights[i] -= dfhist->wl_delta*p_k[i];
412 /* Alternate definition, using logarithms. Shouldn't make very much difference! */
417 di = 1+dfhist->wl_delta*p_k[i];
418 dfhist->sum_weights[i] -= log(di);
424 zero_sum_weights = dfhist->sum_weights[0];
427 dfhist->sum_weights[i] -= zero_sum_weights;
431 if (expand->elamstats==elamstatsBARKER || expand->elamstats==elamstatsMETROPOLIS || expand->elamstats==elamstatsMINVAR) {
433 de_function = 0; /* to get rid of warnings, but this value will not be used because of the logic */
434 maxc = 2*expand->c_range+1;
437 snew(lam_variance,nlim);
439 snew(omegap_array,maxc);
440 snew(weightsp_array,maxc);
441 snew(varp_array,maxc);
442 snew(dwp_array,maxc);
444 snew(omegam_array,maxc);
445 snew(weightsm_array,maxc);
446 snew(varm_array,maxc);
447 snew(dwm_array,maxc);
449 /* unpack the current lambdas -- we will only update 2 of these */
451 for (i=0;i<nlim-1;i++)
452 { /* only through the second to last */
453 lam_dg[i] = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
454 lam_variance[i] = pow(dfhist->sum_variance[i+1],2) - pow(dfhist->sum_variance[i],2);
457 /* accumulate running averages */
458 for (nval = 0; nval<maxc; nval++)
460 /* constants for later use */
461 cnval = (real)(nval-expand->c_range);
462 /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
465 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
466 if (expand->elamstats==elamstatsBARKER || expand->elamstats==elamstatsMINVAR)
468 de_function = 1.0/(1.0+de);
470 else if (expand->elamstats==elamstatsMETROPOLIS)
478 de_function = 1.0/de;
481 dfhist->accum_m[fep_state][nval] += de_function;
482 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
485 if (fep_state < nlim-1)
487 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
488 if (expand->elamstats==elamstatsBARKER || expand->elamstats==elamstatsMINVAR)
490 de_function = 1.0/(1.0+de);
492 else if (expand->elamstats==elamstatsMETROPOLIS)
500 de_function = 1.0/de;
503 dfhist->accum_p[fep_state][nval] += de_function;
504 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
507 /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
509 n0 = dfhist->n_at_lam[fep_state];
510 if (fep_state > 0) {nm1 = dfhist->n_at_lam[fep_state-1];} else {nm1 = 0;}
511 if (fep_state < nlim-1) {np1 = dfhist->n_at_lam[fep_state+1];} else {np1 = 0;}
513 /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
514 chi_m1_0=chi_p1_0=chi_m2_0=chi_p2_0=chi_p1_m1=chi_p2_m1=chi_m1_p1=chi_m2_p1=0;
518 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
519 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
520 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
521 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
524 if ((fep_state > 0 ) && (nm1 > 0))
526 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
527 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
530 if ((fep_state < nlim-1) && (np1 > 0))
532 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
533 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
547 omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
551 omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
553 if ((n0 > 0) && (nm1 > 0))
555 clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
556 clam_varm = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
560 if (fep_state < nlim-1)
564 omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
568 omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
570 if ((n0 > 0) && (np1 > 0))
572 clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
573 clam_varp = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
579 omegam_array[nval] = omega_m1_0;
583 omegam_array[nval] = 0;
585 weightsm_array[nval] = clam_weightsm;
586 varm_array[nval] = clam_varm;
589 dwm_array[nval] = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
593 dwm_array[nval] = fabs( cnval - lam_dg[fep_state-1] );
598 omegap_array[nval] = omega_p1_0;
602 omegap_array[nval] = 0;
604 weightsp_array[nval] = clam_weightsp;
605 varp_array[nval] = clam_varp;
606 if ((np1 > 0) && (n0 > 0))
608 dwp_array[nval] = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
612 dwp_array[nval] = fabs( cnval - lam_dg[fep_state] );
617 /* find the C's closest to the old weights value */
619 min_nvalm = FindMinimum(dwm_array,maxc);
620 omega_m1_0 = omegam_array[min_nvalm];
621 clam_weightsm = weightsm_array[min_nvalm];
622 clam_varm = varm_array[min_nvalm];
624 min_nvalp = FindMinimum(dwp_array,maxc);
625 omega_p1_0 = omegap_array[min_nvalp];
626 clam_weightsp = weightsp_array[min_nvalp];
627 clam_varp = varp_array[min_nvalp];
629 clam_osum = omega_m1_0 + omega_p1_0;
633 clam_minvar = 0.5*log(clam_osum);
638 lam_dg[fep_state-1] = clam_weightsm;
639 lam_variance[fep_state-1] = clam_varm;
642 if (fep_state < nlim-1)
644 lam_dg[fep_state] = clam_weightsp;
645 lam_variance[fep_state] = clam_varp;
648 if (expand->elamstats==elamstatsMINVAR)
650 bSufficientSamples = TRUE;
651 /* make sure they are all past a threshold */
654 if (dfhist->n_at_lam[i] < expand->minvarmin)
656 bSufficientSamples = FALSE;
659 if (bSufficientSamples)
661 dfhist->sum_minvar[fep_state] = clam_minvar;
666 dfhist->sum_minvar[i]+=(expand->minvar_const-clam_minvar);
668 expand->minvar_const = clam_minvar;
669 dfhist->sum_minvar[fep_state] = 0.0;
673 dfhist->sum_minvar[fep_state] -= expand->minvar_const;
678 /* we need to rezero minvar now, since it could change at fep_state = 0 */
679 dfhist->sum_dg[0] = 0.0;
680 dfhist->sum_variance[0] = 0.0;
681 dfhist->sum_weights[0] = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
685 dfhist->sum_dg[i] = lam_dg[i-1] + dfhist->sum_dg[i-1];
686 dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1],2));
687 dfhist->sum_weights[i] = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
694 sfree(weightsm_array);
699 sfree(weightsp_array);
706 static int ChooseNewLambda(FILE *log, int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, real *p_k, gmx_rng_t rng)
708 /* Choose new lambda value, and update transition matrix */
710 int i,ifep,jfep,minfep,maxfep,lamnew,lamtrial,starting_fep_state;
711 real r1,r2,pks,de_old,de_new,de,trialprob,tprob=0;
713 real *propose,*accept,*remainder;
715 gmx_bool bRestricted;
717 starting_fep_state = fep_state;
718 lamnew = fep_state; /* so that there is a default setting -- stays the same */
720 if (!EWL(expand->elamstats)) /* ignore equilibrating the weights if using WL */
722 if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
724 /* Use a marching method to run through the lambdas and get preliminary free energy data,
725 before starting 'free' sampling. We start free sampling when we have enough at each lambda */
727 /* if we have enough at this lambda, move on to the next one */
729 if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
731 lamnew = fep_state+1;
732 if (lamnew == nlim) /* whoops, stepped too far! */
747 snew(remainder,nlim);
749 for (i=0;i<expand->lmc_repeats;i++)
752 for(ifep=0;ifep<nlim;ifep++)
758 if ((expand->elmcmove==elmcmoveGIBBS) || (expand->elmcmove==elmcmoveMETGIBBS))
761 /* use the Gibbs sampler, with restricted range */
762 if (expand->gibbsdeltalam < 0)
770 minfep = fep_state - expand->gibbsdeltalam;
771 maxfep = fep_state + expand->gibbsdeltalam;
782 GenerateGibbsProbabilities(weighted_lamee,p_k,&pks,minfep,maxfep);
784 if (expand->elmcmove == elmcmoveGIBBS)
786 for (ifep=minfep;ifep<=maxfep;ifep++)
788 propose[ifep] = p_k[ifep];
792 r1 = gmx_rng_uniform_real(rng);
793 for (lamnew=minfep;lamnew<=maxfep;lamnew++)
795 if (r1 <= p_k[lamnew])
802 else if (expand->elmcmove==elmcmoveMETGIBBS)
805 /* Metropolized Gibbs sampling */
806 for (ifep=minfep;ifep<=maxfep;ifep++)
808 remainder[ifep] = 1 - p_k[ifep];
811 /* find the proposal probabilities */
813 if (remainder[fep_state] == 0) {
814 /* only the current state has any probability */
815 /* we have to stay at the current state */
818 for (ifep=minfep;ifep<=maxfep;ifep++)
820 if (ifep != fep_state)
822 propose[ifep] = p_k[ifep]/remainder[fep_state];
830 r1 = gmx_rng_uniform_real(rng);
831 for (lamtrial=minfep;lamtrial<=maxfep;lamtrial++)
833 pnorm = p_k[lamtrial]/remainder[fep_state];
834 if (lamtrial!=fep_state)
844 /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
846 /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
847 trialprob = (remainder[fep_state])/(remainder[lamtrial]);
848 if (trialprob < tprob)
852 r2 = gmx_rng_uniform_real(rng);
863 /* now figure out the acceptance probability for each */
864 for (ifep=minfep;ifep<=maxfep;ifep++)
867 if (remainder[ifep] != 0) {
868 trialprob = (remainder[fep_state])/(remainder[ifep]);
872 trialprob = 1.0; /* this state is the only choice! */
874 if (trialprob < tprob)
878 /* probability for fep_state=0, but that's fine, it's never proposed! */
879 accept[ifep] = tprob;
885 /* it's possible some rounding is failing */
886 if (remainder[fep_state] < 2.0e-15)
888 /* probably numerical rounding error -- no state other than the original has weight */
893 /* probably not a numerical issue */
895 int nerror = 200+(maxfep-minfep+1)*60;
897 snew(errorstr,nerror);
898 /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
899 of sum weights. Generated detailed info for failure */
900 loc += sprintf(errorstr,"Something wrong in choosing new lambda state with a Gibbs move -- probably underflow in weight determination.\nDenominator is: %3d%17.10e\n i dE numerator weights\n",0,pks);
901 for (ifep=minfep;ifep<=maxfep;ifep++)
903 loc += sprintf(&errorstr[loc],"%3d %17.10e%17.10e%17.10e\n",ifep,weighted_lamee[ifep],p_k[ifep],dfhist->sum_weights[ifep]);
905 gmx_fatal(FARGS,errorstr);
909 else if ((expand->elmcmove==elmcmoveMETROPOLIS) || (expand->elmcmove==elmcmoveBARKER))
911 /* use the metropolis sampler with trial +/- 1 */
912 r1 = gmx_rng_uniform_real(rng);
915 if (fep_state == 0) {
916 lamtrial = fep_state;
920 lamtrial = fep_state-1;
925 if (fep_state == nlim-1) {
926 lamtrial = fep_state;
930 lamtrial = fep_state+1;
934 de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
935 if (expand->elmcmove==elmcmoveMETROPOLIS)
939 if (trialprob < tprob)
943 propose[fep_state] = 0;
944 propose[lamtrial] = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
945 accept[fep_state] = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
946 accept[lamtrial] = tprob;
949 else if (expand->elmcmove==elmcmoveBARKER)
951 tprob = 1.0/(1.0+exp(-de));
953 propose[fep_state] = (1-tprob);
954 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
955 accept[fep_state] = 1.0;
956 accept[lamtrial] = 1.0;
959 r2 = gmx_rng_uniform_real(rng);
967 for (ifep=0;ifep<nlim;ifep++)
969 dfhist->Tij[fep_state][ifep] += propose[ifep]*accept[ifep];
970 dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
975 dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
984 /* print out the weights to the log, along with current state */
985 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
986 int nlam, int frequency, gmx_large_int_t step)
988 int nlim,i,ifep,jfep;
989 real dw,dg,dv,dm,Tprint;
991 const char *print_names[efptNR] = {" FEPL","MassL","CoulL"," VdwL","BondL","RestT","Temp.(K)"};
992 gmx_bool bSimTemp = FALSE;
994 nlim = fep->n_lambda;
995 if (simtemp != NULL) {
999 if (mod(step,frequency)==0)
1001 fprintf(outfile," MC-lambda information\n");
1002 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) {
1003 fprintf(outfile," Wang-Landau incrementor is: %11.5g\n",dfhist->wl_delta);
1005 fprintf(outfile," N");
1006 for (i=0;i<efptNR;i++)
1008 if (fep->separate_dvdl[i])
1010 fprintf(outfile,"%7s",print_names[i]);
1012 else if ((i == efptTEMPERATURE) && bSimTemp)
1014 fprintf(outfile,"%10s",print_names[i]); /* more space for temperature formats */
1017 fprintf(outfile," Count ");
1018 if (expand->elamstats==elamstatsMINVAR)
1020 fprintf(outfile,"W(in kT) G(in kT) dG(in kT) dV(in kT)\n");
1024 fprintf(outfile,"G(in kT) dG(in kT)\n");
1026 for (ifep=0;ifep<nlim;ifep++)
1037 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1038 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1039 dv = sqrt(pow(dfhist->sum_variance[ifep+1],2) - pow(dfhist->sum_variance[ifep],2));
1040 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1043 fprintf(outfile,"%3d",(ifep+1));
1044 for (i=0;i<efptNR;i++)
1046 if (fep->separate_dvdl[i])
1048 fprintf(outfile,"%7.3f",fep->all_lambda[i][ifep]);
1049 } else if (i == efptTEMPERATURE && bSimTemp)
1051 fprintf(outfile,"%9.3f",simtemp->temperatures[ifep]);
1054 if (EWL(expand->elamstats) && (!(dfhist->bEquil))) /* if performing WL and still haven't equilibrated */
1056 if (expand->elamstats == elamstatsWL)
1058 fprintf(outfile," %8d",(int)dfhist->wl_histo[ifep]);
1060 fprintf(outfile," %8.3f",dfhist->wl_histo[ifep]);
1063 else /* we have equilibrated weights */
1065 fprintf(outfile," %8d",dfhist->n_at_lam[ifep]);
1067 if (expand->elamstats==elamstatsMINVAR)
1069 fprintf(outfile," %10.5f %10.5f %10.5f %10.5f",dfhist->sum_weights[ifep],dfhist->sum_dg[ifep],dg,dv);
1073 fprintf(outfile," %10.5f %10.5f",dfhist->sum_weights[ifep],dw);
1076 fprintf(outfile," <<\n");
1080 fprintf(outfile," \n");
1083 fprintf(outfile,"\n");
1085 if ((mod(step,expand->nstTij)==0) && (expand->nstTij > 0) && (step > 0))
1087 fprintf(outfile," Transition Matrix\n");
1088 for (ifep=0;ifep<nlim;ifep++)
1090 fprintf(outfile,"%12d",(ifep+1));
1092 fprintf(outfile,"\n");
1093 for (ifep=0;ifep<nlim;ifep++)
1095 for (jfep=0;jfep<nlim;jfep++)
1097 if (dfhist->n_at_lam[ifep] > 0)
1099 if (expand->bSymmetrizedTMatrix)
1101 Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1103 Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1110 fprintf(outfile,"%12.8f",Tprint);
1112 fprintf(outfile,"%3d\n",(ifep+1));
1115 fprintf(outfile," Empirical Transition Matrix\n");
1116 for (ifep=0;ifep<nlim;ifep++)
1118 fprintf(outfile,"%12d",(ifep+1));
1120 fprintf(outfile,"\n");
1121 for (ifep=0;ifep<nlim;ifep++)
1123 for (jfep=0;jfep<nlim;jfep++)
1125 if (dfhist->n_at_lam[ifep] > 0)
1127 if (expand->bSymmetrizedTMatrix)
1129 Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1131 Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1138 fprintf(outfile,"%12.8f",Tprint);
1140 fprintf(outfile,"%3d\n",(ifep+1));
1146 extern void get_mc_state(gmx_rng_t rng,t_state *state)
1148 gmx_rng_get_state(rng,state->mc_rng,state->mc_rngi);
1151 extern void set_mc_state(gmx_rng_t rng,t_state *state)
1153 gmx_rng_set_state(rng,state->mc_rng,state->mc_rngi[0]);
1156 extern int ExpandedEnsembleDynamics(FILE *log,t_inputrec *ir, gmx_enerdata_t *enerd,
1157 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
1158 gmx_large_int_t step, gmx_rng_t mcrng,
1159 rvec *v, t_mdatoms *mdatoms)
1161 real *pfep_lamee,*p_k, *scaled_lamee, *weighted_lamee;
1162 int i,nlam,nlim,lamnew,totalsamples;
1163 real oneovert,maxscaled=0,maxweighted=0;
1166 double *temperature_lambdas;
1167 gmx_bool bIfReset,bSwitchtoOneOverT,bDoneEquilibrating=FALSE;
1169 expand = ir->expandedvals;
1170 simtemp = ir->simtempvals;
1171 nlim = ir->fepvals->n_lambda;
1172 nlam = state->fep_state;
1174 snew(scaled_lamee,nlim);
1175 snew(weighted_lamee,nlim);
1176 snew(pfep_lamee,nlim);
1179 if (expand->bInit_weights) /* if initialized weights, we need to fill them in */
1181 dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
1182 for (i=0;i<nlim;i++) {
1183 dfhist->sum_weights[i] = expand->init_lambda_weights[i];
1184 dfhist->sum_dg[i] = expand->init_lambda_weights[i];
1186 expand->bInit_weights = FALSE;
1189 /* update the count at the current lambda*/
1190 dfhist->n_at_lam[nlam]++;
1192 /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1193 pressure controlled.*/
1196 where does this PV term go?
1197 for (i=0;i<nlim;i++)
1199 fep_lamee[i] += pVTerm;
1203 /* determine the minimum value to avoid overflow. Probably a better way to do this */
1204 /* we don't need to include the pressure term, since the volume is the same between the two.
1205 is there some term we are neglecting, however? */
1207 if (ir->efep != efepNO)
1209 for (i=0;i<nlim;i++)
1213 /* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
1214 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1215 + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
1219 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1220 /* mc_temp is currently set to the system reft unless otherwise defined */
1223 /* save these energies for printing, so they don't get overwritten by the next step */
1224 /* they aren't overwritten in the non-free energy case, but we always print with these
1229 for (i=0;i<nlim;i++) {
1230 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
1235 for (i=0;i<nlim;i++) {
1236 pfep_lamee[i] = scaled_lamee[i];
1238 weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1241 maxscaled = scaled_lamee[i];
1242 maxweighted = weighted_lamee[i];
1246 if (scaled_lamee[i] > maxscaled)
1248 maxscaled = scaled_lamee[i];
1250 if (weighted_lamee[i] > maxweighted)
1252 maxweighted = weighted_lamee[i];
1257 for (i=0;i<nlim;i++)
1259 scaled_lamee[i] -= maxscaled;
1260 weighted_lamee[i] -= maxweighted;
1263 /* update weights - we decide whether or not to actually do this inside */
1265 bDoneEquilibrating = UpdateWeights(nlim,expand,dfhist,nlam,scaled_lamee,weighted_lamee,step);
1266 if (bDoneEquilibrating)
1269 fprintf(log,"\nStep %d: Weights have equilibrated, using criteria: %s\n",(int)step,elmceq_names[expand->elmceq]);
1273 lamnew = ChooseNewLambda(log,nlim,expand,dfhist,nlam,weighted_lamee,p_k,mcrng);
1274 /* if using simulated tempering, we need to adjust the temperatures */
1275 if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
1280 int nstart, nend, gt;
1282 snew(buf_ngtc,ir->opts.ngtc);
1284 for (i=0;i<ir->opts.ngtc;i++) {
1285 if (ir->opts.ref_t[i] > 0) {
1286 told = ir->opts.ref_t[i];
1287 ir->opts.ref_t[i] = simtemp->temperatures[lamnew];
1288 buf_ngtc[i] = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1292 /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1294 nstart = mdatoms->start;
1295 nend = nstart + mdatoms->homenr;
1296 for (n=nstart; n<nend; n++)
1301 gt = mdatoms->cTC[n];
1303 for(d=0; d<DIM; d++)
1305 v[n][d] *= buf_ngtc[gt];
1309 if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir)) {
1310 /* we need to recalculate the masses if the temperature has changed */
1311 init_npt_masses(ir,state,MassQ,FALSE);
1312 for (i=0;i<state->nnhpres;i++)
1314 for (j=0;j<ir->opts.nhchainlength;j++)
1316 state->nhpres_vxi[i+j] *= buf_ngtc[i];
1319 for (i=0;i<ir->opts.ngtc;i++)
1321 for (j=0;j<ir->opts.nhchainlength;j++)
1323 state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1330 /* now check on the Wang-Landau updating critera */
1332 if (EWL(expand->elamstats))
1334 bSwitchtoOneOverT = FALSE;
1335 if (expand->bWLoneovert) {
1337 for (i=0;i<nlim;i++)
1339 totalsamples += dfhist->n_at_lam[i];
1341 oneovert = (1.0*nlim)/totalsamples;
1342 /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1343 /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1344 if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1345 (dfhist->wl_delta < expand->init_wl_delta))
1347 bSwitchtoOneOverT = TRUE;
1350 if (bSwitchtoOneOverT) {
1351 dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1353 bIfReset = CheckHistogramRatios(nlim,dfhist->wl_histo,expand->wl_ratio);
1356 for (i=0;i<nlim;i++)
1358 dfhist->wl_histo[i] = 0;
1360 dfhist->wl_delta *= expand->wl_scale;
1362 fprintf(log,"\nStep %d: weights are now:",(int)step);
1363 for (i=0;i<nlim;i++)
1365 fprintf(log," %.5f",dfhist->sum_weights[i]);
1372 sfree(scaled_lamee);
1373 sfree(weighted_lamee);