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