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