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