Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / mdlib / expanded.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
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.
14  *
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROwing Monsters And Cloning Shrimps
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #ifdef GMX_CRAY_XT3
40 #include <catamount/dclock.h>
41 #endif
42
43
44 #include <stdio.h>
45 #include <time.h>
46 #ifdef HAVE_SYS_TIME_H
47 #include <sys/time.h>
48 #endif
49 #include <math.h>
50 #include "typedefs.h"
51 #include "string2.h"
52 #include "gmxfio.h"
53 #include "smalloc.h"
54 #include "names.h"
55 #include "confio.h"
56 #include "mvdata.h"
57 #include "txtdump.h"
58 #include "pbc.h"
59 #include "chargegroup.h"
60 #include "vec.h"
61 #include "nrnb.h"
62 #include "mshift.h"
63 #include "mdrun.h"
64 #include "update.h"
65 #include "physics.h"
66 #include "main.h"
67 #include "mdatoms.h"
68 #include "force.h"
69 #include "bondf.h"
70 #include "pme.h"
71 #include "disre.h"
72 #include "orires.h"
73 #include "network.h"
74 #include "calcmu.h"
75 #include "constr.h"
76 #include "xvgr.h"
77 #include "trnio.h"
78 #include "xtcio.h"
79 #include "gmx_random.h"
80 #include "domdec.h"
81 #include "partdec.h"
82 #include "gmx_wallcycle.h"
83 #include "macros.h"
84
85 #include "gromacs/utility/gmxmpi.h"
86
87 void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
88 {
89
90     int  i;
91     real maxene;
92
93     *pks   = 0.0;
94     maxene = ene[minfep];
95     /* find the maximum value */
96     for (i = minfep; i <= maxfep; i++)
97     {
98         if (ene[i] > maxene)
99         {
100             maxene = ene[i];
101         }
102     }
103     /* find the denominator */
104     for (i = minfep; i <= maxfep; i++)
105     {
106         *pks += exp(ene[i]-maxene);
107     }
108     /*numerators*/
109     for (i = minfep; i <= maxfep; i++)
110     {
111         p_k[i] = exp(ene[i]-maxene) / *pks;
112     }
113 }
114
115 void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
116 {
117
118     int   i;
119     real  maxene;
120     real *nene;
121     *pks = 0.0;
122
123     snew(nene, nlim);
124     for (i = 0; i < nlim; i++)
125     {
126         if (nvals[i] == 0)
127         {
128             /* add the delta, since we need to make sure it's greater than zero, and
129                we need a non-arbitrary number? */
130             nene[i] = ene[i] + log(nvals[i]+delta);
131         }
132         else
133         {
134             nene[i] = ene[i] + log(nvals[i]);
135         }
136     }
137
138     /* find the maximum value */
139     maxene = nene[0];
140     for (i = 0; i < nlim; i++)
141     {
142         if (nene[i] > maxene)
143         {
144             maxene = nene[i];
145         }
146     }
147
148     /* subtract off the maximum, avoiding overflow */
149     for (i = 0; i < nlim; i++)
150     {
151         nene[i] -= maxene;
152     }
153
154     /* find the denominator */
155     for (i = 0; i < nlim; i++)
156     {
157         *pks += exp(nene[i]);
158     }
159
160     /*numerators*/
161     for (i = 0; i < nlim; i++)
162     {
163         p_k[i] = exp(nene[i]) / *pks;
164     }
165     sfree(nene);
166 }
167
168 real do_logsum(int N, real *a_n)
169 {
170
171     /*     RETURN VALUE */
172     /* log(\sum_{i=0}^(N-1) exp[a_n]) */
173     real maxarg;
174     real sum;
175     int  i;
176     real logsum;
177     /*     compute maximum argument to exp(.) */
178
179     maxarg = a_n[0];
180     for (i = 1; i < N; i++)
181     {
182         maxarg = max(maxarg, a_n[i]);
183     }
184
185     /* compute sum of exp(a_n - maxarg) */
186     sum = 0.0;
187     for (i = 0; i < N; i++)
188     {
189         sum = sum + exp(a_n[i] - maxarg);
190     }
191
192     /*     compute log sum */
193     logsum = log(sum) + maxarg;
194     return logsum;
195 }
196
197 int FindMinimum(real *min_metric, int N)
198 {
199
200     real min_val;
201     int  min_nval, nval;
202
203     min_nval = 0;
204     min_val  = min_metric[0];
205
206     for (nval = 0; nval < N; nval++)
207     {
208         if (min_metric[nval] < min_val)
209         {
210             min_val  = min_metric[nval];
211             min_nval = nval;
212         }
213     }
214     return min_nval;
215 }
216
217 static gmx_bool CheckHistogramRatios(int nhisto, real *histo, real ratio)
218 {
219
220     int      i;
221     real     nmean;
222     gmx_bool bIfFlat;
223
224     nmean = 0;
225     for (i = 0; i < nhisto; i++)
226     {
227         nmean += histo[i];
228     }
229
230     if (nmean == 0)
231     {
232         /* no samples! is bad!*/
233         bIfFlat = FALSE;
234         return bIfFlat;
235     }
236     nmean /= (real)nhisto;
237
238     bIfFlat = TRUE;
239     for (i = 0; i < nhisto; i++)
240     {
241         /* make sure that all points are in the ratio < x <  1/ratio range  */
242         if (!((histo[i]/nmean < 1.0/ratio) && (histo[i]/nmean > ratio)))
243         {
244             bIfFlat = FALSE;
245             break;
246         }
247     }
248     return bIfFlat;
249 }
250
251 static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_large_int_t step)
252 {
253
254     int      i, totalsamples;
255     gmx_bool bDoneEquilibrating = TRUE;
256     gmx_bool bIfFlat;
257
258     /* assume we have equilibrated the weights, then check to see if any of the conditions are not met */
259
260     /* calculate the total number of samples */
261     switch (expand->elmceq)
262     {
263         case elmceqNO:
264             /* We have not equilibrated, and won't, ever. */
265             return FALSE;
266         case elmceqYES:
267             /* we have equilibrated -- we're done */
268             return TRUE;
269         case elmceqSTEPS:
270             /* first, check if we are equilibrating by steps, if we're still under */
271             if (step < expand->equil_steps)
272             {
273                 bDoneEquilibrating = FALSE;
274             }
275             break;
276         case elmceqSAMPLES:
277             totalsamples = 0;
278             for (i = 0; i < nlim; i++)
279             {
280                 totalsamples += dfhist->n_at_lam[i];
281             }
282             if (totalsamples < expand->equil_samples)
283             {
284                 bDoneEquilibrating = FALSE;
285             }
286             break;
287         case elmceqNUMATLAM:
288             for (i = 0; i < nlim; i++)
289             {
290                 if (dfhist->n_at_lam[i] < expand->equil_n_at_lam) /* we are still doing the initial sweep, so we're definitely not
291                                                                      done equilibrating*/
292                 {
293                     bDoneEquilibrating  = FALSE;
294                     break;
295                 }
296             }
297             break;
298         case elmceqWLDELTA:
299             if (EWL(expand->elamstats)) /* This check is in readir as well, but
300                                            just to be sure */
301             {
302                 if (dfhist->wl_delta > expand->equil_wl_delta)
303                 {
304                     bDoneEquilibrating = FALSE;
305                 }
306             }
307             break;
308         case elmceqRATIO:
309             /* we can use the flatness as a judge of good weights, as long as
310                we're not doing minvar, or Wang-Landau.
311                But turn off for now until we figure out exactly how we do this.
312              */
313
314             if (!(EWL(expand->elamstats) || expand->elamstats == elamstatsMINVAR))
315             {
316                 /* we want to use flatness -avoiding- the forced-through samples.  Plus, we need to convert to
317                    floats for this histogram function. */
318
319                 real *modhisto;
320                 snew(modhisto, nlim);
321                 for (i = 0; i < nlim; i++)
322                 {
323                     modhisto[i] = 1.0*(dfhist->n_at_lam[i]-expand->lmc_forced_nstart);
324                 }
325                 bIfFlat = CheckHistogramRatios(nlim, modhisto, expand->equil_ratio);
326                 sfree(modhisto);
327                 if (!bIfFlat)
328                 {
329                     bDoneEquilibrating = FALSE;
330                 }
331             }
332         default:
333             bDoneEquilibrating = TRUE;
334     }
335     /* one last case to go though, if we are doing slow growth to get initial values, we haven't finished equilibrating */
336
337     if (expand->lmc_forced_nstart > 0)
338     {
339         for (i = 0; i < nlim; i++)
340         {
341             if (dfhist->n_at_lam[i] < expand->lmc_forced_nstart) /* we are still doing the initial sweep, so we're definitely not
342                                                                     done equilibrating*/
343             {
344                 bDoneEquilibrating = FALSE;
345                 break;
346             }
347         }
348     }
349     return bDoneEquilibrating;
350 }
351
352 static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
353                               int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_large_int_t step)
354 {
355     real     maxdiff = 0.000000001;
356     gmx_bool bSufficientSamples;
357     int      i, k, n, nz, indexi, indexk, min_n, max_n, nlam, totali;
358     int      n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
359     real     omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
360     real     de, de_function, dr, denom, maxdr;
361     real     min_val, cnval, zero_sum_weights;
362     real    *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
363     real     clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
364     real    *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
365     double  *p_k;
366     double  pks = 0;
367     real    *numweighted_lamee, *logfrac;
368     int     *nonzero;
369     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;
370
371     /* if we have equilibrated the weights, exit now */
372     if (dfhist->bEquil)
373     {
374         return FALSE;
375     }
376
377     if (CheckIfDoneEquilibrating(nlim, expand, dfhist, step))
378     {
379         dfhist->bEquil = TRUE;
380         /* zero out the visited states so we know how many equilibrated states we have
381            from here on out.*/
382         for (i = 0; i < nlim; i++)
383         {
384             dfhist->n_at_lam[i] = 0;
385         }
386         return TRUE;
387     }
388
389     /* If we reached this far, we have not equilibrated yet, keep on
390        going resetting the weights */
391
392     if (EWL(expand->elamstats))
393     {
394         if (expand->elamstats == elamstatsWL)  /* Standard Wang-Landau */
395         {
396             dfhist->sum_weights[fep_state] -= dfhist->wl_delta;
397             dfhist->wl_histo[fep_state]    += 1.0;
398         }
399         else if (expand->elamstats == elamstatsWWL) /* Weighted Wang-Landau */
400         {
401             snew(p_k, nlim);
402
403             /* first increment count */
404             GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
405             for (i = 0; i < nlim; i++)
406             {
407                 dfhist->wl_histo[i] += (real)p_k[i];
408             }
409
410             /* then increment weights (uses count) */
411             pks = 0.0;
412             GenerateWeightedGibbsProbabilities(weighted_lamee, p_k, &pks, nlim, dfhist->wl_histo, dfhist->wl_delta);
413
414             for (i = 0; i < nlim; i++)
415             {
416                 dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
417             }
418             /* Alternate definition, using logarithms. Shouldn't make very much difference! */
419             /*
420                real di;
421                for (i=0;i<nlim;i++)
422                {
423                 di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
424                 dfhist->sum_weights[i] -= log(di);
425                }
426              */
427             sfree(p_k);
428         }
429
430         zero_sum_weights =  dfhist->sum_weights[0];
431         for (i = 0; i < nlim; i++)
432         {
433             dfhist->sum_weights[i] -= zero_sum_weights;
434         }
435     }
436
437     if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMETROPOLIS || expand->elamstats == elamstatsMINVAR)
438     {
439
440         de_function = 0;  /* to get rid of warnings, but this value will not be used because of the logic */
441         maxc        = 2*expand->c_range+1;
442
443         snew(lam_dg, nlim);
444         snew(lam_variance, nlim);
445
446         snew(omegap_array, maxc);
447         snew(weightsp_array, maxc);
448         snew(varp_array, maxc);
449         snew(dwp_array, maxc);
450
451         snew(omegam_array, maxc);
452         snew(weightsm_array, maxc);
453         snew(varm_array, maxc);
454         snew(dwm_array, maxc);
455
456         /* unpack the current lambdas -- we will only update 2 of these */
457
458         for (i = 0; i < nlim-1; i++)
459         {   /* only through the second to last */
460             lam_dg[i]       = dfhist->sum_dg[i+1] - dfhist->sum_dg[i];
461             lam_variance[i] = pow(dfhist->sum_variance[i+1], 2) - pow(dfhist->sum_variance[i], 2);
462         }
463
464         /* accumulate running averages */
465         for (nval = 0; nval < maxc; nval++)
466         {
467             /* constants for later use */
468             cnval = (real)(nval-expand->c_range);
469             /* actually, should be able to rewrite it w/o exponential, for better numerical stability */
470             if (fep_state > 0)
471             {
472                 de = exp(cnval - (scaled_lamee[fep_state]-scaled_lamee[fep_state-1]));
473                 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
474                 {
475                     de_function = 1.0/(1.0+de);
476                 }
477                 else if (expand->elamstats == elamstatsMETROPOLIS)
478                 {
479                     if (de < 1.0)
480                     {
481                         de_function = 1.0;
482                     }
483                     else
484                     {
485                         de_function = 1.0/de;
486                     }
487                 }
488                 dfhist->accum_m[fep_state][nval]  += de_function;
489                 dfhist->accum_m2[fep_state][nval] += de_function*de_function;
490             }
491
492             if (fep_state < nlim-1)
493             {
494                 de = exp(-cnval + (scaled_lamee[fep_state+1]-scaled_lamee[fep_state]));
495                 if (expand->elamstats == elamstatsBARKER || expand->elamstats == elamstatsMINVAR)
496                 {
497                     de_function = 1.0/(1.0+de);
498                 }
499                 else if (expand->elamstats == elamstatsMETROPOLIS)
500                 {
501                     if (de < 1.0)
502                     {
503                         de_function = 1.0;
504                     }
505                     else
506                     {
507                         de_function = 1.0/de;
508                     }
509                 }
510                 dfhist->accum_p[fep_state][nval]  += de_function;
511                 dfhist->accum_p2[fep_state][nval] += de_function*de_function;
512             }
513
514             /* Metropolis transition and Barker transition (unoptimized Bennett) acceptance weight determination */
515
516             n0  = dfhist->n_at_lam[fep_state];
517             if (fep_state > 0)
518             {
519                 nm1 = dfhist->n_at_lam[fep_state-1];
520             }
521             else
522             {
523                 nm1 = 0;
524             }
525             if (fep_state < nlim-1)
526             {
527                 np1 = dfhist->n_at_lam[fep_state+1];
528             }
529             else
530             {
531                 np1 = 0;
532             }
533
534             /* logic SHOULD keep these all set correctly whatever the logic, but apparently it can't figure it out. */
535             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;
536
537             if (n0 > 0)
538             {
539                 chi_m1_0 = dfhist->accum_m[fep_state][nval]/n0;
540                 chi_p1_0 = dfhist->accum_p[fep_state][nval]/n0;
541                 chi_m2_0 = dfhist->accum_m2[fep_state][nval]/n0;
542                 chi_p2_0 = dfhist->accum_p2[fep_state][nval]/n0;
543             }
544
545             if ((fep_state > 0 ) && (nm1 > 0))
546             {
547                 chi_p1_m1 = dfhist->accum_p[fep_state-1][nval]/nm1;
548                 chi_p2_m1 = dfhist->accum_p2[fep_state-1][nval]/nm1;
549             }
550
551             if ((fep_state < nlim-1) && (np1 > 0))
552             {
553                 chi_m1_p1 = dfhist->accum_m[fep_state+1][nval]/np1;
554                 chi_m2_p1 = dfhist->accum_m2[fep_state+1][nval]/np1;
555             }
556
557             omega_m1_0    = 0;
558             omega_p1_0    = 0;
559             clam_weightsm = 0;
560             clam_weightsp = 0;
561             clam_varm     = 0;
562             clam_varp     = 0;
563
564             if (fep_state > 0)
565             {
566                 if (n0 > 0)
567                 {
568                     omega_m1_0 = chi_m2_0/(chi_m1_0*chi_m1_0) - 1.0;
569                 }
570                 if (nm1 > 0)
571                 {
572                     omega_p1_m1 = chi_p2_m1/(chi_p1_m1*chi_p1_m1) - 1.0;
573                 }
574                 if ((n0 > 0) && (nm1 > 0))
575                 {
576                     clam_weightsm = (log(chi_m1_0) - log(chi_p1_m1)) + cnval;
577                     clam_varm     = (1.0/n0)*(omega_m1_0) + (1.0/nm1)*(omega_p1_m1);
578                 }
579             }
580
581             if (fep_state < nlim-1)
582             {
583                 if (n0 > 0)
584                 {
585                     omega_p1_0 = chi_p2_0/(chi_p1_0*chi_p1_0) - 1.0;
586                 }
587                 if (np1 > 0)
588                 {
589                     omega_m1_p1 = chi_m2_p1/(chi_m1_p1*chi_m1_p1) - 1.0;
590                 }
591                 if ((n0 > 0) && (np1 > 0))
592                 {
593                     clam_weightsp = (log(chi_m1_p1) - log(chi_p1_0)) + cnval;
594                     clam_varp     = (1.0/np1)*(omega_m1_p1) + (1.0/n0)*(omega_p1_0);
595                 }
596             }
597
598             if (n0 > 0)
599             {
600                 omegam_array[nval]             = omega_m1_0;
601             }
602             else
603             {
604                 omegam_array[nval]             = 0;
605             }
606             weightsm_array[nval]           = clam_weightsm;
607             varm_array[nval]               = clam_varm;
608             if (nm1 > 0)
609             {
610                 dwm_array[nval]  = fabs( (cnval + log((1.0*n0)/nm1)) - lam_dg[fep_state-1] );
611             }
612             else
613             {
614                 dwm_array[nval]  = fabs( cnval - lam_dg[fep_state-1] );
615             }
616
617             if (n0 > 0)
618             {
619                 omegap_array[nval]             = omega_p1_0;
620             }
621             else
622             {
623                 omegap_array[nval]             = 0;
624             }
625             weightsp_array[nval]           = clam_weightsp;
626             varp_array[nval]               = clam_varp;
627             if ((np1 > 0) && (n0 > 0))
628             {
629                 dwp_array[nval]  = fabs( (cnval + log((1.0*np1)/n0)) - lam_dg[fep_state] );
630             }
631             else
632             {
633                 dwp_array[nval]  = fabs( cnval - lam_dg[fep_state] );
634             }
635
636         }
637
638         /* find the C's closest to the old weights value */
639
640         min_nvalm     = FindMinimum(dwm_array, maxc);
641         omega_m1_0    = omegam_array[min_nvalm];
642         clam_weightsm = weightsm_array[min_nvalm];
643         clam_varm     = varm_array[min_nvalm];
644
645         min_nvalp     = FindMinimum(dwp_array, maxc);
646         omega_p1_0    = omegap_array[min_nvalp];
647         clam_weightsp = weightsp_array[min_nvalp];
648         clam_varp     = varp_array[min_nvalp];
649
650         clam_osum   = omega_m1_0 + omega_p1_0;
651         clam_minvar = 0;
652         if (clam_osum > 0)
653         {
654             clam_minvar = 0.5*log(clam_osum);
655         }
656
657         if (fep_state > 0)
658         {
659             lam_dg[fep_state-1]       = clam_weightsm;
660             lam_variance[fep_state-1] = clam_varm;
661         }
662
663         if (fep_state < nlim-1)
664         {
665             lam_dg[fep_state]       = clam_weightsp;
666             lam_variance[fep_state] = clam_varp;
667         }
668
669         if (expand->elamstats == elamstatsMINVAR)
670         {
671             bSufficientSamples = TRUE;
672             /* make sure they are all past a threshold */
673             for (i = 0; i < nlim; i++)
674             {
675                 if (dfhist->n_at_lam[i] < expand->minvarmin)
676                 {
677                     bSufficientSamples = FALSE;
678                 }
679             }
680             if (bSufficientSamples)
681             {
682                 dfhist->sum_minvar[fep_state] = clam_minvar;
683                 if (fep_state == 0)
684                 {
685                     for (i = 0; i < nlim; i++)
686                     {
687                         dfhist->sum_minvar[i] += (expand->minvar_const-clam_minvar);
688                     }
689                     expand->minvar_const          = clam_minvar;
690                     dfhist->sum_minvar[fep_state] = 0.0;
691                 }
692                 else
693                 {
694                     dfhist->sum_minvar[fep_state] -= expand->minvar_const;
695                 }
696             }
697         }
698
699         /* we need to rezero minvar now, since it could change at fep_state = 0 */
700         dfhist->sum_dg[0]       = 0.0;
701         dfhist->sum_variance[0] = 0.0;
702         dfhist->sum_weights[0]  = dfhist->sum_dg[0] + dfhist->sum_minvar[0]; /* should be zero */
703
704         for (i = 1; i < nlim; i++)
705         {
706             dfhist->sum_dg[i]       = lam_dg[i-1] + dfhist->sum_dg[i-1];
707             dfhist->sum_variance[i] = sqrt(lam_variance[i-1] + pow(dfhist->sum_variance[i-1], 2));
708             dfhist->sum_weights[i]  = dfhist->sum_dg[i] + dfhist->sum_minvar[i];
709         }
710
711         sfree(lam_dg);
712         sfree(lam_variance);
713
714         sfree(omegam_array);
715         sfree(weightsm_array);
716         sfree(varm_array);
717         sfree(dwm_array);
718
719         sfree(omegap_array);
720         sfree(weightsp_array);
721         sfree(varp_array);
722         sfree(dwp_array);
723     }
724     return FALSE;
725 }
726
727 static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k, gmx_rng_t rng)
728 {
729     /* Choose new lambda value, and update transition matrix */
730
731     int      i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
732     real     r1, r2, de_old, de_new, de, trialprob, tprob = 0;
733     real   **Tij;
734     double  *propose, *accept, *remainder;
735     double   pks;
736     real     sum, pnorm;
737     gmx_bool bRestricted;
738
739     starting_fep_state = fep_state;
740     lamnew             = fep_state; /* so that there is a default setting -- stays the same */
741
742     if (!EWL(expand->elamstats))    /* ignore equilibrating the weights if using WL */
743     {
744         if ((expand->lmc_forced_nstart > 0) && (dfhist->n_at_lam[nlim-1] <= expand->lmc_forced_nstart))
745         {
746             /* Use a marching method to run through the lambdas and get preliminary free energy data,
747                before starting 'free' sampling.  We start free sampling when we have enough at each lambda */
748
749             /* if we have enough at this lambda, move on to the next one */
750
751             if (dfhist->n_at_lam[fep_state] == expand->lmc_forced_nstart)
752             {
753                 lamnew = fep_state+1;
754                 if (lamnew == nlim)  /* whoops, stepped too far! */
755                 {
756                     lamnew -= 1;
757                 }
758             }
759             else
760             {
761                 lamnew = fep_state;
762             }
763             return lamnew;
764         }
765     }
766
767     snew(propose, nlim);
768     snew(accept, nlim);
769     snew(remainder, nlim);
770
771     for (i = 0; i < expand->lmc_repeats; i++)
772     {
773
774         for (ifep = 0; ifep < nlim; ifep++)
775         {
776             propose[ifep] = 0;
777             accept[ifep]  = 0;
778         }
779
780         if ((expand->elmcmove == elmcmoveGIBBS) || (expand->elmcmove == elmcmoveMETGIBBS))
781         {
782             bRestricted = TRUE;
783             /* use the Gibbs sampler, with restricted range */
784             if (expand->gibbsdeltalam < 0)
785             {
786                 minfep      = 0;
787                 maxfep      = nlim-1;
788                 bRestricted = FALSE;
789             }
790             else
791             {
792                 minfep = fep_state - expand->gibbsdeltalam;
793                 maxfep = fep_state + expand->gibbsdeltalam;
794                 if (minfep < 0)
795                 {
796                     minfep = 0;
797                 }
798                 if (maxfep > nlim-1)
799                 {
800                     maxfep = nlim-1;
801                 }
802             }
803
804             GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, minfep, maxfep);
805
806             if (expand->elmcmove == elmcmoveGIBBS)
807             {
808                 for (ifep = minfep; ifep <= maxfep; ifep++)
809                 {
810                     propose[ifep] = p_k[ifep];
811                     accept[ifep]  = 1.0;
812                 }
813                 /* Gibbs sampling */
814                 r1 = gmx_rng_uniform_real(rng);
815                 for (lamnew = minfep; lamnew <= maxfep; lamnew++)
816                 {
817                     if (r1 <= p_k[lamnew])
818                     {
819                         break;
820                     }
821                     r1 -= p_k[lamnew];
822                 }
823             }
824             else if (expand->elmcmove == elmcmoveMETGIBBS)
825             {
826
827                 /* Metropolized Gibbs sampling */
828                 for (ifep = minfep; ifep <= maxfep; ifep++)
829                 {
830                     remainder[ifep] = 1 - p_k[ifep];
831                 }
832
833                 /* find the proposal probabilities */
834
835                 if (remainder[fep_state] == 0)
836                 {
837                     /* only the current state has any probability */
838                     /* we have to stay at the current state */
839                     lamnew = fep_state;
840                 }
841                 else
842                 {
843                     for (ifep = minfep; ifep <= maxfep; ifep++)
844                     {
845                         if (ifep != fep_state)
846                         {
847                             propose[ifep] = p_k[ifep]/remainder[fep_state];
848                         }
849                         else
850                         {
851                             propose[ifep] = 0;
852                         }
853                     }
854
855                     r1 = gmx_rng_uniform_real(rng);
856                     for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
857                     {
858                         pnorm = p_k[lamtrial]/remainder[fep_state];
859                         if (lamtrial != fep_state)
860                         {
861                             if (r1 <= pnorm)
862                             {
863                                 break;
864                             }
865                             r1 -= pnorm;
866                         }
867                     }
868
869                     /* we have now selected lamtrial according to p(lamtrial)/1-p(fep_state) */
870                     tprob = 1.0;
871                     /* trial probability is min{1,\frac{1 - p(old)}{1-p(new)} MRS 1/8/2008 */
872                     trialprob = (remainder[fep_state])/(remainder[lamtrial]);
873                     if (trialprob < tprob)
874                     {
875                         tprob = trialprob;
876                     }
877                     r2 = gmx_rng_uniform_real(rng);
878                     if (r2 < tprob)
879                     {
880                         lamnew = lamtrial;
881                     }
882                     else
883                     {
884                         lamnew = fep_state;
885                     }
886                 }
887
888                 /* now figure out the acceptance probability for each */
889                 for (ifep = minfep; ifep <= maxfep; ifep++)
890                 {
891                     tprob = 1.0;
892                     if (remainder[ifep] != 0)
893                     {
894                         trialprob = (remainder[fep_state])/(remainder[ifep]);
895                     }
896                     else
897                     {
898                         trialprob = 1.0; /* this state is the only choice! */
899                     }
900                     if (trialprob < tprob)
901                     {
902                         tprob = trialprob;
903                     }
904                     /* probability for fep_state=0, but that's fine, it's never proposed! */
905                     accept[ifep] = tprob;
906                 }
907             }
908
909             if (lamnew > maxfep)
910             {
911                 /* it's possible some rounding is failing */
912                 if (gmx_within_tol(remainder[fep_state],0,50*GMX_DOUBLE_EPS))
913                 {
914                     /* numerical rounding error -- no state other than the original has weight */
915                     lamnew = fep_state;
916                 }
917                 else
918                 {
919                     /* probably not a numerical issue */
920                     int   loc    = 0;
921                     int   nerror = 200+(maxfep-minfep+1)*60;
922                     char *errorstr;
923                     snew(errorstr, nerror);
924                     /* if its greater than maxfep, then something went wrong -- probably underflow in the calculation
925                        of sum weights. Generated detailed info for failure */
926                     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);
927                     for (ifep = minfep; ifep <= maxfep; ifep++)
928                     {
929                         loc += sprintf(&errorstr[loc], "%3d %17.10e%17.10e%17.10e\n", ifep, weighted_lamee[ifep], p_k[ifep], dfhist->sum_weights[ifep]);
930                     }
931                     gmx_fatal(FARGS, errorstr);
932                 }
933             }
934         }
935         else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
936         {
937             /* use the metropolis sampler with trial +/- 1 */
938             r1 = gmx_rng_uniform_real(rng);
939             if (r1 < 0.5)
940             {
941                 if (fep_state == 0)
942                 {
943                     lamtrial = fep_state;
944                 }
945                 else
946                 {
947                     lamtrial = fep_state-1;
948                 }
949             }
950             else
951             {
952                 if (fep_state == nlim-1)
953                 {
954                     lamtrial = fep_state;
955                 }
956                 else
957                 {
958                     lamtrial = fep_state+1;
959                 }
960             }
961
962             de = weighted_lamee[lamtrial] - weighted_lamee[fep_state];
963             if (expand->elmcmove == elmcmoveMETROPOLIS)
964             {
965                 tprob     = 1.0;
966                 trialprob = exp(de);
967                 if (trialprob < tprob)
968                 {
969                     tprob = trialprob;
970                 }
971                 propose[fep_state] = 0;
972                 propose[lamtrial]  = 1.0; /* note that this overwrites the above line if fep_state = ntrial, which only occurs at the ends */
973                 accept[fep_state]  = 1.0; /* doesn't actually matter, never proposed unless fep_state = ntrial, in which case it's 1.0 anyway */
974                 accept[lamtrial]   = tprob;
975
976             }
977             else if (expand->elmcmove == elmcmoveBARKER)
978             {
979                 tprob = 1.0/(1.0+exp(-de));
980
981                 propose[fep_state] = (1-tprob);
982                 propose[lamtrial] += tprob; /* we add, to account for the fact that at the end, they might be the same point */
983                 accept[fep_state]  = 1.0;
984                 accept[lamtrial]   = 1.0;
985             }
986
987             r2 = gmx_rng_uniform_real(rng);
988             if (r2 < tprob)
989             {
990                 lamnew = lamtrial;
991             }
992             else
993             {
994                 lamnew = fep_state;
995             }
996         }
997
998         for (ifep = 0; ifep < nlim; ifep++)
999         {
1000             dfhist->Tij[fep_state][ifep]      += propose[ifep]*accept[ifep];
1001             dfhist->Tij[fep_state][fep_state] += propose[ifep]*(1.0-accept[ifep]);
1002         }
1003         fep_state = lamnew;
1004     }
1005
1006     dfhist->Tij_empirical[starting_fep_state][lamnew] += 1.0;
1007
1008     sfree(propose);
1009     sfree(accept);
1010     sfree(remainder);
1011
1012     return lamnew;
1013 }
1014
1015 /* print out the weights to the log, along with current state */
1016 extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
1017                                       int nlam, int frequency, gmx_large_int_t step)
1018 {
1019     int         nlim, i, ifep, jfep;
1020     real        dw, dg, dv, dm, Tprint;
1021     real       *temps;
1022     const char *print_names[efptNR] = {" FEPL", "MassL", "CoulL", " VdwL", "BondL", "RestT", "Temp.(K)"};
1023     gmx_bool    bSimTemp            = FALSE;
1024
1025     nlim = fep->n_lambda;
1026     if (simtemp != NULL)
1027     {
1028         bSimTemp = TRUE;
1029     }
1030
1031     if (mod(step, frequency) == 0)
1032     {
1033         fprintf(outfile, "             MC-lambda information\n");
1034         if (EWL(expand->elamstats) && (!(dfhist->bEquil)))
1035         {
1036             fprintf(outfile, "  Wang-Landau incrementor is: %11.5g\n", dfhist->wl_delta);
1037         }
1038         fprintf(outfile, "  N");
1039         for (i = 0; i < efptNR; i++)
1040         {
1041             if (fep->separate_dvdl[i])
1042             {
1043                 fprintf(outfile, "%7s", print_names[i]);
1044             }
1045             else if ((i == efptTEMPERATURE) && bSimTemp)
1046             {
1047                 fprintf(outfile, "%10s", print_names[i]); /* more space for temperature formats */
1048             }
1049         }
1050         fprintf(outfile, "    Count   ");
1051         if (expand->elamstats == elamstatsMINVAR)
1052         {
1053             fprintf(outfile, "W(in kT)   G(in kT)  dG(in kT)  dV(in kT)\n");
1054         }
1055         else
1056         {
1057             fprintf(outfile, "G(in kT)  dG(in kT)\n");
1058         }
1059         for (ifep = 0; ifep < nlim; ifep++)
1060         {
1061             if (ifep == nlim-1)
1062             {
1063                 dw = 0.0;
1064                 dg = 0.0;
1065                 dv = 0.0;
1066                 dm = 0.0;
1067             }
1068             else
1069             {
1070                 dw = dfhist->sum_weights[ifep+1] - dfhist->sum_weights[ifep];
1071                 dg = dfhist->sum_dg[ifep+1] - dfhist->sum_dg[ifep];
1072                 dv = sqrt(pow(dfhist->sum_variance[ifep+1], 2) - pow(dfhist->sum_variance[ifep], 2));
1073                 dm = dfhist->sum_minvar[ifep+1] - dfhist->sum_minvar[ifep];
1074
1075             }
1076             fprintf(outfile, "%3d", (ifep+1));
1077             for (i = 0; i < efptNR; i++)
1078             {
1079                 if (fep->separate_dvdl[i])
1080                 {
1081                     fprintf(outfile, "%7.3f", fep->all_lambda[i][ifep]);
1082                 }
1083                 else if (i == efptTEMPERATURE && bSimTemp)
1084                 {
1085                     fprintf(outfile, "%9.3f", simtemp->temperatures[ifep]);
1086                 }
1087             }
1088             if (EWL(expand->elamstats) && (!(dfhist->bEquil)))  /* if performing WL and still haven't equilibrated */
1089             {
1090                 if (expand->elamstats == elamstatsWL)
1091                 {
1092                     fprintf(outfile, " %8d", (int)dfhist->wl_histo[ifep]);
1093                 }
1094                 else
1095                 {
1096                     fprintf(outfile, " %8.3f", dfhist->wl_histo[ifep]);
1097                 }
1098             }
1099             else   /* we have equilibrated weights */
1100             {
1101                 fprintf(outfile, " %8d", dfhist->n_at_lam[ifep]);
1102             }
1103             if (expand->elamstats == elamstatsMINVAR)
1104             {
1105                 fprintf(outfile, " %10.5f %10.5f %10.5f %10.5f", dfhist->sum_weights[ifep], dfhist->sum_dg[ifep], dg, dv);
1106             }
1107             else
1108             {
1109                 fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
1110             }
1111             if (ifep == nlam)
1112             {
1113                 fprintf(outfile, " <<\n");
1114             }
1115             else
1116             {
1117                 fprintf(outfile, "   \n");
1118             }
1119         }
1120         fprintf(outfile, "\n");
1121
1122         if ((mod(step, expand->nstTij) == 0) && (expand->nstTij > 0) && (step > 0))
1123         {
1124             fprintf(outfile, "                     Transition Matrix\n");
1125             for (ifep = 0; ifep < nlim; ifep++)
1126             {
1127                 fprintf(outfile, "%12d", (ifep+1));
1128             }
1129             fprintf(outfile, "\n");
1130             for (ifep = 0; ifep < nlim; ifep++)
1131             {
1132                 for (jfep = 0; jfep < nlim; jfep++)
1133                 {
1134                     if (dfhist->n_at_lam[ifep] > 0)
1135                     {
1136                         if (expand->bSymmetrizedTMatrix)
1137                         {
1138                             Tprint = (dfhist->Tij[ifep][jfep]+dfhist->Tij[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1139                         }
1140                         else
1141                         {
1142                             Tprint = (dfhist->Tij[ifep][jfep])/(dfhist->n_at_lam[ifep]);
1143                         }
1144                     }
1145                     else
1146                     {
1147                         Tprint = 0.0;
1148                     }
1149                     fprintf(outfile, "%12.8f", Tprint);
1150                 }
1151                 fprintf(outfile, "%3d\n", (ifep+1));
1152             }
1153
1154             fprintf(outfile, "                  Empirical Transition Matrix\n");
1155             for (ifep = 0; ifep < nlim; ifep++)
1156             {
1157                 fprintf(outfile, "%12d", (ifep+1));
1158             }
1159             fprintf(outfile, "\n");
1160             for (ifep = 0; ifep < nlim; ifep++)
1161             {
1162                 for (jfep = 0; jfep < nlim; jfep++)
1163                 {
1164                     if (dfhist->n_at_lam[ifep] > 0)
1165                     {
1166                         if (expand->bSymmetrizedTMatrix)
1167                         {
1168                             Tprint = (dfhist->Tij_empirical[ifep][jfep]+dfhist->Tij_empirical[jfep][ifep])/(dfhist->n_at_lam[ifep]+dfhist->n_at_lam[jfep]);
1169                         }
1170                         else
1171                         {
1172                             Tprint = dfhist->Tij_empirical[ifep][jfep]/(dfhist->n_at_lam[ifep]);
1173                         }
1174                     }
1175                     else
1176                     {
1177                         Tprint = 0.0;
1178                     }
1179                     fprintf(outfile, "%12.8f", Tprint);
1180                 }
1181                 fprintf(outfile, "%3d\n", (ifep+1));
1182             }
1183         }
1184     }
1185 }
1186
1187 extern void get_mc_state(gmx_rng_t rng, t_state *state)
1188 {
1189     gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
1190 }
1191
1192 extern void set_mc_state(gmx_rng_t rng, t_state *state)
1193 {
1194     gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
1195 }
1196
1197 extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
1198                                     t_state *state, t_extmass *MassQ, df_history_t *dfhist,
1199                                     gmx_large_int_t step, gmx_rng_t mcrng,
1200                                     rvec *v, t_mdatoms *mdatoms)
1201 {
1202     real       *pfep_lamee, *scaled_lamee, *weighted_lamee;
1203     double     *p_k;
1204     int         i, nlam, nlim, lamnew, totalsamples;
1205     real        oneovert, maxscaled = 0, maxweighted = 0;
1206     t_expanded *expand;
1207     t_simtemp  *simtemp;
1208     double     *temperature_lambdas;
1209     gmx_bool    bIfReset, bSwitchtoOneOverT, bDoneEquilibrating = FALSE;
1210
1211     expand  = ir->expandedvals;
1212     simtemp = ir->simtempvals;
1213     nlim    = ir->fepvals->n_lambda;
1214     nlam    = state->fep_state;
1215
1216     snew(scaled_lamee, nlim);
1217     snew(weighted_lamee, nlim);
1218     snew(pfep_lamee, nlim);
1219     snew(p_k, nlim);
1220
1221     if (expand->bInit_weights)                    /* if initialized weights, we need to fill them in */
1222     {
1223         dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
1224         for (i = 0; i < nlim; i++)
1225         {
1226             dfhist->sum_weights[i] = expand->init_lambda_weights[i];
1227             dfhist->sum_dg[i]      = expand->init_lambda_weights[i];
1228         }
1229         expand->bInit_weights = FALSE;
1230     }
1231
1232     /* update the count at the current lambda*/
1233     dfhist->n_at_lam[nlam]++;
1234
1235     /* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
1236        pressure controlled.*/
1237     /*
1238        pVTerm = 0;
1239        where does this PV term go?
1240        for (i=0;i<nlim;i++)
1241        {
1242        fep_lamee[i] += pVTerm;
1243        }
1244      */
1245
1246     /* determine the minimum value to avoid overflow.  Probably a better way to do this */
1247     /* we don't need to include the pressure term, since the volume is the same between the two.
1248        is there some term we are neglecting, however? */
1249
1250     if (ir->efep != efepNO)
1251     {
1252         for (i = 0; i < nlim; i++)
1253         {
1254             if (ir->bSimTemp)
1255             {
1256                 /* Note -- this assumes no mass changes, since kinetic energy is not added  . . . */
1257                 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
1258                     + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
1259             }
1260             else
1261             {
1262                 scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(expand->mc_temp*BOLTZ);
1263                 /* mc_temp is currently set to the system reft unless otherwise defined */
1264             }
1265
1266             /* save these energies for printing, so they don't get overwritten by the next step */
1267             /* they aren't overwritten in the non-free energy case, but we always print with these
1268                for simplicity */
1269         }
1270     }
1271     else
1272     {
1273         if (ir->bSimTemp)
1274         {
1275             for (i = 0; i < nlim; i++)
1276             {
1277                 scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
1278             }
1279         }
1280     }
1281
1282     for (i = 0; i < nlim; i++)
1283     {
1284         pfep_lamee[i] = scaled_lamee[i];
1285
1286         weighted_lamee[i] = dfhist->sum_weights[i] - scaled_lamee[i];
1287         if (i == 0)
1288         {
1289             maxscaled   = scaled_lamee[i];
1290             maxweighted = weighted_lamee[i];
1291         }
1292         else
1293         {
1294             if (scaled_lamee[i] > maxscaled)
1295             {
1296                 maxscaled = scaled_lamee[i];
1297             }
1298             if (weighted_lamee[i] > maxweighted)
1299             {
1300                 maxweighted = weighted_lamee[i];
1301             }
1302         }
1303     }
1304
1305     for (i = 0; i < nlim; i++)
1306     {
1307         scaled_lamee[i]   -= maxscaled;
1308         weighted_lamee[i] -= maxweighted;
1309     }
1310
1311     /* update weights - we decide whether or not to actually do this inside */
1312
1313     bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, nlam, scaled_lamee, weighted_lamee, step);
1314     if (bDoneEquilibrating)
1315     {
1316         if (log)
1317         {
1318             fprintf(log, "\nStep %d: Weights have equilibrated, using criteria: %s\n", (int)step, elmceq_names[expand->elmceq]);
1319         }
1320     }
1321
1322     lamnew = ChooseNewLambda(nlim, expand, dfhist, nlam, weighted_lamee, p_k, mcrng);
1323     /* if using simulated tempering, we need to adjust the temperatures */
1324     if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
1325     {
1326         int   i, j, n, d;
1327         real *buf_ngtc;
1328         real  told;
1329         int   nstart, nend, gt;
1330
1331         snew(buf_ngtc, ir->opts.ngtc);
1332
1333         for (i = 0; i < ir->opts.ngtc; i++)
1334         {
1335             if (ir->opts.ref_t[i] > 0)
1336             {
1337                 told              = ir->opts.ref_t[i];
1338                 ir->opts.ref_t[i] =  simtemp->temperatures[lamnew];
1339                 buf_ngtc[i]       = sqrt(ir->opts.ref_t[i]/told); /* using the buffer as temperature scaling */
1340             }
1341         }
1342
1343         /* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
1344
1345         nstart = mdatoms->start;
1346         nend   = nstart + mdatoms->homenr;
1347         for (n = nstart; n < nend; n++)
1348         {
1349             gt = 0;
1350             if (mdatoms->cTC)
1351             {
1352                 gt = mdatoms->cTC[n];
1353             }
1354             for (d = 0; d < DIM; d++)
1355             {
1356                 v[n][d] *= buf_ngtc[gt];
1357             }
1358         }
1359
1360         if (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir) || IR_NVT_TROTTER(ir))
1361         {
1362             /* we need to recalculate the masses if the temperature has changed */
1363             init_npt_masses(ir, state, MassQ, FALSE);
1364             for (i = 0; i < state->nnhpres; i++)
1365             {
1366                 for (j = 0; j < ir->opts.nhchainlength; j++)
1367                 {
1368                     state->nhpres_vxi[i+j] *= buf_ngtc[i];
1369                 }
1370             }
1371             for (i = 0; i < ir->opts.ngtc; i++)
1372             {
1373                 for (j = 0; j < ir->opts.nhchainlength; j++)
1374                 {
1375                     state->nosehoover_vxi[i+j] *= buf_ngtc[i];
1376                 }
1377             }
1378         }
1379         sfree(buf_ngtc);
1380     }
1381
1382     /* now check on the Wang-Landau updating critera */
1383
1384     if (EWL(expand->elamstats))
1385     {
1386         bSwitchtoOneOverT = FALSE;
1387         if (expand->bWLoneovert)
1388         {
1389             totalsamples = 0;
1390             for (i = 0; i < nlim; i++)
1391             {
1392                 totalsamples += dfhist->n_at_lam[i];
1393             }
1394             oneovert = (1.0*nlim)/totalsamples;
1395             /* oneovert has decreasd by a bit since last time, so we actually make sure its within one of this number */
1396             /* switch to 1/t incrementing when wl_delta has decreased at least once, and wl_delta is now less than 1/t */
1397             if ((dfhist->wl_delta <= ((totalsamples)/(totalsamples-1.00001))*oneovert) &&
1398                 (dfhist->wl_delta < expand->init_wl_delta))
1399             {
1400                 bSwitchtoOneOverT = TRUE;
1401             }
1402         }
1403         if (bSwitchtoOneOverT)
1404         {
1405             dfhist->wl_delta = oneovert; /* now we reduce by this each time, instead of only at flatness */
1406         }
1407         else
1408         {
1409             bIfReset = CheckHistogramRatios(nlim, dfhist->wl_histo, expand->wl_ratio);
1410             if (bIfReset)
1411             {
1412                 for (i = 0; i < nlim; i++)
1413                 {
1414                     dfhist->wl_histo[i] = 0;
1415                 }
1416                 dfhist->wl_delta *= expand->wl_scale;
1417                 if (log)
1418                 {
1419                     fprintf(log, "\nStep %d: weights are now:", (int)step);
1420                     for (i = 0; i < nlim; i++)
1421                     {
1422                         fprintf(log, " %.5f", dfhist->sum_weights[i]);
1423                     }
1424                     fprintf(log, "\n");
1425                 }
1426             }
1427         }
1428     }
1429     sfree(pfep_lamee);
1430     sfree(scaled_lamee);
1431     sfree(weighted_lamee);
1432     sfree(p_k);
1433
1434     return lamnew;
1435 }