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