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