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