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