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