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