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