2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,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.
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.
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.
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.
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.
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.
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/fileio/xvgr.h"
52 #include "gromacs/legacyheaders/viewit.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/mdebin.h"
60 /* Structure for the names of lambda vector components */
61 typedef struct lambda_components_t
63 char **names; /* Array of strings with names for the lambda vector
65 int N; /* The number of components */
66 int Nalloc; /* The number of allocated components */
67 } lambda_components_t;
69 /* Structure for a lambda vector or a dhdl derivative direction */
70 typedef struct lambda_vec_t
72 double *val; /* The lambda vector component values. Only valid if
74 int dhdl; /* The coordinate index for the derivative described by this
76 const lambda_components_t *lc; /* the associated lambda_components
78 int index; /* The state number (init-lambda-state) of this lambda
79 vector, if known. If not, it is set to -1 */
82 /* the dhdl.xvg data from a simulation */
86 int ftp; /* file type */
87 int nset; /* number of lambdas, including dhdl */
88 int *np; /* number of data points (du or hists) per lambda */
89 int np_alloc; /* number of points (du or hists) allocated */
90 double temp; /* temperature */
91 lambda_vec_t *lambda; /* the lambdas (of first index for y). */
92 double *t; /* the times (of second index for y) */
93 double **y; /* the dU values. y[0] holds the derivative, while
94 further ones contain the energy differences between
95 the native lambda and the 'foreign' lambdas. */
96 lambda_vec_t native_lambda; /* the native lambda */
98 struct xvg_t *next, *prev; /*location in the global linked list of xvg_ts*/
102 typedef struct hist_t
104 unsigned int *bin[2]; /* the (forward + reverse) histogram values */
105 double dx[2]; /* the histogram spacing. The reverse
106 dx is the negative of the forward dx.*/
107 gmx_int64_t x0[2]; /* the (forward + reverse) histogram start
110 int nbin[2]; /* the (forward+reverse) number of bins */
111 gmx_int64_t sum; /* the total number of counts. Must be
112 the same for forward + reverse. */
113 int nhist; /* number of hist datas (forward or reverse) */
115 double start_time, delta_time; /* start time, end time of histogram */
119 /* an aggregate of samples for partial free energy calculation */
120 typedef struct samples_t
122 lambda_vec_t *native_lambda; /* pointer to native lambda vector */
123 lambda_vec_t *foreign_lambda; /* pointer to foreign lambda vector */
124 double temp; /* the temperature */
125 gmx_bool derivative; /* whether this sample is a derivative */
127 /* The samples come either as either delta U lists: */
128 int ndu; /* the number of delta U samples */
129 double *du; /* the delta u's */
130 double *t; /* the times associated with those samples, or: */
131 double start_time, delta_time; /*start time and delta time for linear time*/
133 /* or as histograms: */
134 hist_t *hist; /* a histogram */
136 /* allocation data: (not NULL for data 'owned' by this struct) */
137 double *du_alloc, *t_alloc; /* allocated delta u arrays */
138 size_t ndu_alloc, nt_alloc; /* pre-allocated sizes */
139 hist_t *hist_alloc; /* allocated hist */
141 gmx_int64_t ntot; /* total number of samples */
142 const char *filename; /* the file name this sample comes from */
145 /* a sample range (start to end for du-style data, or boolean
146 for both du-style data and histograms */
147 typedef struct sample_range_t
149 int start, end; /* start and end index for du style data */
150 gmx_bool use; /* whether to use this sample */
152 samples_t *s; /* the samples this range belongs to */
156 /* a collection of samples for a partial free energy calculation
157 (i.e. the collection of samples from one native lambda to one
159 typedef struct sample_coll_t
161 lambda_vec_t *native_lambda; /* these should be the same for all samples
163 lambda_vec_t *foreign_lambda; /* collection */
164 double temp; /* the temperature */
166 int nsamples; /* the number of samples */
167 samples_t **s; /* the samples themselves */
168 sample_range_t *r; /* the sample ranges */
169 int nsamples_alloc; /* number of allocated samples */
171 gmx_int64_t ntot; /* total number of samples in the ranges of
174 struct sample_coll_t *next, *prev; /* next and previous in the list */
177 /* all the samples associated with a lambda point */
178 typedef struct lambda_data_t
180 lambda_vec_t *lambda; /* the native lambda (at start time if dynamic) */
181 double temp; /* temperature */
183 sample_coll_t *sc; /* the samples */
185 sample_coll_t sc_head; /*the pre-allocated list head for the linked list.*/
187 struct lambda_data_t *next, *prev; /* the next and prev in the list */
190 /* Top-level data structure of simulation data */
191 typedef struct sim_data_t
193 lambda_data_t *lb; /* a lambda data linked list */
194 lambda_data_t lb_head; /* The head element of the linked list */
196 lambda_components_t lc; /* the allowed components of the lambda
200 /* Top-level data structure with calculated values. */
202 sample_coll_t *a, *b; /* the simulation data */
204 double dg; /* the free energy difference */
205 double dg_err; /* the free energy difference */
207 double dg_disc_err; /* discretization error */
208 double dg_histrange_err; /* histogram range error */
210 double sa; /* relative entropy of b in state a */
211 double sa_err; /* error in sa */
212 double sb; /* relative entropy of a in state b */
213 double sb_err; /* error in sb */
215 double dg_stddev; /* expected dg stddev per sample */
216 double dg_stddev_err; /* error in dg_stddev */
220 /* Initialize a lambda_components structure */
221 static void lambda_components_init(lambda_components_t *lc)
225 snew(lc->names, lc->Nalloc);
228 /* Add a component to a lambda_components structure */
229 static void lambda_components_add(lambda_components_t *lc,
230 const char *name, size_t name_length)
232 while (lc->N + 1 > lc->Nalloc)
234 lc->Nalloc = (lc->Nalloc == 0) ? 2 : 2*lc->Nalloc;
235 srenew(lc->names, lc->Nalloc);
237 snew(lc->names[lc->N], name_length+1);
238 strncpy(lc->names[lc->N], name, name_length);
242 /* check whether a component with index 'index' matches the given name, or
243 is also NULL. Returns TRUE if this is the case.
244 the string name does not need to end */
245 static gmx_bool lambda_components_check(const lambda_components_t *lc,
255 if (name == NULL && lc->names[index] == NULL)
259 if ((name == NULL) != (lc->names[index] == NULL))
263 len = strlen(lc->names[index]);
264 if (len != name_length)
268 if (strncmp(lc->names[index], name, name_length) == 0)
275 /* Find the index of a given lambda component name, or -1 if not found */
276 static int lambda_components_find(const lambda_components_t *lc,
282 for (i = 0; i < lc->N; i++)
284 if (strncmp(lc->names[i], name, name_length) == 0)
294 /* initialize a lambda vector */
295 static void lambda_vec_init(lambda_vec_t *lv, const lambda_components_t *lc)
297 snew(lv->val, lc->N);
303 static void lambda_vec_destroy(lambda_vec_t *lv)
308 static void lambda_vec_copy(lambda_vec_t *lv, const lambda_vec_t *orig)
312 lambda_vec_init(lv, orig->lc);
313 lv->dhdl = orig->dhdl;
314 lv->index = orig->index;
315 for (i = 0; i < lv->lc->N; i++)
317 lv->val[i] = orig->val[i];
321 /* write a lambda vec to a preallocated string */
322 static void lambda_vec_print(const lambda_vec_t *lv, char *str, gmx_bool named)
327 str[0] = 0; /* reset the string */
332 str += sprintf(str, "delta H to ");
336 str += sprintf(str, "(");
338 for (i = 0; i < lv->lc->N; i++)
340 str += sprintf(str, "%g", lv->val[i]);
343 str += sprintf(str, ", ");
348 str += sprintf(str, ")");
353 /* this lambda vector describes a derivative */
354 str += sprintf(str, "dH/dl");
355 if (strlen(lv->lc->names[lv->dhdl]) > 0)
357 sprintf(str, " (%s)", lv->lc->names[lv->dhdl]);
362 /* write a shortened version of the lambda vec to a preallocated string */
363 static void lambda_vec_print_short(const lambda_vec_t *lv, char *str)
370 sprintf(str, "%6d", lv->index);
376 sprintf(str, "%6.3f", lv->val[0]);
380 sprintf(str, "dH/dl[%d]", lv->dhdl);
385 /* write an intermediate version of two lambda vecs to a preallocated string */
386 static void lambda_vec_print_intermediate(const lambda_vec_t *a,
387 const lambda_vec_t *b, char *str)
393 if ( (a->index >= 0) && (b->index >= 0) )
395 sprintf(str, "%6.3f", ((double)a->index+(double)b->index)/2.);
399 if ( (a->dhdl < 0) && (b->dhdl < 0) )
401 sprintf(str, "%6.3f", (a->val[0]+b->val[0])/2.);
408 /* calculate the difference in lambda vectors: c = a-b.
409 c must be initialized already, and a and b must describe non-derivative
411 static void lambda_vec_diff(const lambda_vec_t *a, const lambda_vec_t *b,
416 if ( (a->dhdl > 0) || (b->dhdl > 0) )
419 "Trying to calculate the difference between derivatives instead of lambda points");
421 if ((a->lc != b->lc) || (a->lc != c->lc) )
424 "Trying to calculate the difference lambdas with differing basis set");
426 for (i = 0; i < a->lc->N; i++)
428 c->val[i] = a->val[i] - b->val[i];
432 /* calculate and return the absolute difference in lambda vectors: c = |a-b|.
433 a and b must describe non-derivative lambda points */
434 static double lambda_vec_abs_diff(const lambda_vec_t *a, const lambda_vec_t *b)
439 if ( (a->dhdl > 0) || (b->dhdl > 0) )
442 "Trying to calculate the difference between derivatives instead of lambda points");
447 "Trying to calculate the difference lambdas with differing basis set");
449 for (i = 0; i < a->lc->N; i++)
451 double df = a->val[i] - b->val[i];
458 /* check whether two lambda vectors are the same */
459 static gmx_bool lambda_vec_same(const lambda_vec_t *a, const lambda_vec_t *b)
469 for (i = 0; i < a->lc->N; i++)
471 if (!gmx_within_tol(a->val[i], b->val[i], 10*GMX_REAL_EPS))
480 /* they're derivatives, so we check whether the indices match */
481 return (a->dhdl == b->dhdl);
485 /* Compare the sort order of two foreign lambda vectors
487 returns 1 if a is 'bigger' than b,
488 returns 0 if they're the same,
489 returns -1 if a is 'smaller' than b.*/
490 static gmx_bool lambda_vec_cmp_foreign(const lambda_vec_t *a,
491 const lambda_vec_t *b)
494 double norm_a = 0, norm_b = 0;
495 gmx_bool different = FALSE;
499 gmx_fatal(FARGS, "Can't compare lambdas with differing basis sets");
501 /* if either one has an index we sort based on that */
502 if ((a->index >= 0) || (b->index >= 0))
504 if (a->index == b->index)
508 return (a->index > b->index) ? 1 : -1;
510 if (a->dhdl >= 0 || b->dhdl >= 0)
512 /* lambda vectors that are derivatives always sort higher than those
513 without derivatives */
514 if ((a->dhdl >= 0) != (b->dhdl >= 0) )
516 return (a->dhdl >= 0) ? 1 : -1;
518 return a->dhdl > b->dhdl;
521 /* neither has an index, so we can only sort on the lambda components,
522 which is only valid if there is one component */
523 for (i = 0; i < a->lc->N; i++)
525 if (!gmx_within_tol(a->val[i], b->val[i], 10*GMX_REAL_EPS))
529 norm_a += a->val[i]*a->val[i];
530 norm_b += b->val[i]*b->val[i];
536 return norm_a > norm_b;
539 /* Compare the sort order of two native lambda vectors
541 returns 1 if a is 'bigger' than b,
542 returns 0 if they're the same,
543 returns -1 if a is 'smaller' than b.*/
544 static gmx_bool lambda_vec_cmp_native(const lambda_vec_t *a,
545 const lambda_vec_t *b)
551 gmx_fatal(FARGS, "Can't compare lambdas with differing basis sets");
553 /* if either one has an index we sort based on that */
554 if ((a->index >= 0) || (b->index >= 0))
556 if (a->index == b->index)
560 return (a->index > b->index) ? 1 : -1;
562 /* neither has an index, so we can only sort on the lambda components,
563 which is only valid if there is one component */
567 "Can't compare lambdas with no index and > 1 component");
569 if (a->dhdl >= 0 || b->dhdl >= 0)
572 "Can't compare native lambdas that are derivatives");
574 if (gmx_within_tol(a->val[0], b->val[0], 10*GMX_REAL_EPS))
578 return a->val[0] > b->val[0] ? 1 : -1;
584 static void hist_init(hist_t *h, int nhist, int *nbin)
589 gmx_fatal(FARGS, "histogram with more than two sets of data!");
591 for (i = 0; i < nhist; i++)
593 snew(h->bin[i], nbin[i]);
595 h->nbin[i] = nbin[i];
596 h->start_time = h->delta_time = 0;
603 static void hist_destroy(hist_t *h)
609 static void xvg_init(xvg_t *ba)
618 static void samples_init(samples_t *s, lambda_vec_t *native_lambda,
619 lambda_vec_t *foreign_lambda, double temp,
620 gmx_bool derivative, const char *filename)
622 s->native_lambda = native_lambda;
623 s->foreign_lambda = foreign_lambda;
625 s->derivative = derivative;
630 s->start_time = s->delta_time = 0;
634 s->hist_alloc = NULL;
639 s->filename = filename;
642 static void sample_range_init(sample_range_t *r, samples_t *s)
650 static void sample_coll_init(sample_coll_t *sc, lambda_vec_t *native_lambda,
651 lambda_vec_t *foreign_lambda, double temp)
653 sc->native_lambda = native_lambda;
654 sc->foreign_lambda = foreign_lambda;
660 sc->nsamples_alloc = 0;
663 sc->next = sc->prev = NULL;
666 static void sample_coll_destroy(sample_coll_t *sc)
668 /* don't free the samples themselves */
674 static void lambda_data_init(lambda_data_t *l, lambda_vec_t *native_lambda,
677 l->lambda = native_lambda;
683 l->sc = &(l->sc_head);
685 sample_coll_init(l->sc, native_lambda, NULL, 0.);
690 static void barres_init(barres_t *br)
699 br->dg_stddev_err = 0;
706 /* calculate the total number of samples in a sample collection */
707 static void sample_coll_calc_ntot(sample_coll_t *sc)
712 for (i = 0; i < sc->nsamples; i++)
718 sc->ntot += sc->s[i]->ntot;
722 sc->ntot += sc->r[i].end - sc->r[i].start;
729 /* find the barsamples_t associated with a lambda that corresponds to
730 a specific foreign lambda */
731 static sample_coll_t *lambda_data_find_sample_coll(lambda_data_t *l,
732 lambda_vec_t *foreign_lambda)
734 sample_coll_t *sc = l->sc->next;
738 if (lambda_vec_same(sc->foreign_lambda, foreign_lambda))
748 /* insert li into an ordered list of lambda_colls */
749 static void lambda_data_insert_sample_coll(lambda_data_t *l, sample_coll_t *sc)
751 sample_coll_t *scn = l->sc->next;
752 while ( (scn != l->sc) )
754 if (lambda_vec_cmp_foreign(scn->foreign_lambda, sc->foreign_lambda) > 0)
760 /* now insert it before the found scn */
762 sc->prev = scn->prev;
763 scn->prev->next = sc;
767 /* insert li into an ordered list of lambdas */
768 static void lambda_data_insert_lambda(lambda_data_t *head, lambda_data_t *li)
770 lambda_data_t *lc = head->next;
773 if (lambda_vec_cmp_native(lc->lambda, li->lambda) > 0)
779 /* now insert ourselves before the found lc */
786 /* insert a sample and a sample_range into a sample_coll. The
787 samples are stored as a pointer, the range is copied. */
788 static void sample_coll_insert_sample(sample_coll_t *sc, samples_t *s,
791 /* first check if it belongs here */
792 if (sc->temp != s->temp)
794 gmx_fatal(FARGS, "Temperatures in files %s and %s are not the same!",
795 s->filename, sc->next->s[0]->filename);
797 if (!lambda_vec_same(sc->native_lambda, s->native_lambda))
799 gmx_fatal(FARGS, "Native lambda in files %s and %s are not the same (and they should be)!",
800 s->filename, sc->next->s[0]->filename);
802 if (!lambda_vec_same(sc->foreign_lambda, s->foreign_lambda))
804 gmx_fatal(FARGS, "Foreign lambda in files %s and %s are not the same (and they should be)!",
805 s->filename, sc->next->s[0]->filename);
808 /* check if there's room */
809 if ( (sc->nsamples + 1) > sc->nsamples_alloc)
811 sc->nsamples_alloc = max(2*sc->nsamples_alloc, 2);
812 srenew(sc->s, sc->nsamples_alloc);
813 srenew(sc->r, sc->nsamples_alloc);
815 sc->s[sc->nsamples] = s;
816 sc->r[sc->nsamples] = *r;
819 sample_coll_calc_ntot(sc);
822 /* insert a sample into a lambda_list, creating the right sample_coll if
824 static void lambda_data_list_insert_sample(lambda_data_t *head, samples_t *s)
826 gmx_bool found = FALSE;
830 lambda_data_t *l = head->next;
832 /* first search for the right lambda_data_t */
835 if (lambda_vec_same(l->lambda, s->native_lambda) )
845 snew(l, 1); /* allocate a new one */
846 lambda_data_init(l, s->native_lambda, s->temp); /* initialize it */
847 lambda_data_insert_lambda(head, l); /* add it to the list */
850 /* now look for a sample collection */
851 sc = lambda_data_find_sample_coll(l, s->foreign_lambda);
854 snew(sc, 1); /* allocate a new one */
855 sample_coll_init(sc, s->native_lambda, s->foreign_lambda, s->temp);
856 lambda_data_insert_sample_coll(l, sc);
859 /* now insert the samples into the sample coll */
860 sample_range_init(&r, s);
861 sample_coll_insert_sample(sc, s, &r);
865 /* make a histogram out of a sample collection */
866 static void sample_coll_make_hist(sample_coll_t *sc, int **bin,
867 int *nbin_alloc, int *nbin,
868 double *dx, double *xmin, int nbin_default)
871 gmx_bool dx_set = FALSE;
872 gmx_bool xmin_set = FALSE;
874 gmx_bool xmax_set = FALSE;
875 gmx_bool xmax_set_hard = FALSE; /* whether the xmax is bounded by the
876 limits of a histogram */
879 /* first determine dx and xmin; try the histograms */
880 for (i = 0; i < sc->nsamples; i++)
884 hist_t *hist = sc->s[i]->hist;
885 for (k = 0; k < hist->nhist; k++)
887 double hdx = hist->dx[k];
888 double xmax_now = (hist->x0[k]+hist->nbin[k])*hdx;
890 /* we use the biggest dx*/
891 if ( (!dx_set) || hist->dx[0] > *dx)
896 if ( (!xmin_set) || (hist->x0[k]*hdx) < *xmin)
899 *xmin = (hist->x0[k]*hdx);
902 if ( (!xmax_set) || (xmax_now > xmax && !xmax_set_hard) )
906 if (hist->bin[k][hist->nbin[k]-1] != 0)
908 xmax_set_hard = TRUE;
911 if (hist->bin[k][hist->nbin[k]-1] != 0 && (xmax_now < xmax) )
913 xmax_set_hard = TRUE;
919 /* and the delta us */
920 for (i = 0; i < sc->nsamples; i++)
922 if (sc->s[i]->ndu > 0)
924 /* determine min and max */
925 int starti = sc->r[i].start;
926 int endi = sc->r[i].end;
927 double du_xmin = sc->s[i]->du[starti];
928 double du_xmax = sc->s[i]->du[starti];
929 for (j = starti+1; j < endi; j++)
931 if (sc->s[i]->du[j] < du_xmin)
933 du_xmin = sc->s[i]->du[j];
935 if (sc->s[i]->du[j] > du_xmax)
937 du_xmax = sc->s[i]->du[j];
941 /* and now change the limits */
942 if ( (!xmin_set) || (du_xmin < *xmin) )
947 if ( (!xmax_set) || ((du_xmax > xmax) && !xmax_set_hard) )
955 if (!xmax_set || !xmin_set)
964 *nbin = nbin_default;
965 *dx = (xmax-(*xmin))/((*nbin)-2); /* -2 because we want the last bin to
966 be 0, and we count from 0 */
970 *nbin = (xmax-(*xmin))/(*dx);
973 if (*nbin > *nbin_alloc)
976 srenew(*bin, *nbin_alloc);
979 /* reset the histogram */
980 for (i = 0; i < (*nbin); i++)
985 /* now add the actual data */
986 for (i = 0; i < sc->nsamples; i++)
990 hist_t *hist = sc->s[i]->hist;
991 for (k = 0; k < hist->nhist; k++)
993 double hdx = hist->dx[k];
994 double xmin_hist = hist->x0[k]*hdx;
995 for (j = 0; j < hist->nbin[k]; j++)
997 /* calculate the bin corresponding to the middle of the
999 double x = hdx*(j+0.5) + xmin_hist;
1000 int binnr = (int)((x-(*xmin))/(*dx));
1002 if (binnr >= *nbin || binnr < 0)
1007 (*bin)[binnr] += hist->bin[k][j];
1013 int starti = sc->r[i].start;
1014 int endi = sc->r[i].end;
1015 for (j = starti; j < endi; j++)
1017 int binnr = (int)((sc->s[i]->du[j]-(*xmin))/(*dx));
1018 if (binnr >= *nbin || binnr < 0)
1029 /* write a collection of histograms to a file */
1030 void sim_data_histogram(sim_data_t *sd, const char *filename,
1031 int nbin_default, const output_env_t oenv)
1033 char label_x[STRLEN];
1034 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda";
1035 const char *title = "N(\\DeltaH)";
1036 const char *label_y = "Samples";
1040 char **setnames = NULL;
1041 gmx_bool first_set = FALSE;
1042 /* histogram data: */
1049 lambda_data_t *bl_head = sd->lb;
1051 printf("\nWriting histogram to %s\n", filename);
1052 sprintf(label_x, "\\DeltaH (%s)", unit_energy);
1054 fp = xvgropen_type(filename, title, label_x, label_y, exvggtXNY, oenv);
1056 /* first get all the set names */
1058 /* iterate over all lambdas */
1059 while (bl != bl_head)
1061 sample_coll_t *sc = bl->sc->next;
1063 /* iterate over all samples */
1064 while (sc != bl->sc)
1066 char buf[STRLEN], buf2[STRLEN];
1069 srenew(setnames, nsets);
1070 snew(setnames[nsets-1], STRLEN);
1071 if (sc->foreign_lambda->dhdl < 0)
1073 lambda_vec_print(sc->native_lambda, buf, FALSE);
1074 lambda_vec_print(sc->foreign_lambda, buf2, FALSE);
1075 sprintf(setnames[nsets-1], "N(%s(%s=%s) | %s=%s)",
1076 deltag, lambda, buf2, lambda, buf);
1080 lambda_vec_print(sc->native_lambda, buf, FALSE);
1081 sprintf(setnames[nsets-1], "N(%s | %s=%s)",
1089 xvgr_legend(fp, nsets, (const char**)setnames, oenv);
1092 /* now make the histograms */
1094 /* iterate over all lambdas */
1095 while (bl != bl_head)
1097 sample_coll_t *sc = bl->sc->next;
1099 /* iterate over all samples */
1100 while (sc != bl->sc)
1104 xvgr_new_dataset(fp, 0, 0, NULL, oenv);
1107 sample_coll_make_hist(sc, &hist, &nbin_alloc, &nbin, &dx, &min,
1110 for (i = 0; i < nbin; i++)
1112 double xmin = i*dx + min;
1113 double xmax = (i+1)*dx + min;
1115 fprintf(fp, "%g %d\n%g %d\n", xmin, hist[i], xmax, hist[i]);
1133 /* create a collection (array) of barres_t object given a ordered linked list
1134 of barlamda_t sample collections */
1135 static barres_t *barres_list_create(sim_data_t *sd, int *nres,
1142 gmx_bool dhdl = FALSE;
1143 gmx_bool first = TRUE;
1144 lambda_data_t *bl_head = sd->lb;
1146 /* first count the lambdas */
1148 while (bl != bl_head)
1153 snew(res, nlambda-1);
1155 /* next put the right samples in the res */
1157 bl = bl_head->next->next; /* we start with the second one. */
1158 while (bl != bl_head)
1160 sample_coll_t *sc, *scprev;
1161 barres_t *br = &(res[*nres]);
1162 /* there is always a previous one. we search for that as a foreign
1164 scprev = lambda_data_find_sample_coll(bl->prev, bl->lambda);
1165 sc = lambda_data_find_sample_coll(bl, bl->prev->lambda);
1173 scprev = lambda_data_find_sample_coll(bl->prev, bl->prev->lambda);
1174 sc = lambda_data_find_sample_coll(bl, bl->lambda);
1178 printf("\nWARNING: Using the derivative data (dH/dlambda) to extrapolate delta H values.\nThis will only work if the Hamiltonian is linear in lambda.\n");
1183 gmx_fatal(FARGS, "Some dhdl files contain only one value (dH/dl), while others \ncontain multiple values (dH/dl and/or Delta H), will not proceed \nbecause of possible inconsistencies.\n");
1186 else if (!scprev && !sc)
1188 gmx_fatal(FARGS, "There is no path from lambda=%f -> %f that is covered by foreign lambdas:\ncannot proceed with BAR.\nUse thermodynamic integration of dH/dl by calculating the averages of dH/dl\nwith g_analyze and integrating them.\nAlternatively, use the -extp option if (and only if) the Hamiltonian\ndepends linearly on lambda, which is NOT normally the case.\n", bl->prev->lambda, bl->lambda);
1191 /* normal delta H */
1194 gmx_fatal(FARGS, "Could not find a set for foreign lambda = %f\nin the files for lambda = %f", bl->lambda, bl->prev->lambda);
1198 gmx_fatal(FARGS, "Could not find a set for foreign lambda = %f\nin the files for lambda = %f", bl->prev->lambda, bl->lambda);
1210 /* estimate the maximum discretization error */
1211 static double barres_list_max_disc_err(barres_t *res, int nres)
1214 double disc_err = 0.;
1215 double delta_lambda;
1217 for (i = 0; i < nres; i++)
1219 barres_t *br = &(res[i]);
1221 delta_lambda = lambda_vec_abs_diff(br->b->native_lambda,
1222 br->a->native_lambda);
1224 for (j = 0; j < br->a->nsamples; j++)
1226 if (br->a->s[j]->hist)
1229 if (br->a->s[j]->derivative)
1231 Wfac = delta_lambda;
1234 disc_err = max(disc_err, Wfac*br->a->s[j]->hist->dx[0]);
1237 for (j = 0; j < br->b->nsamples; j++)
1239 if (br->b->s[j]->hist)
1242 if (br->b->s[j]->derivative)
1244 Wfac = delta_lambda;
1246 disc_err = max(disc_err, Wfac*br->b->s[j]->hist->dx[0]);
1254 /* impose start and end times on a sample collection, updating sample_ranges */
1255 static void sample_coll_impose_times(sample_coll_t *sc, double begin_t,
1259 for (i = 0; i < sc->nsamples; i++)
1261 samples_t *s = sc->s[i];
1262 sample_range_t *r = &(sc->r[i]);
1265 double end_time = s->hist->delta_time*s->hist->sum +
1266 s->hist->start_time;
1267 if (s->hist->start_time < begin_t || end_time > end_t)
1277 if (s->start_time < begin_t)
1279 r->start = (int)((begin_t - s->start_time)/s->delta_time);
1281 end_time = s->delta_time*s->ndu + s->start_time;
1282 if (end_time > end_t)
1284 r->end = (int)((end_t - s->start_time)/s->delta_time);
1290 for (j = 0; j < s->ndu; j++)
1292 if (s->t[j] < begin_t)
1297 if (s->t[j] >= end_t)
1304 if (r->start > r->end)
1310 sample_coll_calc_ntot(sc);
1313 static void sim_data_impose_times(sim_data_t *sd, double begin, double end)
1315 double first_t, last_t;
1316 double begin_t, end_t;
1318 lambda_data_t *head = sd->lb;
1321 if (begin <= 0 && end < 0)
1326 /* first determine the global start and end times */
1332 sample_coll_t *sc = lc->sc->next;
1333 while (sc != lc->sc)
1335 for (j = 0; j < sc->nsamples; j++)
1337 double start_t, end_t;
1339 start_t = sc->s[j]->start_time;
1340 end_t = sc->s[j]->start_time;
1343 end_t += sc->s[j]->delta_time*sc->s[j]->hist->sum;
1349 end_t = sc->s[j]->t[sc->s[j]->ndu-1];
1353 end_t += sc->s[j]->delta_time*sc->s[j]->ndu;
1357 if (start_t < first_t || first_t < 0)
1371 /* calculate the actual times */
1389 printf("\n Samples in time interval: %.3f - %.3f\n", first_t, last_t);
1391 if (begin_t > end_t)
1395 printf("Removing samples outside of: %.3f - %.3f\n", begin_t, end_t);
1397 /* then impose them */
1401 sample_coll_t *sc = lc->sc->next;
1402 while (sc != lc->sc)
1404 sample_coll_impose_times(sc, begin_t, end_t);
1412 /* create subsample i out of ni from an existing sample_coll */
1413 static gmx_bool sample_coll_create_subsample(sample_coll_t *sc,
1414 sample_coll_t *sc_orig,
1418 int hist_start, hist_end;
1420 gmx_int64_t ntot_start;
1421 gmx_int64_t ntot_end;
1422 gmx_int64_t ntot_so_far;
1424 *sc = *sc_orig; /* just copy all fields */
1426 /* allocate proprietary memory */
1427 snew(sc->s, sc_orig->nsamples);
1428 snew(sc->r, sc_orig->nsamples);
1430 /* copy the samples */
1431 for (j = 0; j < sc_orig->nsamples; j++)
1433 sc->s[j] = sc_orig->s[j];
1434 sc->r[j] = sc_orig->r[j]; /* copy the ranges too */
1437 /* now fix start and end fields */
1438 /* the casts avoid possible overflows */
1439 ntot_start = (gmx_int64_t)(sc_orig->ntot*(double)i/(double)ni);
1440 ntot_end = (gmx_int64_t)(sc_orig->ntot*(double)(i+1)/(double)ni);
1442 for (j = 0; j < sc->nsamples; j++)
1444 gmx_int64_t ntot_add;
1445 gmx_int64_t new_start, new_end;
1451 ntot_add = sc->s[j]->hist->sum;
1455 ntot_add = sc->r[j].end - sc->r[j].start;
1463 if (!sc->s[j]->hist)
1465 if (ntot_so_far < ntot_start)
1467 /* adjust starting point */
1468 new_start = sc->r[j].start + (ntot_start - ntot_so_far);
1472 new_start = sc->r[j].start;
1474 /* adjust end point */
1475 new_end = sc->r[j].start + (ntot_end - ntot_so_far);
1476 if (new_end > sc->r[j].end)
1478 new_end = sc->r[j].end;
1481 /* check if we're in range at all */
1482 if ( (new_end < new_start) || (new_start > sc->r[j].end) )
1487 /* and write the new range */
1488 sc->r[j].start = (int)new_start;
1489 sc->r[j].end = (int)new_end;
1496 double ntot_start_norm, ntot_end_norm;
1497 /* calculate the amount of overlap of the
1498 desired range (ntot_start -- ntot_end) onto
1499 the histogram range (ntot_so_far -- ntot_so_far+ntot_add)*/
1501 /* first calculate normalized bounds
1502 (where 0 is the start of the hist range, and 1 the end) */
1503 ntot_start_norm = (ntot_start-ntot_so_far)/(double)ntot_add;
1504 ntot_end_norm = (ntot_end-ntot_so_far)/(double)ntot_add;
1506 /* now fix the boundaries */
1507 ntot_start_norm = min(1, max(0., ntot_start_norm));
1508 ntot_end_norm = max(0, min(1., ntot_end_norm));
1510 /* and calculate the overlap */
1511 overlap = ntot_end_norm - ntot_start_norm;
1513 if (overlap > 0.95) /* we allow for 5% slack */
1515 sc->r[j].use = TRUE;
1517 else if (overlap < 0.05)
1519 sc->r[j].use = FALSE;
1527 ntot_so_far += ntot_add;
1529 sample_coll_calc_ntot(sc);
1534 /* calculate minimum and maximum work values in sample collection */
1535 static void sample_coll_min_max(sample_coll_t *sc, double Wfac,
1536 double *Wmin, double *Wmax)
1543 for (i = 0; i < sc->nsamples; i++)
1545 samples_t *s = sc->s[i];
1546 sample_range_t *r = &(sc->r[i]);
1551 for (j = r->start; j < r->end; j++)
1553 *Wmin = min(*Wmin, s->du[j]*Wfac);
1554 *Wmax = max(*Wmax, s->du[j]*Wfac);
1559 int hd = 0; /* determine the histogram direction: */
1561 if ( (s->hist->nhist > 1) && (Wfac < 0) )
1565 dx = s->hist->dx[hd];
1567 for (j = s->hist->nbin[hd]-1; j >= 0; j--)
1569 *Wmin = min(*Wmin, Wfac*(s->hist->x0[hd])*dx);
1570 *Wmax = max(*Wmax, Wfac*(s->hist->x0[hd])*dx);
1571 /* look for the highest value bin with values */
1572 if (s->hist->bin[hd][j] > 0)
1574 *Wmin = min(*Wmin, Wfac*(j+s->hist->x0[hd]+1)*dx);
1575 *Wmax = max(*Wmax, Wfac*(j+s->hist->x0[hd]+1)*dx);
1584 /* Initialize a sim_data structure */
1585 static void sim_data_init(sim_data_t *sd)
1587 /* make linked list */
1588 sd->lb = &(sd->lb_head);
1589 sd->lb->next = sd->lb;
1590 sd->lb->prev = sd->lb;
1592 lambda_components_init(&(sd->lc));
1596 static double calc_bar_sum(int n, const double *W, double Wfac, double sbMmDG)
1603 for (i = 0; i < n; i++)
1605 sum += 1./(1. + exp(Wfac*W[i] + sbMmDG));
1611 /* calculate the BAR average given a histogram
1613 if type== 0, calculate the best estimate for the average,
1614 if type==-1, calculate the minimum possible value given the histogram
1615 if type== 1, calculate the maximum possible value given the histogram */
1616 static double calc_bar_sum_hist(const hist_t *hist, double Wfac, double sbMmDG,
1622 /* normalization factor multiplied with bin width and
1623 number of samples (we normalize through M): */
1625 int hd = 0; /* determine the histogram direction: */
1628 if ( (hist->nhist > 1) && (Wfac < 0) )
1633 max = hist->nbin[hd]-1;
1636 max = hist->nbin[hd]; /* we also add whatever was out of range */
1639 for (i = 0; i < max; i++)
1641 double x = Wfac*((i+hist->x0[hd])+0.5)*dx; /* bin middle */
1642 double pxdx = hist->bin[0][i]*normdx; /* p(x)dx */
1644 sum += pxdx/(1. + exp(x + sbMmDG));
1650 static double calc_bar_lowlevel(sample_coll_t *ca, sample_coll_t *cb,
1651 double temp, double tol, int type)
1656 double Wfac1, Wfac2, Wmin, Wmax;
1657 double DG0, DG1, DG2, dDG1;
1659 double n1, n2; /* numbers of samples as doubles */
1664 /* count the numbers of samples */
1670 /*if (!lambda_vec_same(ca->native_lambda, ca->foreign_lambda))*/
1671 if (ca->foreign_lambda->dhdl < 0)
1673 /* this is the case when the delta U were calculated directly
1674 (i.e. we're not scaling dhdl) */
1680 /* we're using dhdl, so delta_lambda needs to be a
1681 multiplication factor. */
1682 /*double delta_lambda=cb->native_lambda-ca->native_lambda;*/
1683 double delta_lambda = lambda_vec_abs_diff(cb->native_lambda,
1685 if (cb->native_lambda->lc->N > 1)
1688 "Can't (yet) do multi-component dhdl interpolation");
1691 Wfac1 = beta*delta_lambda;
1692 Wfac2 = -beta*delta_lambda;
1697 /* We print the output both in kT and kJ/mol.
1698 * Here we determine DG in kT, so when beta < 1
1699 * the precision has to be increased.
1704 /* Calculate minimum and maximum work to give an initial estimate of
1705 * delta G as their average.
1708 double Wmin1, Wmin2, Wmax1, Wmax2;
1709 sample_coll_min_max(ca, Wfac1, &Wmin1, &Wmax1);
1710 sample_coll_min_max(cb, Wfac2, &Wmin2, &Wmax2);
1712 Wmin = min(Wmin1, Wmin2);
1713 Wmax = max(Wmax1, Wmax2);
1721 fprintf(debug, "DG %9.5f %9.5f\n", DG0, DG2);
1723 /* We approximate by bisection: given our initial estimates
1724 we keep checking whether the halfway point is greater or
1725 smaller than what we get out of the BAR averages.
1727 For the comparison we can use twice the tolerance. */
1728 while (DG2 - DG0 > 2*tol)
1730 DG1 = 0.5*(DG0 + DG2);
1732 /* calculate the BAR averages */
1735 for (i = 0; i < ca->nsamples; i++)
1737 samples_t *s = ca->s[i];
1738 sample_range_t *r = &(ca->r[i]);
1743 dDG1 += calc_bar_sum_hist(s->hist, Wfac1, (M-DG1), type);
1747 dDG1 += calc_bar_sum(r->end - r->start, s->du + r->start,
1752 for (i = 0; i < cb->nsamples; i++)
1754 samples_t *s = cb->s[i];
1755 sample_range_t *r = &(cb->r[i]);
1760 dDG1 -= calc_bar_sum_hist(s->hist, Wfac2, -(M-DG1), type);
1764 dDG1 -= calc_bar_sum(r->end - r->start, s->du + r->start,
1780 fprintf(debug, "DG %9.5f %9.5f\n", DG0, DG2);
1784 return 0.5*(DG0 + DG2);
1787 static void calc_rel_entropy(sample_coll_t *ca, sample_coll_t *cb,
1788 double temp, double dg, double *sa, double *sb)
1794 double Wfac1, Wfac2;
1800 /* count the numbers of samples */
1804 /* to ensure the work values are the same as during the delta_G */
1805 /*if (!lambda_vec_same(ca->native_lambda, ca->foreign_lambda))*/
1806 if (ca->foreign_lambda->dhdl < 0)
1808 /* this is the case when the delta U were calculated directly
1809 (i.e. we're not scaling dhdl) */
1815 /* we're using dhdl, so delta_lambda needs to be a
1816 multiplication factor. */
1817 double delta_lambda = lambda_vec_abs_diff(cb->native_lambda,
1819 Wfac1 = beta*delta_lambda;
1820 Wfac2 = -beta*delta_lambda;
1823 /* first calculate the average work in both directions */
1824 for (i = 0; i < ca->nsamples; i++)
1826 samples_t *s = ca->s[i];
1827 sample_range_t *r = &(ca->r[i]);
1832 for (j = r->start; j < r->end; j++)
1834 W_ab += Wfac1*s->du[j];
1839 /* normalization factor multiplied with bin width and
1840 number of samples (we normalize through M): */
1843 int hd = 0; /* histogram direction */
1844 if ( (s->hist->nhist > 1) && (Wfac1 < 0) )
1848 dx = s->hist->dx[hd];
1850 for (j = 0; j < s->hist->nbin[0]; j++)
1852 double x = Wfac1*((j+s->hist->x0[0])+0.5)*dx; /*bin ctr*/
1853 double pxdx = s->hist->bin[0][j]*normdx; /* p(x)dx */
1861 for (i = 0; i < cb->nsamples; i++)
1863 samples_t *s = cb->s[i];
1864 sample_range_t *r = &(cb->r[i]);
1869 for (j = r->start; j < r->end; j++)
1871 W_ba += Wfac1*s->du[j];
1876 /* normalization factor multiplied with bin width and
1877 number of samples (we normalize through M): */
1880 int hd = 0; /* histogram direction */
1881 if ( (s->hist->nhist > 1) && (Wfac2 < 0) )
1885 dx = s->hist->dx[hd];
1887 for (j = 0; j < s->hist->nbin[0]; j++)
1889 double x = Wfac1*((j+s->hist->x0[0])+0.5)*dx; /*bin ctr*/
1890 double pxdx = s->hist->bin[0][j]*normdx; /* p(x)dx */
1898 /* then calculate the relative entropies */
1903 static void calc_dg_stddev(sample_coll_t *ca, sample_coll_t *cb,
1904 double temp, double dg, double *stddev)
1908 double sigmafact = 0.;
1910 double Wfac1, Wfac2;
1916 /* count the numbers of samples */
1920 /* to ensure the work values are the same as during the delta_G */
1921 /*if (!lambda_vec_same(ca->native_lambda, ca->foreign_lambda))*/
1922 if (ca->foreign_lambda->dhdl < 0)
1924 /* this is the case when the delta U were calculated directly
1925 (i.e. we're not scaling dhdl) */
1931 /* we're using dhdl, so delta_lambda needs to be a
1932 multiplication factor. */
1933 double delta_lambda = lambda_vec_abs_diff(cb->native_lambda,
1935 Wfac1 = beta*delta_lambda;
1936 Wfac2 = -beta*delta_lambda;
1942 /* calculate average in both directions */
1943 for (i = 0; i < ca->nsamples; i++)
1945 samples_t *s = ca->s[i];
1946 sample_range_t *r = &(ca->r[i]);
1951 for (j = r->start; j < r->end; j++)
1953 sigmafact += 1./(2. + 2.*cosh((M + Wfac1*s->du[j] - dg)));
1958 /* normalization factor multiplied with bin width and
1959 number of samples (we normalize through M): */
1962 int hd = 0; /* histogram direction */
1963 if ( (s->hist->nhist > 1) && (Wfac1 < 0) )
1967 dx = s->hist->dx[hd];
1969 for (j = 0; j < s->hist->nbin[0]; j++)
1971 double x = Wfac1*((j+s->hist->x0[0])+0.5)*dx; /*bin ctr*/
1972 double pxdx = s->hist->bin[0][j]*normdx; /* p(x)dx */
1974 sigmafact += pxdx/(2. + 2.*cosh((M + x - dg)));
1979 for (i = 0; i < cb->nsamples; i++)
1981 samples_t *s = cb->s[i];
1982 sample_range_t *r = &(cb->r[i]);
1987 for (j = r->start; j < r->end; j++)
1989 sigmafact += 1./(2. + 2.*cosh((M - Wfac2*s->du[j] - dg)));
1994 /* normalization factor multiplied with bin width and
1995 number of samples (we normalize through M): */
1998 int hd = 0; /* histogram direction */
1999 if ( (s->hist->nhist > 1) && (Wfac2 < 0) )
2003 dx = s->hist->dx[hd];
2005 for (j = 0; j < s->hist->nbin[0]; j++)
2007 double x = Wfac2*((j+s->hist->x0[0])+0.5)*dx; /*bin ctr*/
2008 double pxdx = s->hist->bin[0][j]*normdx; /* p(x)dx */
2010 sigmafact += pxdx/(2. + 2.*cosh((M - x - dg)));
2016 sigmafact /= (n1 + n2);
2020 Shirts, Bair, Hooker & Pande, Phys. Rev. Lett 91, 140601 (2003): */
2021 *stddev = sqrt(((1./sigmafact) - ( (n1+n2)/n1 + (n1+n2)/n2 )));
2026 static void calc_bar(barres_t *br, double tol,
2027 int npee_min, int npee_max, gmx_bool *bEE,
2031 double dg_sig2, sa_sig2, sb_sig2, stddev_sig2; /* intermediate variance values
2032 for calculated quantities */
2033 int nsample1, nsample2;
2034 double temp = br->a->temp;
2036 double dg_min, dg_max;
2037 gmx_bool have_hist = FALSE;
2039 br->dg = calc_bar_lowlevel(br->a, br->b, temp, tol, 0);
2041 br->dg_disc_err = 0.;
2042 br->dg_histrange_err = 0.;
2044 /* check if there are histograms */
2045 for (i = 0; i < br->a->nsamples; i++)
2047 if (br->a->r[i].use && br->a->s[i]->hist)
2055 for (i = 0; i < br->b->nsamples; i++)
2057 if (br->b->r[i].use && br->b->s[i]->hist)
2065 /* calculate histogram-specific errors */
2068 dg_min = calc_bar_lowlevel(br->a, br->b, temp, tol, -1);
2069 dg_max = calc_bar_lowlevel(br->a, br->b, temp, tol, 1);
2071 if (fabs(dg_max - dg_min) > GMX_REAL_EPS*10)
2073 /* the histogram range error is the biggest of the differences
2074 between the best estimate and the extremes */
2075 br->dg_histrange_err = fabs(dg_max - dg_min);
2077 br->dg_disc_err = 0.;
2078 for (i = 0; i < br->a->nsamples; i++)
2080 if (br->a->s[i]->hist)
2082 br->dg_disc_err = max(br->dg_disc_err, br->a->s[i]->hist->dx[0]);
2085 for (i = 0; i < br->b->nsamples; i++)
2087 if (br->b->s[i]->hist)
2089 br->dg_disc_err = max(br->dg_disc_err, br->b->s[i]->hist->dx[0]);
2093 calc_rel_entropy(br->a, br->b, temp, br->dg, &(br->sa), &(br->sb));
2095 calc_dg_stddev(br->a, br->b, temp, br->dg, &(br->dg_stddev) );
2104 sample_coll_t ca, cb;
2106 /* initialize the samples */
2107 sample_coll_init(&ca, br->a->native_lambda, br->a->foreign_lambda,
2109 sample_coll_init(&cb, br->b->native_lambda, br->b->foreign_lambda,
2112 for (npee = npee_min; npee <= npee_max; npee++)
2121 double dstddev2 = 0;
2124 for (p = 0; p < npee; p++)
2131 cac = sample_coll_create_subsample(&ca, br->a, p, npee);
2132 cbc = sample_coll_create_subsample(&cb, br->b, p, npee);
2136 printf("WARNING: histogram number incompatible with block number for averaging: can't do error estimate\n");
2140 sample_coll_destroy(&ca);
2144 sample_coll_destroy(&cb);
2149 dgp = calc_bar_lowlevel(&ca, &cb, temp, tol, 0);
2153 partsum[npee*(npee_max+1)+p] += dgp;
2155 calc_rel_entropy(&ca, &cb, temp, dgp, &sac, &sbc);
2160 calc_dg_stddev(&ca, &cb, temp, dgp, &stddevc );
2163 dstddev2 += stddevc*stddevc;
2165 sample_coll_destroy(&ca);
2166 sample_coll_destroy(&cb);
2170 dg_sig2 += (dgs2-dgs*dgs)/(npee-1);
2176 sa_sig2 += (dsa2-dsa*dsa)/(npee-1);
2177 sb_sig2 += (dsb2-dsb*dsb)/(npee-1);
2181 stddev_sig2 += (dstddev2-dstddev*dstddev)/(npee-1);
2183 br->dg_err = sqrt(dg_sig2/(npee_max - npee_min + 1));
2184 br->sa_err = sqrt(sa_sig2/(npee_max - npee_min + 1));
2185 br->sb_err = sqrt(sb_sig2/(npee_max - npee_min + 1));
2186 br->dg_stddev_err = sqrt(stddev_sig2/(npee_max - npee_min + 1));
2191 static double bar_err(int nbmin, int nbmax, const double *partsum)
2194 double svar, s, s2, dg;
2197 for (nb = nbmin; nb <= nbmax; nb++)
2201 for (b = 0; b < nb; b++)
2203 dg = partsum[nb*(nbmax+1)+b];
2209 svar += (s2 - s*s)/(nb - 1);
2212 return sqrt(svar/(nbmax + 1 - nbmin));
2216 /* Seek the end of an identifier (consecutive non-spaces), followed by
2217 an optional number of spaces or '='-signs. Returns a pointer to the
2218 first non-space value found after that. Returns NULL if the string
2221 static const char *find_value(const char *str)
2223 gmx_bool name_end_found = FALSE;
2225 /* if the string is a NULL pointer, return a NULL pointer. */
2230 while (*str != '\0')
2232 /* first find the end of the name */
2233 if (!name_end_found)
2235 if (isspace(*str) || (*str == '=') )
2237 name_end_found = TRUE;
2242 if (!( isspace(*str) || (*str == '=') ))
2254 /* read a vector-notation description of a lambda vector */
2255 static gmx_bool read_lambda_compvec(const char *str,
2257 const lambda_components_t *lc_in,
2258 lambda_components_t *lc_out,
2262 gmx_bool initialize_lc = FALSE; /* whether to initialize the lambda
2263 components, or to check them */
2264 gmx_bool start_reached = FALSE; /* whether the start of component names
2266 gmx_bool vector = FALSE; /* whether there are multiple components */
2267 int n = 0; /* current component number */
2268 const char *val_start = NULL; /* start of the component name, or NULL
2269 if not in a value */
2279 if (lc_out && lc_out->N == 0)
2281 initialize_lc = TRUE;
2296 start_reached = TRUE;
2299 else if (*str == '(')
2302 start_reached = TRUE;
2304 else if (!isspace(*str))
2306 gmx_fatal(FARGS, "Error in lambda components in %s", fn);
2313 if (isspace(*str) || *str == ')' || *str == ',' || *str == '\0')
2320 lambda_components_add(lc_out, val_start,
2325 if (!lambda_components_check(lc_out, n, val_start,
2334 /* add a vector component to lv */
2335 lv->val[n] = strtod(val_start, &strtod_end);
2336 if (val_start == strtod_end)
2339 "Error reading lambda vector in %s", fn);
2342 /* reset for the next identifier */
2351 else if (isalnum(*str))
2364 gmx_fatal(FARGS, "Error in lambda components in %s", fn);
2372 else if (lv == NULL)
2378 gmx_fatal(FARGS, "Incomplete lambda vector data in %s",
2398 gmx_fatal(FARGS, "Incomplete lambda components data in %s", fn);
2404 /* read and check the component names from a string */
2405 static gmx_bool read_lambda_components(const char *str,
2406 lambda_components_t *lc,
2410 return read_lambda_compvec(str, NULL, NULL, lc, end, fn);
2413 /* read an initialized lambda vector from a string */
2414 static gmx_bool read_lambda_vector(const char *str,
2419 return read_lambda_compvec(str, lv, lv->lc, NULL, end, fn);
2424 /* deduce lambda value from legend.
2426 legend = the legend string
2428 lam = the initialized lambda vector
2429 returns whether to use the data in this set.
2431 static gmx_bool legend2lambda(const char *fn,
2436 const char *ptr = NULL, *ptr2 = NULL;
2437 gmx_bool ok = FALSE;
2438 gmx_bool bdhdl = FALSE;
2439 const char *tostr = " to ";
2443 gmx_fatal(FARGS, "There is no legend in file '%s', can not deduce lambda", fn);
2446 /* look for the last 'to': */
2450 ptr2 = strstr(ptr2, tostr);
2457 while (ptr2 != NULL && *ptr2 != '\0');
2461 ptr += strlen(tostr)-1; /* and advance past that 'to' */
2465 /* look for the = sign */
2466 ptr = strrchr(legend, '=');
2469 /* otherwise look for the last space */
2470 ptr = strrchr(legend, ' ');
2474 if (strstr(legend, "dH"))
2479 else if (strchr(legend, 'D') != NULL && strchr(legend, 'H') != NULL)
2484 else /*if (strstr(legend, "pV"))*/
2495 gmx_fatal(FARGS, "There is no proper lambda legend in file '%s', can not deduce lambda", fn);
2499 ptr = find_value(ptr);
2500 if (!ptr || !read_lambda_vector(ptr, lam, NULL, fn))
2502 gmx_fatal(FARGS, "lambda vector '%s' %s faulty", legend, fn);
2511 ptr = strrchr(legend, '=');
2515 /* there must be a component name */
2519 gmx_fatal(FARGS, "dhdl legend '%s' %s faulty", legend, fn);
2521 /* now backtrack to the start of the identifier */
2522 while (isspace(*ptr))
2528 gmx_fatal(FARGS, "dhdl legend '%s' %s faulty", legend, fn);
2531 while (!isspace(*ptr))
2536 gmx_fatal(FARGS, "dhdl legend '%s' %s faulty", legend, fn);
2540 strncpy(buf, ptr, (end-ptr));
2541 buf[(end-ptr)] = '\0';
2542 dhdl_index = lambda_components_find(lam->lc, ptr, (end-ptr));
2546 strncpy(buf, ptr, (end-ptr));
2547 buf[(end-ptr)] = '\0';
2549 "Did not find lambda component for '%s' in %s",
2558 "dhdl without component name with >1 lambda component in %s",
2563 lam->dhdl = dhdl_index;
2568 static gmx_bool subtitle2lambda(const char *subtitle, xvg_t *ba, const char *fn,
2569 lambda_components_t *lc)
2574 double native_lambda;
2578 /* first check for a state string */
2579 ptr = strstr(subtitle, "state");
2583 const char *val_end;
2585 /* the new 4.6 style lambda vectors */
2586 ptr = find_value(ptr);
2589 index = strtol(ptr, &end, 10);
2592 gmx_fatal(FARGS, "Incomplete state data in %s", fn);
2599 gmx_fatal(FARGS, "Incomplete state data in %s", fn);
2602 /* now find the lambda vector component names */
2603 while (*ptr != '(' && !isalnum(*ptr))
2609 "Incomplete lambda vector component data in %s", fn);
2614 if (!read_lambda_components(ptr, lc, &val_end, fn))
2617 "lambda vector components in %s don't match those previously read",
2620 ptr = find_value(val_end);
2623 gmx_fatal(FARGS, "Incomplete state data in %s", fn);
2626 lambda_vec_init(&(ba->native_lambda), lc);
2627 if (!read_lambda_vector(ptr, &(ba->native_lambda), NULL, fn))
2629 gmx_fatal(FARGS, "lambda vector in %s faulty", fn);
2631 ba->native_lambda.index = index;
2636 /* compatibility mode: check for lambda in other ways. */
2637 /* plain text lambda string */
2638 ptr = strstr(subtitle, "lambda");
2641 /* xmgrace formatted lambda string */
2642 ptr = strstr(subtitle, "\\xl\\f{}");
2646 /* xmgr formatted lambda string */
2647 ptr = strstr(subtitle, "\\8l\\4");
2651 ptr = strstr(ptr, "=");
2655 bFound = (sscanf(ptr+1, "%lf", &(native_lambda)) == 1);
2656 /* add the lambda component name as an empty string */
2659 if (!lambda_components_check(lc, 0, "", 0))
2662 "lambda vector components in %s don't match those previously read",
2668 lambda_components_add(lc, "", 0);
2670 lambda_vec_init(&(ba->native_lambda), lc);
2671 ba->native_lambda.val[0] = native_lambda;
2678 static void filename2lambda(const char *fn)
2681 const char *ptr, *digitptr;
2685 /* go to the end of the path string and search backward to find the last
2686 directory in the path which has to contain the value of lambda
2688 while (ptr[1] != '\0')
2692 /* searching backward to find the second directory separator */
2697 if (ptr[0] != DIR_SEPARATOR && ptr[1] == DIR_SEPARATOR)
2705 /* save the last position of a digit between the last two
2706 separators = in the last dirname */
2707 if (dirsep > 0 && isdigit(*ptr))
2715 gmx_fatal(FARGS, "While trying to read the lambda value from the file path:"
2716 " last directory in the path '%s' does not contain a number", fn);
2718 if (digitptr[-1] == '-')
2722 lambda = strtod(digitptr, &endptr);
2723 if (endptr == digitptr)
2725 gmx_fatal(FARGS, "Malformed number in file path '%s'", fn);
2729 static void read_bar_xvg_lowlevel(const char *fn, real *temp, xvg_t *ba,
2730 lambda_components_t *lc)
2733 char *subtitle, **legend, *ptr;
2735 gmx_bool native_lambda_read = FALSE;
2743 np = read_xvg_legend(fn, &ba->y, &ba->nset, &subtitle, &legend);
2746 gmx_fatal(FARGS, "File %s contains no usable data.", fn);
2748 /* Reorder the data */
2750 for (i = 1; i < ba->nset; i++)
2752 ba->y[i-1] = ba->y[i];
2756 snew(ba->np, ba->nset);
2757 for (i = 0; i < ba->nset; i++)
2763 if (subtitle != NULL)
2765 /* try to extract temperature */
2766 ptr = strstr(subtitle, "T =");
2770 if (sscanf(ptr, "%lf", &ba->temp) == 1)
2774 gmx_fatal(FARGS, "Found temperature of %f in file '%s'",
2784 gmx_fatal(FARGS, "Did not find a temperature in the subtitle in file '%s', use the -temp option of [TT]gmx bar[tt]", fn);
2789 /* Try to deduce lambda from the subtitle */
2792 if (subtitle2lambda(subtitle, ba, fn, lc))
2794 native_lambda_read = TRUE;
2797 snew(ba->lambda, ba->nset);
2800 /* Check if we have a single set, no legend, nset=1 means t and dH/dl */
2803 if (!native_lambda_read)
2805 /* Deduce lambda from the file name */
2806 filename2lambda(fn);
2807 native_lambda_read = TRUE;
2809 ba->lambda[0] = ba->native_lambda;
2813 gmx_fatal(FARGS, "File %s contains multiple sets but no legends, can not determine the lambda values", fn);
2818 for (i = 0; i < ba->nset; )
2820 gmx_bool use = FALSE;
2821 /* Read lambda from the legend */
2822 lambda_vec_init( &(ba->lambda[i]), lc );
2823 lambda_vec_copy( &(ba->lambda[i]), &(ba->native_lambda));
2824 use = legend2lambda(fn, legend[i], &(ba->lambda[i]));
2827 lambda_vec_print(&(ba->lambda[i]), buf, FALSE);
2833 printf("%s: Ignoring set '%s'.\n", fn, legend[i]);
2834 for (j = i+1; j < ba->nset; j++)
2836 ba->y[j-1] = ba->y[j];
2837 legend[j-1] = legend[j];
2844 if (!native_lambda_read)
2846 gmx_fatal(FARGS, "File %s contains multiple sets but no indication of the native lambda", fn);
2851 for (i = 0; i < ba->nset-1; i++)
2859 static void read_bar_xvg(char *fn, real *temp, sim_data_t *sd)
2868 read_bar_xvg_lowlevel(fn, temp, barsim, &(sd->lc));
2870 if (barsim->nset < 1)
2872 gmx_fatal(FARGS, "File '%s' contains fewer than two columns", fn);
2875 if (!gmx_within_tol(*temp, barsim->temp, GMX_FLOAT_EPS) && (*temp > 0) )
2877 gmx_fatal(FARGS, "Temperature in file %s different from earlier files or setting\n", fn);
2879 *temp = barsim->temp;
2881 /* now create a series of samples_t */
2882 snew(s, barsim->nset);
2883 for (i = 0; i < barsim->nset; i++)
2885 samples_init(s+i, &(barsim->native_lambda), &(barsim->lambda[i]),
2886 barsim->temp, lambda_vec_same(&(barsim->native_lambda),
2887 &(barsim->lambda[i])),
2889 s[i].du = barsim->y[i];
2890 s[i].ndu = barsim->np[i];
2893 lambda_data_list_insert_sample(sd->lb, s+i);
2898 lambda_vec_print(s[0].native_lambda, buf, FALSE);
2899 printf("%s: %.1f - %.1f; lambda = %s\n dH/dl & foreign lambdas:\n",
2900 fn, s[0].t[0], s[0].t[s[0].ndu-1], buf);
2901 for (i = 0; i < barsim->nset; i++)
2903 lambda_vec_print(s[i].foreign_lambda, buf, TRUE);
2904 printf(" %s (%d pts)\n", buf, s[i].ndu);
2910 static void read_edr_rawdh_block(samples_t **smp, int *ndu, t_enxblock *blk,
2911 double start_time, double delta_time,
2912 lambda_vec_t *native_lambda, double temp,
2913 double *last_t, const char *filename)
2917 double old_foreign_lambda;
2918 lambda_vec_t *foreign_lambda;
2920 samples_t *s; /* convenience pointer */
2923 /* check the block types etc. */
2924 if ( (blk->nsub < 3) ||
2925 (blk->sub[0].type != xdr_datatype_int) ||
2926 (blk->sub[1].type != xdr_datatype_double) ||
2928 (blk->sub[2].type != xdr_datatype_float) &&
2929 (blk->sub[2].type != xdr_datatype_double)
2931 (blk->sub[0].nr < 1) ||
2932 (blk->sub[1].nr < 1) )
2935 "Unexpected/corrupted block data in file %s around time %f.",
2936 filename, start_time);
2939 snew(foreign_lambda, 1);
2940 lambda_vec_init(foreign_lambda, native_lambda->lc);
2941 lambda_vec_copy(foreign_lambda, native_lambda);
2942 type = blk->sub[0].ival[0];
2945 for (i = 0; i < native_lambda->lc->N; i++)
2947 foreign_lambda->val[i] = blk->sub[1].dval[i];
2952 if (blk->sub[0].nr > 1)
2954 foreign_lambda->dhdl = blk->sub[0].ival[1];
2958 foreign_lambda->dhdl = 0;
2964 /* initialize the samples structure if it's empty. */
2966 samples_init(*smp, native_lambda, foreign_lambda, temp,
2967 type == dhbtDHDL, filename);
2968 (*smp)->start_time = start_time;
2969 (*smp)->delta_time = delta_time;
2972 /* set convenience pointer */
2975 /* now double check */
2976 if (!lambda_vec_same(s->foreign_lambda, foreign_lambda) )
2978 char buf[STRLEN], buf2[STRLEN];
2979 lambda_vec_print(foreign_lambda, buf, FALSE);
2980 lambda_vec_print(s->foreign_lambda, buf2, FALSE);
2981 fprintf(stderr, "Got foreign lambda=%s, expected: %s\n", buf, buf2);
2982 gmx_fatal(FARGS, "Corrupted data in file %s around t=%f.",
2983 filename, start_time);
2986 /* make room for the data */
2987 if (s->ndu_alloc < (size_t)(s->ndu + blk->sub[2].nr) )
2989 s->ndu_alloc += (s->ndu_alloc < (size_t)blk->sub[2].nr) ?
2990 blk->sub[2].nr*2 : s->ndu_alloc;
2991 srenew(s->du_alloc, s->ndu_alloc);
2992 s->du = s->du_alloc;
2995 s->ndu += blk->sub[2].nr;
2996 s->ntot += blk->sub[2].nr;
2997 *ndu = blk->sub[2].nr;
2999 /* and copy the data*/
3000 for (j = 0; j < blk->sub[2].nr; j++)
3002 if (blk->sub[2].type == xdr_datatype_float)
3004 s->du[startj+j] = blk->sub[2].fval[j];
3008 s->du[startj+j] = blk->sub[2].dval[j];
3011 if (start_time + blk->sub[2].nr*delta_time > *last_t)
3013 *last_t = start_time + blk->sub[2].nr*delta_time;
3017 static samples_t *read_edr_hist_block(int *nsamples, t_enxblock *blk,
3018 double start_time, double delta_time,
3019 lambda_vec_t *native_lambda, double temp,
3020 double *last_t, const char *filename)
3025 double old_foreign_lambda;
3026 lambda_vec_t *foreign_lambda;
3030 /* check the block types etc. */
3031 if ( (blk->nsub < 2) ||
3032 (blk->sub[0].type != xdr_datatype_double) ||
3033 (blk->sub[1].type != xdr_datatype_int64) ||
3034 (blk->sub[0].nr < 2) ||
3035 (blk->sub[1].nr < 2) )
3038 "Unexpected/corrupted block data in file %s around time %f",
3039 filename, start_time);
3042 nhist = blk->nsub-2;
3050 "Unexpected/corrupted block data in file %s around time %f",
3051 filename, start_time);
3057 snew(foreign_lambda, 1);
3058 lambda_vec_init(foreign_lambda, native_lambda->lc);
3059 lambda_vec_copy(foreign_lambda, native_lambda);
3060 type = (int)(blk->sub[1].lval[1]);
3063 double old_foreign_lambda;
3065 old_foreign_lambda = blk->sub[0].dval[0];
3066 if (old_foreign_lambda >= 0)
3068 foreign_lambda->val[0] = old_foreign_lambda;
3069 if (foreign_lambda->lc->N > 1)
3072 "Single-component lambda in multi-component file %s",
3078 for (i = 0; i < native_lambda->lc->N; i++)
3080 foreign_lambda->val[i] = blk->sub[0].dval[i+2];
3086 if (foreign_lambda->lc->N > 1)
3088 if (blk->sub[1].nr < 3 + nhist)
3091 "Missing derivative coord in multi-component file %s",
3094 foreign_lambda->dhdl = blk->sub[1].lval[2 + nhist];
3098 foreign_lambda->dhdl = 0;
3102 samples_init(s, native_lambda, foreign_lambda, temp, type == dhbtDHDL,
3106 for (i = 0; i < nhist; i++)
3108 nbins[i] = blk->sub[i+2].nr;
3111 hist_init(s->hist, nhist, nbins);
3113 for (i = 0; i < nhist; i++)
3115 s->hist->x0[i] = blk->sub[1].lval[2+i];
3116 s->hist->dx[i] = blk->sub[0].dval[1];
3119 s->hist->dx[i] = -s->hist->dx[i];
3123 s->hist->start_time = start_time;
3124 s->hist->delta_time = delta_time;
3125 s->start_time = start_time;
3126 s->delta_time = delta_time;
3128 for (i = 0; i < nhist; i++)
3131 gmx_int64_t sum = 0;
3133 for (j = 0; j < s->hist->nbin[i]; j++)
3135 int binv = (int)(blk->sub[i+2].ival[j]);
3137 s->hist->bin[i][j] = binv;
3150 gmx_fatal(FARGS, "Histogram counts don't match in %s",
3156 if (start_time + s->hist->sum*delta_time > *last_t)
3158 *last_t = start_time + s->hist->sum*delta_time;
3164 static void read_barsim_edr(char *fn, real *temp, sim_data_t *sd)
3170 gmx_enxnm_t *enm = NULL;
3171 double first_t = -1;
3173 samples_t **samples_rawdh = NULL; /* contains samples for raw delta_h */
3174 int *nhists = NULL; /* array to keep count & print at end */
3175 int *npts = NULL; /* array to keep count & print at end */
3176 lambda_vec_t **lambdas = NULL; /* array to keep count & print at end */
3177 lambda_vec_t *native_lambda;
3178 double end_time; /* the end time of the last batch of samples */
3180 lambda_vec_t start_lambda;
3182 fp = open_enx(fn, "r");
3183 do_enxnms(fp, &nre, &enm);
3186 snew(native_lambda, 1);
3187 start_lambda.lc = NULL;
3189 while (do_enx(fp, fr))
3191 /* count the data blocks */
3192 int nblocks_raw = 0;
3193 int nblocks_hist = 0;
3196 /* DHCOLL block information: */
3197 double start_time = 0, delta_time = 0, old_start_lambda = 0, delta_lambda = 0;
3200 /* count the blocks and handle collection information: */
3201 for (i = 0; i < fr->nblock; i++)
3203 if (fr->block[i].id == enxDHHIST)
3207 if (fr->block[i].id == enxDH)
3211 if (fr->block[i].id == enxDHCOLL)
3214 if ( (fr->block[i].nsub < 1) ||
3215 (fr->block[i].sub[0].type != xdr_datatype_double) ||
3216 (fr->block[i].sub[0].nr < 5))
3218 gmx_fatal(FARGS, "Unexpected block data in file %s", fn);
3221 /* read the data from the DHCOLL block */
3222 rtemp = fr->block[i].sub[0].dval[0];
3223 start_time = fr->block[i].sub[0].dval[1];
3224 delta_time = fr->block[i].sub[0].dval[2];
3225 old_start_lambda = fr->block[i].sub[0].dval[3];
3226 delta_lambda = fr->block[i].sub[0].dval[4];
3228 if (delta_lambda > 0)
3230 gmx_fatal(FARGS, "Lambda values not constant in %s: can't apply BAR method", fn);
3232 if ( ( *temp != rtemp) && (*temp > 0) )
3234 gmx_fatal(FARGS, "Temperature in file %s different from earlier files or setting\n", fn);
3238 if (old_start_lambda >= 0)
3242 if (!lambda_components_check(&(sd->lc), 0, "", 0))
3245 "lambda vector components in %s don't match those previously read",
3251 lambda_components_add(&(sd->lc), "", 0);
3253 if (!start_lambda.lc)
3255 lambda_vec_init(&start_lambda, &(sd->lc));
3257 start_lambda.val[0] = old_start_lambda;
3261 /* read lambda vector */
3263 gmx_bool check = (sd->lc.N > 0);
3264 if (fr->block[i].nsub < 2)
3267 "No lambda vector, but start_lambda=%f\n",
3270 n_lambda_vec = fr->block[i].sub[1].ival[1];
3271 for (j = 0; j < n_lambda_vec; j++)
3274 efpt_singular_names[fr->block[i].sub[1].ival[1+j]];
3277 /* check the components */
3278 lambda_components_check(&(sd->lc), j, name,
3283 lambda_components_add(&(sd->lc), name,
3287 lambda_vec_init(&start_lambda, &(sd->lc));
3288 start_lambda.index = fr->block[i].sub[1].ival[0];
3289 for (j = 0; j < n_lambda_vec; j++)
3291 start_lambda.val[j] = fr->block[i].sub[0].dval[5+j];
3296 first_t = start_time;
3303 gmx_fatal(FARGS, "Did not find delta H information in file %s", fn);
3305 if (nblocks_raw > 0 && nblocks_hist > 0)
3307 gmx_fatal(FARGS, "Can't handle both raw delta U data and histograms in the same file %s", fn);
3312 /* check the native lambda */
3313 if (!lambda_vec_same(&start_lambda, native_lambda) )
3315 gmx_fatal(FARGS, "Native lambda not constant in file %s: started at %f, and becomes %f at time %f",
3316 fn, native_lambda, start_lambda, start_time);
3318 /* check the number of samples against the previous number */
3319 if ( ((nblocks_raw+nblocks_hist) != nsamples) || (nlam != 1 ) )
3321 gmx_fatal(FARGS, "Unexpected block count in %s: was %d, now %d\n",
3322 fn, nsamples+1, nblocks_raw+nblocks_hist+nlam);
3324 /* check whether last iterations's end time matches with
3325 the currrent start time */
3326 if ( (fabs(last_t - start_time) > 2*delta_time) && last_t >= 0)
3328 /* it didn't. We need to store our samples and reallocate */
3329 for (i = 0; i < nsamples; i++)
3331 if (samples_rawdh[i])
3333 /* insert it into the existing list */
3334 lambda_data_list_insert_sample(sd->lb,
3336 /* and make sure we'll allocate a new one this time
3338 samples_rawdh[i] = NULL;
3345 /* this is the first round; allocate the associated data
3347 /*native_lambda=start_lambda;*/
3348 lambda_vec_init(native_lambda, &(sd->lc));
3349 lambda_vec_copy(native_lambda, &start_lambda);
3350 nsamples = nblocks_raw+nblocks_hist;
3351 snew(nhists, nsamples);
3352 snew(npts, nsamples);
3353 snew(lambdas, nsamples);
3354 snew(samples_rawdh, nsamples);
3355 for (i = 0; i < nsamples; i++)
3360 samples_rawdh[i] = NULL; /* init to NULL so we know which
3361 ones contain values */
3366 k = 0; /* counter for the lambdas, etc. arrays */
3367 for (i = 0; i < fr->nblock; i++)
3369 if (fr->block[i].id == enxDH)
3371 int type = (fr->block[i].sub[0].ival[0]);
3372 if (type == dhbtDH || type == dhbtDHDL)
3375 read_edr_rawdh_block(&(samples_rawdh[k]),
3378 start_time, delta_time,
3379 native_lambda, rtemp,
3382 if (samples_rawdh[k])
3384 lambdas[k] = samples_rawdh[k]->foreign_lambda;
3389 else if (fr->block[i].id == enxDHHIST)
3391 int type = (int)(fr->block[i].sub[1].lval[1]);
3392 if (type == dhbtDH || type == dhbtDHDL)
3396 samples_t *s; /* this is where the data will go */
3397 s = read_edr_hist_block(&nb, &(fr->block[i]),
3398 start_time, delta_time,
3399 native_lambda, rtemp,
3404 lambdas[k] = s->foreign_lambda;
3407 /* and insert the new sample immediately */
3408 for (j = 0; j < nb; j++)
3410 lambda_data_list_insert_sample(sd->lb, s+j);
3416 /* Now store all our extant sample collections */
3417 for (i = 0; i < nsamples; i++)
3419 if (samples_rawdh[i])
3421 /* insert it into the existing list */
3422 lambda_data_list_insert_sample(sd->lb, samples_rawdh[i]);
3430 lambda_vec_print(native_lambda, buf, FALSE);
3431 printf("%s: %.1f - %.1f; lambda = %s\n foreign lambdas:\n",
3432 fn, first_t, last_t, buf);
3433 for (i = 0; i < nsamples; i++)
3437 lambda_vec_print(lambdas[i], buf, TRUE);
3440 printf(" %s (%d hists)\n", buf, nhists[i]);
3444 printf(" %s (%d pts)\n", buf, npts[i]);
3456 int gmx_bar(int argc, char *argv[])
3458 static const char *desc[] = {
3459 "[THISMODULE] calculates free energy difference estimates through ",
3460 "Bennett's acceptance ratio method (BAR). It also automatically",
3461 "adds series of individual free energies obtained with BAR into",
3462 "a combined free energy estimate.[PAR]",
3464 "Every individual BAR free energy difference relies on two ",
3465 "simulations at different states: say state A and state B, as",
3466 "controlled by a parameter, [GRK]lambda[grk] (see the [TT].mdp[tt] parameter",
3467 "[TT]init_lambda[tt]). The BAR method calculates a ratio of weighted",
3468 "average of the Hamiltonian difference of state B given state A and",
3470 "The energy differences to the other state must be calculated",
3471 "explicitly during the simulation. This can be done with",
3472 "the [TT].mdp[tt] option [TT]foreign_lambda[tt].[PAR]",
3474 "Input option [TT]-f[tt] expects multiple [TT]dhdl.xvg[tt] files. ",
3475 "Two types of input files are supported:[BR]",
3476 "[TT]*[tt] Files with more than one [IT]y[it]-value. ",
3477 "The files should have columns ",
3478 "with dH/d[GRK]lambda[grk] and [GRK]Delta[grk][GRK]lambda[grk]. ",
3479 "The [GRK]lambda[grk] values are inferred ",
3480 "from the legends: [GRK]lambda[grk] of the simulation from the legend of ",
3481 "dH/d[GRK]lambda[grk] and the foreign [GRK]lambda[grk] values from the ",
3482 "legends of Delta H",
3484 "[TT]*[tt] Files with only one [IT]y[it]-value. Using the",
3485 "[TT]-extp[tt] option for these files, it is assumed",
3486 "that the [IT]y[it]-value is dH/d[GRK]lambda[grk] and that the ",
3487 "Hamiltonian depends linearly on [GRK]lambda[grk]. ",
3488 "The [GRK]lambda[grk] value of the simulation is inferred from the ",
3489 "subtitle (if present), otherwise from a number in the subdirectory ",
3490 "in the file name.[PAR]",
3492 "The [GRK]lambda[grk] of the simulation is parsed from ",
3493 "[TT]dhdl.xvg[tt] file's legend containing the string 'dH', the ",
3494 "foreign [GRK]lambda[grk] values from the legend containing the ",
3495 "capitalized letters 'D' and 'H'. The temperature is parsed from ",
3496 "the legend line containing 'T ='.[PAR]",
3498 "The input option [TT]-g[tt] expects multiple [TT].edr[tt] files. ",
3499 "These can contain either lists of energy differences (see the ",
3500 "[TT].mdp[tt] option [TT]separate_dhdl_file[tt]), or a series of ",
3501 "histograms (see the [TT].mdp[tt] options [TT]dh_hist_size[tt] and ",
3502 "[TT]dh_hist_spacing[tt]).", "The temperature and [GRK]lambda[grk] ",
3503 "values are automatically deduced from the [TT]ener.edr[tt] file.[PAR]",
3505 "In addition to the [TT].mdp[tt] option [TT]foreign_lambda[tt], ",
3506 "the energy difference can also be extrapolated from the ",
3507 "dH/d[GRK]lambda[grk] values. This is done with the[TT]-extp[tt]",
3508 "option, which assumes that the system's Hamiltonian depends linearly",
3509 "on [GRK]lambda[grk], which is not normally the case.[PAR]",
3511 "The free energy estimates are determined using BAR with bisection, ",
3512 "with the precision of the output set with [TT]-prec[tt]. ",
3513 "An error estimate taking into account time correlations ",
3514 "is made by splitting the data into blocks and determining ",
3515 "the free energy differences over those blocks and assuming ",
3516 "the blocks are independent. ",
3517 "The final error estimate is determined from the average variance ",
3518 "over 5 blocks. A range of block numbers for error estimation can ",
3519 "be provided with the options [TT]-nbmin[tt] and [TT]-nbmax[tt].[PAR]",
3521 "[THISMODULE] tries to aggregate samples with the same 'native' and ",
3522 "'foreign' [GRK]lambda[grk] values, but always assumes independent ",
3523 "samples. [BB]Note[bb] that when aggregating energy ",
3524 "differences/derivatives with different sampling intervals, this is ",
3525 "almost certainly not correct. Usually subsequent energies are ",
3526 "correlated and different time intervals mean different degrees ",
3527 "of correlation between samples.[PAR]",
3529 "The results are split in two parts: the last part contains the final ",
3530 "results in kJ/mol, together with the error estimate for each part ",
3531 "and the total. The first part contains detailed free energy ",
3532 "difference estimates and phase space overlap measures in units of ",
3533 "kT (together with their computed error estimate). The printed ",
3535 "[TT]*[tt] lam_A: the [GRK]lambda[grk] values for point A.[BR]",
3536 "[TT]*[tt] lam_B: the [GRK]lambda[grk] values for point B.[BR]",
3537 "[TT]*[tt] DG: the free energy estimate.[BR]",
3538 "[TT]*[tt] s_A: an estimate of the relative entropy of B in A.[BR]",
3539 "[TT]*[tt] s_B: an estimate of the relative entropy of A in B.[BR]",
3540 "[TT]*[tt] stdev: an estimate expected per-sample standard deviation.[PAR]",
3542 "The relative entropy of both states in each other's ensemble can be ",
3543 "interpreted as a measure of phase space overlap: ",
3544 "the relative entropy s_A of the work samples of lambda_B in the ",
3545 "ensemble of lambda_A (and vice versa for s_B), is a ",
3546 "measure of the 'distance' between Boltzmann distributions of ",
3547 "the two states, that goes to zero for identical distributions. See ",
3548 "Wu & Kofke, J. Chem. Phys. 123 084109 (2005) for more information.",
3550 "The estimate of the expected per-sample standard deviation, as given ",
3551 "in Bennett's original BAR paper: Bennett, J. Comp. Phys. 22, p 245 (1976).",
3552 "Eq. 10 therein gives an estimate of the quality of sampling (not directly",
3553 "of the actual statistical error, because it assumes independent samples).[PAR]",
3555 "To get a visual estimate of the phase space overlap, use the ",
3556 "[TT]-oh[tt] option to write series of histograms, together with the ",
3557 "[TT]-nbin[tt] option.[PAR]"
3559 static real begin = 0, end = -1, temp = -1;
3560 int nd = 2, nbmin = 5, nbmax = 5;
3562 gmx_bool use_dhdl = FALSE;
3563 gmx_bool calc_s, calc_v;
3565 { "-b", FALSE, etREAL, {&begin}, "Begin time for BAR" },
3566 { "-e", FALSE, etREAL, {&end}, "End time for BAR" },
3567 { "-temp", FALSE, etREAL, {&temp}, "Temperature (K)" },
3568 { "-prec", FALSE, etINT, {&nd}, "The number of digits after the decimal point" },
3569 { "-nbmin", FALSE, etINT, {&nbmin}, "Minimum number of blocks for error estimation" },
3570 { "-nbmax", FALSE, etINT, {&nbmax}, "Maximum number of blocks for error estimation" },
3571 { "-nbin", FALSE, etINT, {&nbin}, "Number of bins for histogram output"},
3572 { "-extp", FALSE, etBOOL, {&use_dhdl}, "Whether to linearly extrapolate dH/dl values to use as energies"}
3576 { efXVG, "-f", "dhdl", ffOPTRDMULT },
3577 { efEDR, "-g", "ener", ffOPTRDMULT },
3578 { efXVG, "-o", "bar", ffOPTWR },
3579 { efXVG, "-oi", "barint", ffOPTWR },
3580 { efXVG, "-oh", "histogram", ffOPTWR }
3582 #define NFILE asize(fnm)
3585 int nf = 0; /* file counter */
3587 int nfile_tot; /* total number of input files */
3592 sim_data_t sim_data; /* the simulation data */
3593 barres_t *results; /* the results */
3594 int nresults; /* number of results in results array */
3597 double prec, dg_tot, dg, sig, dg_tot_max, dg_tot_min;
3599 char dgformat[20], xvg2format[STRLEN], xvg3format[STRLEN];
3600 char buf[STRLEN], buf2[STRLEN];
3601 char ktformat[STRLEN], sktformat[STRLEN];
3602 char kteformat[STRLEN], skteformat[STRLEN];
3605 gmx_bool result_OK = TRUE, bEE = TRUE;
3607 gmx_bool disc_err = FALSE;
3608 double sum_disc_err = 0.; /* discretization error */
3609 gmx_bool histrange_err = FALSE;
3610 double sum_histrange_err = 0.; /* histogram range error */
3611 double stat_err = 0.; /* statistical error */
3613 if (!parse_common_args(&argc, argv,
3615 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
3620 if (opt2bSet("-f", NFILE, fnm))
3622 nxvgfile = opt2fns(&fxvgnms, "-f", NFILE, fnm);
3624 if (opt2bSet("-g", NFILE, fnm))
3626 nedrfile = opt2fns(&fedrnms, "-g", NFILE, fnm);
3629 sim_data_init(&sim_data);
3631 /* make linked list */
3633 lambda_data_init(lb, 0, 0);
3639 nfile_tot = nxvgfile + nedrfile;
3643 gmx_fatal(FARGS, "No input files!");
3648 gmx_fatal(FARGS, "Can not have negative number of digits");
3650 prec = pow(10, -nd);
3652 snew(partsum, (nbmax+1)*(nbmax+1));
3655 /* read in all files. First xvg files */
3656 for (f = 0; f < nxvgfile; f++)
3658 read_bar_xvg(fxvgnms[f], &temp, &sim_data);
3661 /* then .edr files */
3662 for (f = 0; f < nedrfile; f++)
3664 read_barsim_edr(fedrnms[f], &temp, &sim_data);;
3668 /* fix the times to allow for equilibration */
3669 sim_data_impose_times(&sim_data, begin, end);
3671 if (opt2bSet("-oh", NFILE, fnm))
3673 sim_data_histogram(&sim_data, opt2fn("-oh", NFILE, fnm), nbin, oenv);
3676 /* assemble the output structures from the lambdas */
3677 results = barres_list_create(&sim_data, &nresults, use_dhdl);
3679 sum_disc_err = barres_list_max_disc_err(results, nresults);
3683 printf("\nNo results to calculate.\n");
3687 if (sum_disc_err > prec)
3689 prec = sum_disc_err;
3690 nd = ceil(-log10(prec));
3691 printf("WARNING: setting the precision to %g because that is the minimum\n reasonable number, given the expected discretization error.\n", prec);
3695 /*sprintf(lamformat,"%%6.3f");*/
3696 sprintf( dgformat, "%%%d.%df", 3+nd, nd);
3697 /* the format strings of the results in kT */
3698 sprintf( ktformat, "%%%d.%df", 5+nd, nd);
3699 sprintf( sktformat, "%%%ds", 6+nd);
3700 /* the format strings of the errors in kT */
3701 sprintf( kteformat, "%%%d.%df", 3+nd, nd);
3702 sprintf( skteformat, "%%%ds", 4+nd);
3703 sprintf(xvg2format, "%s %s\n", "%s", dgformat);
3704 sprintf(xvg3format, "%s %s %s\n", "%s", dgformat, dgformat);
3709 if (opt2bSet("-o", NFILE, fnm))
3711 sprintf(buf, "%s (%s)", "\\DeltaG", "kT");
3712 fpb = xvgropen_type(opt2fn("-o", NFILE, fnm), "Free energy differences",
3713 "\\lambda", buf, exvggtXYDY, oenv);
3717 if (opt2bSet("-oi", NFILE, fnm))
3719 sprintf(buf, "%s (%s)", "\\DeltaG", "kT");
3720 fpi = xvgropen(opt2fn("-oi", NFILE, fnm), "Free energy integral",
3721 "\\lambda", buf, oenv);
3731 /* first calculate results */
3734 for (f = 0; f < nresults; f++)
3736 /* Determine the free energy difference with a factor of 10
3737 * more accuracy than requested for printing.
3739 calc_bar(&(results[f]), 0.1*prec, nbmin, nbmax,
3742 if (results[f].dg_disc_err > prec/10.)
3746 if (results[f].dg_histrange_err > prec/10.)
3748 histrange_err = TRUE;
3752 /* print results in kT */
3756 printf("\nTemperature: %g K\n", temp);
3758 printf("\nDetailed results in kT (see help for explanation):\n\n");
3759 printf("%6s ", " lam_A");
3760 printf("%6s ", " lam_B");
3761 printf(sktformat, "DG ");
3764 printf(skteformat, "+/- ");
3768 printf(skteformat, "disc ");
3772 printf(skteformat, "range ");
3774 printf(sktformat, "s_A ");
3777 printf(skteformat, "+/- " );
3779 printf(sktformat, "s_B ");
3782 printf(skteformat, "+/- " );
3784 printf(sktformat, "stdev ");
3787 printf(skteformat, "+/- ");
3790 for (f = 0; f < nresults; f++)
3792 lambda_vec_print_short(results[f].a->native_lambda, buf);
3794 lambda_vec_print_short(results[f].b->native_lambda, buf);
3796 printf(ktformat, results[f].dg);
3800 printf(kteformat, results[f].dg_err);
3805 printf(kteformat, results[f].dg_disc_err);
3810 printf(kteformat, results[f].dg_histrange_err);
3813 printf(ktformat, results[f].sa);
3817 printf(kteformat, results[f].sa_err);
3820 printf(ktformat, results[f].sb);
3824 printf(kteformat, results[f].sb_err);
3827 printf(ktformat, results[f].dg_stddev);
3831 printf(kteformat, results[f].dg_stddev_err);
3835 /* Check for negative relative entropy with a 95% certainty. */
3836 if (results[f].sa < -2*results[f].sa_err ||
3837 results[f].sb < -2*results[f].sb_err)
3845 printf("\nWARNING: Some of these results violate the Second Law of "
3846 "Thermodynamics: \n"
3847 " This is can be the result of severe undersampling, or "
3849 " there is something wrong with the simulations.\n");
3853 /* final results in kJ/mol */
3854 printf("\n\nFinal results in kJ/mol:\n\n");
3856 for (f = 0; f < nresults; f++)
3861 lambda_vec_print_short(results[f].a->native_lambda, buf);
3862 fprintf(fpi, xvg2format, buf, dg_tot);
3868 lambda_vec_print_intermediate(results[f].a->native_lambda,
3869 results[f].b->native_lambda,
3872 fprintf(fpb, xvg3format, buf, results[f].dg, results[f].dg_err);
3876 lambda_vec_print_short(results[f].a->native_lambda, buf);
3877 lambda_vec_print_short(results[f].b->native_lambda, buf2);
3878 printf("%s - %s", buf, buf2);
3881 printf(dgformat, results[f].dg*kT);
3885 printf(dgformat, results[f].dg_err*kT);
3889 printf(" (max. range err. = ");
3890 printf(dgformat, results[f].dg_histrange_err*kT);
3892 sum_histrange_err += results[f].dg_histrange_err*kT;
3896 dg_tot += results[f].dg;
3900 lambda_vec_print_short(results[0].a->native_lambda, buf);
3901 lambda_vec_print_short(results[nresults-1].b->native_lambda, buf2);
3902 printf("%s - %s", buf, buf2);
3905 printf(dgformat, dg_tot*kT);
3908 stat_err = bar_err(nbmin, nbmax, partsum)*kT;
3910 printf(dgformat, max(max(stat_err, sum_disc_err), sum_histrange_err));
3915 printf("\nmaximum discretization error = ");
3916 printf(dgformat, sum_disc_err);
3917 if (bEE && stat_err < sum_disc_err)
3919 printf("WARNING: discretization error (%g) is larger than statistical error.\n Decrease histogram spacing for more accurate results\n", stat_err);
3924 printf("\nmaximum histogram range error = ");
3925 printf(dgformat, sum_histrange_err);
3926 if (bEE && stat_err < sum_histrange_err)
3928 printf("WARNING: histogram range error (%g) is larger than statistical error.\n Increase histogram range for more accurate results\n", stat_err);
3937 lambda_vec_print_short(results[nresults-1].b->native_lambda, buf);
3938 fprintf(fpi, xvg2format, buf, dg_tot);
3942 do_view(oenv, opt2fn_null("-o", NFILE, fnm), "-xydy");
3943 do_view(oenv, opt2fn_null("-oi", NFILE, fnm), "-xydy");