Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_bar.c
blob501c717b8c0f49dfa4349ef54fe156e4e2fbc4a8
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Green Red Orange Magenta Azure Cyan Skyblue
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39 #include <math.h>
40 #include <string.h>
41 #include <ctype.h>
42 #include <math.h>
44 #include "sysstuff.h"
45 #include "typedefs.h"
46 #include "smalloc.h"
47 #include "futil.h"
48 #include "statutil.h"
49 #include "copyrite.h"
50 #include "macros.h"
51 #include "physics.h"
52 #include "gmx_fatal.h"
53 #include "xvgr.h"
54 #include "gmx_ana.h"
55 #include "maths.h"
57 typedef struct {
58 char *filename;
59 int nset;
60 int np;
61 int begin;
62 int end;
63 double temp;
64 double *lambda;
65 double *t;
66 double **y;
67 } barsim_t;
69 /* calculated values */
70 typedef struct {
71 barsim_t *a, *b; /* the simulation data */
73 double lambda_a, lambda_b; /* the lambda values at a and b */
75 double dg; /* the free energy difference */
76 double dg_err; /* the free energy difference */
78 double sa; /* relative entropy of b in state a */
79 double sa_err; /* error in sa */
80 double sb; /* relative entropy of a in state b */
81 double sb_err; /* error in sb */
83 double dg_stddev; /* expected dg stddev per sample */
84 double dg_stddev_err; /* error in dg_stddev */
85 } barres_t;
88 static double calc_bar_sum(int n,double *W,double Wfac,double sbMmDG)
90 int i;
91 double sum;
93 sum = 0;
95 for(i=0; i<n; i++)
97 sum += 1./(1. + exp(Wfac*W[i] + sbMmDG));
100 return sum;
103 static double calc_bar_lowlevel(int n1,double *W1,int n2,double *W2,
104 double delta_lambda,double temp,double tol)
106 double kT,beta,beta_dl,M;
107 double DG;
108 int i;
109 double Wfac1,Wfac2,Wmin,Wmax;
110 double DG0,DG1,DG2,dDG1;
111 double sum1,sum2;
113 kT = BOLTZ*temp;
114 beta = 1/kT;
116 M = log((double)n1/(double)n2);
118 if (delta_lambda == 0)
120 Wfac1 = beta;
121 Wfac2 = beta;
123 else
125 Wfac1 = beta*delta_lambda;
126 Wfac2 = -beta*delta_lambda;
129 if (beta < 1)
131 /* We print the output both in kT and kJ/mol.
132 * Here we determine DG in kT, so when beta < 1
133 * the precision has to be increased.
135 tol *= beta;
138 /* Calculate minimum and maximum work to give an initial estimate of
139 * delta G as their average.
141 Wmin = W1[0];
142 Wmax = W1[0];
143 for(i=0; i<n1; i++)
145 Wmin = min(Wmin,W1[i]*Wfac1);
146 Wmax = max(Wmax,W1[i]*Wfac1);
148 for(i=0; i<n2; i++)
150 Wmin = min(Wmin,-W2[i]*Wfac2);
151 Wmax = max(Wmax,-W2[i]*Wfac2);
153 DG0 = Wmin;
154 DG2 = Wmax;
156 /* For the comparison we can use twice the tolerance. */
157 if (debug)
159 fprintf(debug,"DG %9.5f %9.5f\n",DG0,DG2);
161 while (DG2 - DG0 > 2*tol)
163 DG1 = 0.5*(DG0 + DG2);
165 /*printf("Wfac1=%g, Wfac2=%g, beta=%g, DG1=%g\n",Wfac1,Wfac2,beta,
166 DG1);*/
167 dDG1 =
168 calc_bar_sum(n1,W1,Wfac1, (M-DG1)) -
169 calc_bar_sum(n2,W2,Wfac2,-(M-DG1));
171 if (dDG1 < 0)
173 DG0 = DG1;
175 else
177 DG2 = DG1;
179 if (debug)
181 fprintf(debug,"DG %9.5f %9.5f\n",DG0,DG2);
185 return 0.5*(DG0 + DG2);
188 static void calc_rel_entropy(int n1,double *W1,int n2,double *W2,
189 double delta_lambda, double temp,
190 double dg, double *sa, double *sb)
192 int i;
193 double W_ab=0.;
194 double W_ba=0.;
195 double kT, beta;
196 double Wfac1, Wfac2;
198 kT = BOLTZ*temp;
199 beta = 1/kT;
201 /* to ensure the work values are the same as during the delta_G */
202 if (delta_lambda == 0)
204 Wfac1 = beta;
205 Wfac2 = beta;
207 else
209 Wfac1 = beta*delta_lambda;
210 Wfac2 = -beta*delta_lambda;
213 /* first calculate the average work in both directions */
214 for(i=0;i<n1;i++)
216 W_ab += Wfac1*W1[i];
218 W_ab/=n1;
219 for(i=0;i<n2;i++)
221 W_ba += Wfac2*W2[i];
223 W_ba/=n2;
225 /* then calculate the relative entropies */
226 *sa = (W_ab - dg);
227 *sb = (W_ba + dg);
230 static void calc_dg_stddev(int n1, double *W1, int n2, double *W2,
231 double delta_lambda, double temp,
232 double dg, double *stddev)
234 int i;
235 double M;
236 double sigmafact=0.;
237 double kT, beta;
238 double Wfac1, Wfac2;
240 double nn1=n1; /* this makes the fraction in the *stddev eq below nicer */
241 double nn2=n2;
243 kT = BOLTZ*temp;
244 beta = 1/kT;
246 /* to ensure the work values are the same as during the delta_G */
247 if (delta_lambda == 0)
249 Wfac1 = beta;
250 Wfac2 = beta;
252 else
254 Wfac1 = beta*delta_lambda;
255 Wfac2 = -beta*delta_lambda;
258 M = log(nn1/nn2);
260 /* calculate average in both directions */
261 for(i=0;i<n1;i++)
263 sigmafact += 1./(2. + 2.*cosh((M + Wfac1*W1[i] - dg)));
265 for(i=0;i<n2;i++)
267 sigmafact += 1./(2. + 2.*cosh((M - Wfac2*W2[i] - dg)));
269 sigmafact /= (n1 + n2);
271 /* Eq. 10 from
272 Shirts, Bair, Hooker & Pande, Phys. Rev. Lett 91, 140601 (2003): */
273 *stddev = sqrt(((1./sigmafact) - ( (nn1+nn2)/nn1 + (nn1+nn2)/nn2 )));
280 static void get_begin_end(barsim_t *ba,real begin,real end,int *b,int *e)
282 int i;
284 i = 0;
285 while (i + 1 < ba->np && ba->t[i] < begin)
287 i++;
289 if (i >= ba->np)
291 gmx_fatal(FARGS,"Some data end before the start time %g",begin);
293 *b = i;
295 i = ba->np;
296 if (end >= begin)
298 while (i > *b && ba->t[i-1] > end)
300 i--;
303 *e = i;
306 static int get_lam_set(barsim_t *ba,double lambda)
308 int i;
310 i = 1;
311 while (i < ba->nset &&
312 !gmx_within_tol(ba->lambda[i],lambda,10*GMX_REAL_EPS))
314 i++;
316 if (i == ba->nset)
318 gmx_fatal(FARGS,"Could not find a set for lambda = %g in the file '%s' of lambda = %g",lambda,ba->filename,ba->lambda[0]);
321 return i;
324 static void calc_bar(barsim_t *ba1,barsim_t *ba2,bool bUsedhdl,
325 double tol, int npee_min,int npee_max,
326 barres_t *br, bool *bEE, double *partsum)
328 int np1,np2,s1,s2,npee,p;
329 double delta_lambda;
330 double dg_sig2,sa_sig2,sb_sig2,stddev_sig2; /* intermediate variance values
331 for calculated quantities */
332 br->a = ba1;
333 br->b = ba2;
334 br->lambda_a = ba1->lambda[0];
335 br->lambda_b = ba2->lambda[0];
337 if (bUsedhdl)
339 s1 = 0;
340 s2 = 0;
342 delta_lambda = ba2->lambda[0] - ba1->lambda[0];
344 else
346 s1 = get_lam_set(ba1,ba2->lambda[0]);
347 s2 = get_lam_set(ba2,ba1->lambda[0]);
349 delta_lambda = 0;
352 np1 = ba1->end - ba1->begin;
353 np2 = ba2->end - ba2->begin;
355 br->dg = calc_bar_lowlevel(np1,ba1->y[s1]+ba1->begin,
356 np2,ba2->y[s2]+ba2->begin,
357 delta_lambda,ba1->temp,tol);
360 calc_rel_entropy(np1, ba1->y[s1]+ba1->begin,
361 np2, ba2->y[s2]+ba2->begin,
362 delta_lambda, ba1->temp, br->dg, &(br->sa), &(br->sb));
363 calc_dg_stddev(np1, ba1->y[s1]+ba1->begin,
364 np2, ba2->y[s2]+ba2->begin,
365 delta_lambda, ba1->temp, br->dg, &(br->dg_stddev) );
368 dg_sig2 = 0;
369 sa_sig2 = 0;
370 sb_sig2 = 0;
371 stddev_sig2 = 0;
372 if (np1 >= npee_max && np2 >= npee_max)
374 for(npee=npee_min; npee<=npee_max; npee++)
376 double dgs = 0;
377 double dgs2 = 0;
378 double dsa = 0;
379 double dsb = 0;
380 double dsa2 = 0;
381 double dsb2 = 0;
382 double dstddev = 0;
383 double dstddev2 = 0;
386 for(p=0; p<npee; p++)
388 double dgp;
389 double stddevc;
390 double sac, sbc;
391 dgp = calc_bar_lowlevel(np1/npee,
392 ba1->y[s1]+ba1->begin+p*(np1/npee),
393 np2/npee,
394 ba2->y[s2]+ba2->begin+p*(np2/npee),
395 delta_lambda,ba1->temp,tol);
396 dgs += dgp;
397 dgs2 += dgp*dgp;
399 partsum[npee*(npee_max+1)+p] += dgp;
401 calc_rel_entropy(np1/npee,
402 ba1->y[s1]+ba1->begin+p*(np1/npee),
403 np2/npee,
404 ba2->y[s2]+ba2->begin+p*(np2/npee),
405 delta_lambda, ba1->temp, dgp, &sac, &sbc);
406 dsa += sac;
407 dsa2 += sac*sac;
408 dsb += sbc;
409 dsb2 += sbc*sbc;
410 calc_dg_stddev(np1/npee,
411 ba1->y[s1]+ba1->begin+p*(np1/npee),
412 np2/npee,
413 ba2->y[s2]+ba2->begin+p*(np2/npee),
414 delta_lambda, ba1->temp, dgp, &stddevc );
416 dstddev += stddevc;
417 dstddev2 += stddevc*stddevc;
419 dgs /= npee;
420 dgs2 /= npee;
421 dg_sig2 += (dgs2-dgs*dgs)/(npee-1);
423 dsa /= npee;
424 dsa2 /= npee;
425 dsb /= npee;
426 dsb2 /= npee;
427 sa_sig2 += (dsa2-dsa*dsa)/(npee-1);
428 sb_sig2 += (dsb2-dsb*dsb)/(npee-1);
430 dstddev /= npee;
431 dstddev2 /= npee;
432 stddev_sig2 += (dstddev2-dstddev*dstddev)/(npee-1);
434 br->dg_err = sqrt(dg_sig2/(npee_max - npee_min + 1));
435 br->sa_err = sqrt(sa_sig2/(npee_max - npee_min + 1));
436 br->sb_err = sqrt(sb_sig2/(npee_max - npee_min + 1));
437 br->dg_stddev_err = sqrt(stddev_sig2/(npee_max - npee_min + 1));
439 else
441 *bEE = FALSE;
446 static double bar_err(int nbmin, int nbmax, const double *partsum)
448 int nb,b;
449 double svar,s,s2,dg;
451 svar = 0;
452 for(nb=nbmin; nb<=nbmax; nb++)
454 s = 0;
455 s2 = 0;
456 for(b=0; b<nb; b++)
458 dg = partsum[nb*(nbmax+1)+b];
459 s += dg;
460 s2 += dg*dg;
462 s /= nb;
463 s2 /= nb;
464 svar += (s2 - s*s)/(nb - 1);
467 return sqrt(svar/(nbmax + 1 - nbmin));
471 static double legend2lambda(char *fn,const char *legend,bool bdhdl)
473 double lambda=0;
474 const char *ptr;
476 if (legend == NULL)
478 gmx_fatal(FARGS,"There is no legend in file '%s', can not deduce lambda",fn);
480 ptr = strrchr(legend,' ');
481 if (( bdhdl && strstr(legend,"dH") == NULL) ||
482 (!bdhdl && (strchr(legend,'D') == NULL ||
483 strchr(legend,'H') == NULL)) ||
484 ptr == NULL)
486 gmx_fatal(FARGS,"There is no proper lambda legend in file '%s', can not deduce lambda",fn);
488 if (sscanf(ptr,"%lf",&lambda) != 1)
490 gmx_fatal(FARGS,"There is no proper lambda legend in file '%s', can not deduce lambda",fn);
493 return lambda;
496 static bool subtitle2lambda(const char *subtitle,double *lambda)
498 bool bFound;
499 char *ptr;
501 bFound = FALSE;
503 /* plain text lambda string */
504 ptr = strstr(subtitle,"lambda");
505 if (ptr == NULL)
507 /* xmgrace formatted lambda string */
508 ptr = strstr(subtitle,"\\xl\\f{}");
510 if (ptr == NULL)
512 /* xmgr formatted lambda string */
513 ptr = strstr(subtitle,"\\8l\\4");
515 if (ptr != NULL)
517 ptr = strstr(ptr,"=");
519 if (ptr != NULL)
521 bFound = (sscanf(ptr+1,"%lf",lambda) == 1);
524 return bFound;
527 static double filename2lambda(char *fn)
529 double lambda;
530 char *ptr,*endptr,*digitptr;
531 int dirsep;
532 ptr = fn;
533 /* go to the end of the path string and search backward to find the last
534 directory in the path which has to contain the value of lambda
536 while (ptr[1] != '\0')
538 ptr++;
540 /* searching backward to find the second directory separator */
541 dirsep = 0;
542 digitptr = NULL;
543 while (ptr >= fn)
545 if (ptr[0] != DIR_SEPARATOR && ptr[1] == DIR_SEPARATOR)
547 if (dirsep == 1) break;
548 dirsep++;
550 /* save the last position of a digit between the last two
551 separators = in the last dirname */
552 if (dirsep > 0 && isdigit(*ptr))
554 digitptr = ptr;
556 ptr--;
558 if (!digitptr)
560 gmx_fatal(FARGS,"While trying to read the lambda value from the file path:"
561 " last directory in the path '%s' does not contain a number",fn);
563 if (digitptr[-1] == '-')
565 digitptr--;
567 lambda = strtod(digitptr,&endptr);
568 if (endptr == digitptr)
570 gmx_fatal(FARGS,"Malformed number in file path '%s'",fn);
573 return lambda;
576 static void read_barsim(char *fn,double begin,double end,real temp,
577 barsim_t *ba)
579 int i;
580 char *subtitle,**legend,*ptr;
582 ba->filename = fn;
584 printf("'%s' ",ba->filename);
586 ba->np = read_xvg_legend(fn,&ba->y,&ba->nset,&subtitle,&legend);
587 if (!ba->y)
589 gmx_fatal(FARGS,"File %s contains no usable data.",fn);
591 ba->t = ba->y[0];
593 get_begin_end(ba,begin,end,&ba->begin,&ba->end);
594 printf("%.1f - %.1f, %6d points, lam:",
595 ba->t[ba->begin],ba->t[ba->end-1],ba->end-ba->begin);
597 ba->temp = -1;
598 if (subtitle != NULL)
600 ptr = strstr(subtitle,"T =");
601 if (ptr != NULL)
603 ptr += 3;
604 if (sscanf(ptr,"%lf",&ba->temp) == 1)
606 if (ba->temp <= 0)
608 gmx_fatal(FARGS,"Found temperature of %g in file '%s'",
609 ba->temp,fn);
614 if (ba->temp < 0)
616 if (temp <= 0)
618 gmx_fatal(FARGS,"Did not find a temperature in the subtitle in file '%s', use the -temp option of g_bar",fn);
620 ba->temp = temp;
623 snew(ba->lambda,ba->nset-1);
624 if (legend == NULL)
626 /* Check if we have a single set, nset=2 means t and dH/dl */
627 if (ba->nset == 2)
629 /* Try to deduce lambda from the subtitle */
630 if (subtitle != NULL &&
631 !subtitle2lambda(subtitle,&ba->lambda[0]))
633 /* Deduce lambda from the file name */
634 ba->lambda[0] = filename2lambda(fn);
636 printf(" %g",ba->lambda[0]);
638 else
640 gmx_fatal(FARGS,"File %s contains multiple sets but no legends, can not determine the lambda values",fn);
643 else
645 for(i=0; i<ba->nset-1; i++)
647 /* Read lambda from the legend */
648 ba->lambda[i] = legend2lambda(fn,legend[i],i==0);
649 printf(" %g",ba->lambda[i]);
652 printf("\n");
654 /* Reorder the data */
655 for(i=1; i<ba->nset; i++)
657 ba->y[i-1] = ba->y[i];
659 if (legend != NULL)
661 for(i=0; i<ba->nset-1; i++)
663 sfree(legend[i]);
665 sfree(legend);
667 ba->nset--;
670 int gmx_bar(int argc,char *argv[])
672 static const char *desc[] = {
673 "g_bar calculates free energy difference estimates through ",
674 "Bennett's acceptance ratio method. ",
675 "Input option [TT]-f[tt] expects multiple dhdl files. ",
676 "Two types of input files are supported:[BR]",
677 "* Files with only one y-value, for such files it is assumed ",
678 "that the y-value is dH/dlambda and that the Hamiltonian depends ",
679 "linearly on lambda. The lambda value of the simulation is inferred ",
680 "from the subtitle if present, otherwise from a number in the",
681 "subdirectory in the file name.",
682 "[BR]",
683 "* Files with more than one y-value. The files should have columns ",
684 "with dH/dlambda and Delta lambda. The lambda values are inferred ",
685 "from the legends: ",
686 "lambda of the simulation from the legend of dH/dlambda ",
687 "and the foreign lambda's from the legends of Delta H.[PAR]",
689 "The lambda of the simulation is parsed from dhdl.xvg file's legend ",
690 "containing the string 'dH', the foreign lambda's from the legend ",
691 "containing the capitalized letters 'D' and 'H'. The temperature ",
692 "is parsed from the legend line containing 'T ='.[PAR]",
694 "The free energy estimates are determined using BAR with bisection, ",
695 "the precision of the output is set with [TT]-prec[tt]. ",
696 "An error estimate taking into account time correlations ",
697 "is made by splitting the data into blocks and determining ",
698 "the free energy differences over those blocks and assuming ",
699 "the blocks are independent. ",
700 "The final error estimate is determined from the average variance ",
701 "over 5 blocks. A range of blocks numbers for error estimation can ",
702 "be provided with the options [TT]-nbmin[tt] and [TT]-nbmax[tt].[PAR]",
704 "The results are split in two parts: the last part contains the final ",
705 "results in kJ/mol, together with the error estimate for each part ",
706 "and the total. The first part contains detailed free energy ",
707 "difference estimates and phase space overlap measures in units of ",
708 "kT (together with their computed error estimate). The printed ",
709 "values are:[BR]",
710 "* lam_A: the lambda values for point A.[BR]",
711 "* lam_B: the lambda values for point B.[BR]",
712 "* DG: the free energy estimate.[BR]",
713 "* s_A: an estimate of the relative entropy of B in A.[BR]",
714 "* s_A: an estimate of the relative entropy of A in B.[BR]",
715 "* stdev: an estimate expected per-sample standard deviation.[PAR]",
717 "The relative entropy of both states in each other's ensemble can be ",
718 "interpreted as a measure of phase space overlap: ",
719 "the relative entropy s_A of the work samples of lambda_B in the ",
720 "ensemble of lambda_A (and vice versa for s_B), is a ",
721 "measure of the 'distance' between Boltzmann distributions of ",
722 "the two states, that goes to zero for identical distributions. See ",
723 "Wu & Kofke, J. Chem. Phys. 123 084109 (2009) for more information.",
724 "[PAR]",
725 "The estimate of the expected per-sample standard deviation, as given ",
726 "in Bennett's original BAR paper: ",
727 "Bennett, J. Comp. Phys. 22, p 245 (1976), Eq. 10 gives an estimate ",
728 "of the quality of sampling (not directly of the actual statistical ",
729 "error, because it assumes independent samples).[PAR]",
732 static real begin=0,end=-1,temp=-1;
733 static int nd=2,nbmin=5,nbmax=5;
734 bool calc_s,calc_v;
735 t_pargs pa[] = {
736 { "-b", FALSE, etREAL, {&begin}, "Begin time for BAR" },
737 { "-e", FALSE, etREAL, {&end}, "End time for BAR" },
738 { "-temp", FALSE, etREAL, {&temp}, "Temperature (K)" },
739 { "-prec", FALSE, etINT, {&nd}, "The number of digits after the decimal point" },
740 { "-nbmin", FALSE, etINT, {&nbmin}, "Minimum number of blocks for error estimation" },
741 { "-nbmax", FALSE, etINT, {&nbmax}, "Maximum number of blocks for error estimation" }
744 t_filenm fnm[] = {
745 { efXVG, "-f", "dhdl", ffRDMULT },
746 { efXVG, "-o", "bar", ffOPTWR },
747 { efXVG, "-oi", "barint", ffOPTWR }
749 #define NFILE asize(fnm)
751 int nfile,f,f2,fm,n1,nm;
752 char **fnms;
753 barsim_t *ba,ba_tmp;
754 barres_t *results;
755 double *partsum;
756 double prec,dg_tot,dg,sig;
757 FILE *fpb,*fpi;
758 char lamformat[20];
759 char dgformat[20],xvg2format[STRLEN],xvg3format[STRLEN],buf[STRLEN];
760 char ktformat[STRLEN], sktformat[STRLEN];
761 char kteformat[STRLEN], skteformat[STRLEN];
762 output_env_t oenv;
763 double kT, beta;
764 bool result_OK=TRUE,bEE=TRUE;
766 CopyRight(stderr,argv[0]);
767 parse_common_args(&argc,argv,
768 PCA_CAN_VIEW,
769 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
771 nfile = opt2fns(&fnms,"-f",NFILE,fnm);
772 if (nfile == 0)
774 gmx_fatal(FARGS,"No input files!");
777 if (nd < 0)
779 gmx_fatal(FARGS,"Can not have negative number of digits");
781 prec = pow(10,-nd);
782 sprintf(lamformat,"%%6.3f");
783 sprintf( dgformat,"%%%d.%df",3+nd,nd);
784 /* the format strings of the results in kT */
785 sprintf( ktformat,"%%%d.%df",5+nd,nd);
786 sprintf( sktformat,"%%%ds",6+nd);
787 /* the format strings of the errors in kT */
788 sprintf( kteformat,"%%%d.%df",3+nd,nd);
789 sprintf( skteformat,"%%%ds",4+nd);
790 sprintf(xvg2format,"%s %s\n","%g",dgformat);
791 sprintf(xvg3format,"%s %s %s\n","%g",dgformat,dgformat);
794 snew(ba,nfile);
795 snew(results,nfile-1);
796 snew(partsum,(nbmax+1)*(nbmax+1));
797 n1 = 0;
798 nm = 0;
799 for(f=0; f<nfile; f++)
801 read_barsim(fnms[f],begin,end,temp,&ba[f]);
802 if (f > 0 && ba[f].temp != ba[0].temp)
804 printf("\nWARNING: temperature for file '%s' (%g) is not equal to that of file '%s' (%g)\n\n",fnms[f],ba[f].temp,fnms[0],ba[0].temp);
807 if (ba[f].nset == 0)
809 gmx_fatal(FARGS,"File '%s' contains less than two columns",fnms[f]);
811 else if (ba[f].nset == 1)
813 n1++;
815 else
817 nm++;
820 printf("\n");
822 if (n1 > 0 && nm > 0)
824 gmx_fatal(FARGS,"Some dhdl files contain only one value (assuming dH/dl), while others contain multiple values (assuming dH/dl and Delta H), will not proceed because of possible inconsistencies");
827 /* Sort the data sets on lambda */
828 for(f=0; f<nfile-1; f++)
830 fm = f;
831 for(f2=f+1; f2<nfile; f2++)
833 if (ba[f2].lambda[0] == ba[fm].lambda[0])
835 gmx_fatal(FARGS,"There are multiple files with lambda = %g",
836 ba[fm].lambda[0]);
838 else if (ba[f2].lambda[0] < ba[fm].lambda[0])
840 fm = f2;
843 ba_tmp = ba[f];
844 ba[f] = ba[fm];
845 ba[fm] = ba_tmp;
848 if (n1 > 0)
850 printf("Only one y value in all files,\n"
851 "assuming the Hamiltonian depends linearly on lambda\n\n");
854 fpb = NULL;
855 if (opt2bSet("-o",NFILE,fnm))
857 sprintf(buf,"%s (%s)","\\DeltaG",unit_energy);
858 fpb = xvgropen_type(opt2fn("-o",NFILE,fnm),"Free energy differences",
859 "\\lambda",buf,exvggtXYDY,oenv);
862 fpi = NULL;
863 if (opt2bSet("-oi",NFILE,fnm))
865 sprintf(buf,"%s (%s)","\\DeltaG",unit_energy);
866 fpi = xvgropen(opt2fn("-oi",NFILE,fnm),"Free energy integral",
867 "\\lambda",buf,oenv);
870 /* first calculate results */
871 bEE = TRUE;
872 for(f=0; f<nfile-1; f++)
874 /* Determine the free energy difference with a factor of 10
875 * more accuracy than requested for printing.
877 calc_bar(&ba[f], &ba[f+1], n1>0, 0.1*prec, nbmin, nbmax,
878 &(results[f]), &bEE, partsum);
881 /* print results in kT */
882 kT = BOLTZ*ba[0].temp;
883 beta = 1/kT;
885 printf("\nTemperature: %g K\n", ba[0].temp);
887 printf("\nDetailed results in kT (see help for explanation):\n\n");
888 printf("%6s ", " lam_A");
889 printf("%6s ", " lam_B");
890 printf(sktformat, "DG ");
891 printf(skteformat, "+/- ");
892 printf(sktformat, "s_A ");
893 printf(skteformat, "+/- " );
894 printf(sktformat, "s_B ");
895 printf(skteformat, "+/- " );
896 printf(sktformat, "stdev ");
897 printf(skteformat, "+/- ");
898 printf("\n");
899 for(f=0; f<nfile-1; f++)
901 printf(lamformat, results[f].lambda_a);
902 printf(" ");
903 printf(lamformat, results[f].lambda_b);
904 printf(" ");
905 printf(ktformat, results[f].dg);
906 printf(" ");
907 printf(kteformat, results[f].dg_err);
908 printf(" ");
909 printf(ktformat, results[f].sa);
910 printf(" ");
911 printf(kteformat, results[f].sa_err);
912 printf(" ");
913 printf(ktformat, results[f].sb);
914 printf(" ");
915 printf(kteformat, results[f].sb_err);
916 printf(" ");
917 printf(ktformat, results[f].dg_stddev);
918 printf(" ");
919 printf(kteformat, results[f].dg_stddev_err);
920 printf("\n");
922 /* Check for negative relative entropy with a 95% certainty. */
923 if (results[f].sa < -2*results[f].sa_err ||
924 results[f].sb < -2*results[f].sb_err)
926 result_OK=FALSE;
930 if (!result_OK)
932 printf("\nWARNING: Some of these results violate the Second Law of "
933 "Thermodynamics: \n"
934 " This is can be the result of severe undersampling, or "
935 "(more likely)\n"
936 " there is something wrong with the simulations.\n");
940 /* final results in kJ/mol */
941 printf("\n\nFinal results in kJ/mol:\n\n");
942 dg_tot = 0;
943 for(f=0; f<nfile-1; f++)
946 if (fpi != NULL)
948 fprintf(fpi, xvg2format, ba[f].lambda[0], dg_tot);
952 if (fpb != NULL)
954 fprintf(fpb, xvg3format,
955 0.5*(ba[f].lambda[0] + ba[f+1].lambda[0]),
956 results[f].dg,results[f].dg_err);
959 /*printf("lambda %4.2f - %4.2f, DG ", results[f].lambda_a,
960 results[f].lambda_b);*/
961 printf("lambda ");
962 printf(lamformat, results[f].lambda_a);
963 printf(" - ");
964 printf(lamformat, results[f].lambda_b);
965 printf(", DG ");
967 printf(dgformat,results[f].dg*kT);
968 printf(" +/- ");
969 printf(dgformat,results[f].dg_err*kT);
971 printf("\n");
972 dg_tot += results[f].dg;
974 printf("\n");
975 printf("total ");
976 printf(lamformat, ba[0].lambda[0]);
977 printf(" - ");
978 printf(lamformat, ba[nfile-1].lambda[0]);
979 printf(", DG ");
981 printf(dgformat,dg_tot*kT);
982 if (bEE)
984 printf(" +/- ");
985 printf(dgformat,bar_err(nbmin,nbmax,partsum)*kT);
987 printf("\n");
989 if (fpi != NULL)
991 fprintf(fpi, xvg2format,
992 ba[nfile-1].lambda[0], dg_tot);
993 ffclose(fpi);
996 do_view(oenv,opt2fn_null("-o",NFILE,fnm),"-xydy");
997 do_view(oenv,opt2fn_null("-oi",NFILE,fnm),"-xydy");
999 thanx(stderr);
1001 return 0;