2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief Implementation of the Weighted Histogram Analysis Method (WHAM)
41 * \author Jochen Hub <jhub@gwdg.de>
56 #include "gromacs/commandline/pargs.h"
57 #include "gromacs/fileio/tpxio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxana/gmx_ana.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/units.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/pull-params.h"
66 #include "gromacs/pulling/pull.h"
67 #include "gromacs/random/tabulatednormaldistribution.h"
68 #include "gromacs/random/threefry.h"
69 #include "gromacs/random/uniformintdistribution.h"
70 #include "gromacs/random/uniformrealdistribution.h"
71 #include "gromacs/utility/arraysize.h"
72 #include "gromacs/utility/cstringutil.h"
73 #include "gromacs/utility/exceptions.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/futil.h"
76 #include "gromacs/utility/gmxomp.h"
77 #include "gromacs/utility/pleasecite.h"
78 #include "gromacs/utility/smalloc.h"
80 //! longest file names allowed in input files
81 #define WHAM_MAXFILELEN 2048
84 * enum for energy units
87 enSel, en_kJ, en_kCal, en_kT, enNr
90 * enum for type of input files (pdos, tpr, or pullf)
93 whamin_unknown, whamin_tpr, whamin_pullxf, whamin_pdo
96 /*! \brief enum for bootstrapping method
98 * These bootstrap methods are supported:
99 * - bootstrap complete histograms with continuous weights (Bayesian bootstrap)
100 * (bsMethod_BayesianHist)
101 * - bootstrap complete histograms (bsMethod_hist)
102 * - bootstrap trajectories from given umbrella histograms. This generates new
103 * "synthetic" histograms (bsMethod_traj)
104 * - bootstrap trajectories from Gaussian with mu/sigma computed from
105 * the respective histogram (bsMethod_trajGauss). This gives very similar
106 * results compared to bsMethod_traj.
108 * ********************************************************************
109 * FOR MORE DETAILS ON THE BOOTSTRAP METHODS (INCLUDING EXAMPLES), SEE
110 * JS Hub, BL de Groot, D van der Spoel
111 * g_wham - A free weighted histogram analysis implementation including
112 * robust error and autocorrelation estimates,
113 * J Chem Theory Comput, 6(12), 3713-3720 (2010)
114 * ********************************************************************
117 bsMethod_unknown, bsMethod_BayesianHist, bsMethod_hist,
118 bsMethod_traj, bsMethod_trajGauss
121 //! Parameters of one pull coodinate
124 int pull_type; //!< such as constraint, umbrella, ...
125 int geometry; //!< such as distance, direction, cylinder
126 int ngroup; //!< the number of pull groups involved
127 ivec dim; //!< pull dimension with geometry distance
128 int ndim; //!< nr of pull_dim != 0
129 real k; //!< force constants in tpr file
130 real init_dist; //!< reference displacement
131 char coord_unit[256]; //!< unit of the displacement
134 //! Parameters of the umbrella potentials
138 * \name Using umbrella pull code since gromacs 4.x
141 int npullcrds; //!< nr of umbrella pull coordinates for reading
142 t_pullcoord *pcrd; //!< the pull coordinates
143 gmx_bool bPrintCOM; //!< COMs of pull groups writtn in pullx.xvg
144 gmx_bool bPrintRefValue; //!< Reference value for the coordinate written in pullx.xvg
145 gmx_bool bPrintComp; //!< Components of pull distance written to pullx.xvg ?
149 * \name Using PDO files common until gromacs 3.x
157 char PullName[4][256];
159 double UmbCons[4][3];
163 //! Data in the umbrella histograms
166 int nPull; //!< nr of pull groups in this pdo or pullf/x file
167 double **Histo; //!< nPull histograms
168 double **cum; //!< nPull cumulative distribution functions
169 int nBin; //!< nr of bins. identical to opt->bins
170 double *k; //!< force constants for the nPull coords
171 double *pos; //!< umbrella positions for the nPull coords
172 double *z; //!< z=(-Fi/kT) for the nPull coords. These values are iteratively computed during wham
173 int *N; //!< nr of data points in nPull histograms.
174 int *Ntot; //!< also nr of data points. N and Ntot only differ if bHistEq==TRUE
176 /*! \brief g = 1 + 2*tau[int]/dt where tau is the integrated autocorrelation time.
178 * Compare, e.g. Ferrenberg/Swendsen, PRL 63:1195 (1989),
179 * Kumar et al, J Comp Chem 13, 1011-1021 (1992), eq. 28
182 double *tau; //!< intetrated autocorrelation time (IACT)
183 double *tausmooth; //!< smoothed IACT
185 double dt; //!< timestep in the input data. Can be adapted with gmx wham option -dt
187 /*! \brief TRUE, if any data point of the histogram is within min and max, otherwise FALSE */
189 real **ztime; //!< input data z(t) as a function of time. Required to compute ACTs
191 /*! \brief average force estimated from average displacement, fAv=dzAv*k
193 * Used for integration to guess the potential.
196 real *aver; //!< average of histograms
197 real *sigma; //!< stddev of histograms
198 double *bsWeight; //!< for bootstrapping complete histograms with continuous weights
201 //! Selection of pull coordinates to be used in WHAM (one structure for each tpr file)
204 int n; //!< total nr of pull coords in this tpr file
205 int nUse; //!< nr of pull coords used
206 gmx_bool *bUse; //!< boolean array of size n. =1 if used, =0 if not
209 //! Parameters of WHAM
216 const char *fnTpr, *fnPullf, *fnCoordSel;
217 const char *fnPdo, *fnPullx; //!< file names of input
218 gmx_bool bTpr, bPullf, bPdo, bPullx; //!< input file types given?
219 real tmin, tmax, dt; //!< only read input within tmin and tmax with dt
221 gmx_bool bInitPotByIntegration; //!< before WHAM, guess potential by force integration. Yields 1.5 to 2 times faster convergence
222 int stepUpdateContrib; //!< update contribution table every ... iterations. Accelerates WHAM.
223 int nCoordsel; //!< if >0: use only certain group in WHAM, if ==0: use all groups
224 t_coordselection *coordsel; //!< for each tpr file: which pull coordinates to use in WHAM?
227 * \name Basic WHAM options
230 int bins; //!< nr of bins, min, max, and dz of profile
232 real Temperature, Tolerance; //!< temperature, converged when probability changes less than Tolerance
233 gmx_bool bCycl; //!< generate cyclic (periodic) PMF
236 * \name Output control
239 gmx_bool bLog; //!< energy output (instead of probability) for profile
240 int unit; //!< unit for PMF output kJ/mol or kT or kCal/mol
241 gmx_bool bSym; //!< symmetrize PMF around z=0 after WHAM, useful for membranes etc.
242 /*! \brief after wham, set prof to zero at this z-position.
243 * When bootstrapping, set zProf0 to a "stable" reference position.
246 gmx_bool bProf0Set; //!< setting profile to 0 at zProf0?
248 gmx_bool bBoundsOnly, bHistOnly; //!< determine min and max, or write histograms and exit
249 gmx_bool bAuto; //!< determine min and max automatically but do not exit
251 gmx_bool verbose; //!< more noisy wham mode
252 int stepchange; //!< print maximum change in prof after how many interations
253 gmx_output_env_t *oenv; //!< xvgr options
256 * \name Autocorrelation stuff
259 gmx_bool bTauIntGiven, bCalcTauInt; //!< IACT given or should be calculated?
260 real sigSmoothIact; //!< sigma of Gaussian to smooth ACTs
261 gmx_bool bAllowReduceIact; //!< Allow to reduce ACTs during smoothing. Otherwise ACT are only increased during smoothing
262 real acTrestart; //!< when computing ACT, time between restarting points
264 /* \brief Enforce the same weight for each umbella window, that is
265 * calculate with the same number of data points for
266 * each window. That can be reasonable, if the histograms
267 * have different length, but due to autocorrelation,
268 * a longer simulation should not have larger weightin wham.
274 * \name Bootstrapping stuff
277 int nBootStrap; //!< nr of bootstraps (50 is usually enough)
279 /* \brief bootstrap method
281 * if == bsMethod_hist, consider complete histograms as independent
282 * data points and, hence, only mix complete histograms.
283 * if == bsMethod_BayesianHist, consider complete histograms
284 * as independent data points, but assign random weights
285 * to the histograms during the bootstrapping ("Bayesian bootstrap")
286 * In case of long correlations (e.g., inside a channel), these
287 * will yield a more realistic error.
288 * if == bsMethod_traj(Gauss), generate synthetic histograms
290 * histogram by generating an autocorrelated random sequence
291 * that is distributed according to the respective given
292 * histogram. With bsMethod_trajGauss, bootstrap from a Gaussian
293 * (instead of from the umbrella histogram) to generate a new
298 /* \brief autocorrelation time (ACT) used to generate synthetic histograms. If ==0, use calculated ACF */
301 /* \brief when mixing histograms, mix only histograms withing blocks
302 long the reaction coordinate xi. Avoids gaps along xi. */
303 int histBootStrapBlockLength;
305 int bsSeed; //!< random seed for bootstrapping
307 /* \brief Write cumulative distribution functions (CDFs) of histograms
308 and write the generated histograms for each bootstrap */
312 * \name tabulated umbrella potential stuff
316 double *tabX, *tabY, tabMin, tabMax, tabDz;
319 gmx::DefaultRandomEngine rng; //!< gromacs random number generator
320 gmx::TabulatedNormalDistribution<> normalDistribution; //!< Uses default: real output, 14-bit table
323 //! Make an umbrella window (may contain several histograms)
324 t_UmbrellaWindow * initUmbrellaWindows(int nwin)
326 t_UmbrellaWindow *win;
329 for (i = 0; i < nwin; i++)
331 win[i].Histo = win[i].cum = 0;
332 win[i].k = win[i].pos = win[i].z = 0;
333 win[i].N = win[i].Ntot = 0;
334 win[i].g = win[i].tau = win[i].tausmooth = 0;
338 win[i].aver = win[i].sigma = 0;
344 //! Delete an umbrella window (may contain several histograms)
345 void freeUmbrellaWindows(t_UmbrellaWindow *win, int nwin)
348 for (i = 0; i < nwin; i++)
352 for (j = 0; j < win[i].nPull; j++)
354 sfree(win[i].Histo[j]);
359 for (j = 0; j < win[i].nPull; j++)
361 sfree(win[i].cum[j]);
366 for (j = 0; j < win[i].nPull; j++)
368 sfree(win[i].bContrib[j]);
380 sfree(win[i].tausmooth);
381 sfree(win[i].bContrib);
383 sfree(win[i].forceAv);
386 sfree(win[i].bsWeight);
392 * Read and setup tabulated umbrella potential
394 void setup_tab(const char *fn, t_UmbrellaOptions *opt)
399 printf("Setting up tabulated potential from file %s\n", fn);
400 nl = read_xvg(fn, &y, &ny);
404 gmx_fatal(FARGS, "Found %d columns in %s. Expected 2.\n", ny, fn);
406 opt->tabMin = y[0][0];
407 opt->tabMax = y[0][nl-1];
408 opt->tabDz = (opt->tabMax-opt->tabMin)/(nl-1);
411 gmx_fatal(FARGS, "The tabulated potential in %s must be provided in \n"
412 "ascending z-direction", fn);
414 for (i = 0; i < nl-1; i++)
416 if (std::abs(y[0][i+1]-y[0][i]-opt->tabDz) > opt->tabDz/1e6)
418 gmx_fatal(FARGS, "z-values in %s are not equally spaced.\n", ny, fn);
423 for (i = 0; i < nl; i++)
425 opt->tabX[i] = y[0][i];
426 opt->tabY[i] = y[1][i];
428 printf("Found equally spaced tabulated potential from %g to %g, spacing %g\n",
429 opt->tabMin, opt->tabMax, opt->tabDz);
432 //! Read the header of an PDO file (position, force const, nr of groups)
433 void read_pdo_header(FILE * file, t_UmbrellaHeader * header, t_UmbrellaOptions *opt)
436 char Buffer0[256], Buffer1[256], Buffer2[256], Buffer3[256], Buffer4[256];
438 std::istringstream ist;
441 if (fgets(line, 2048, file) == NULL)
443 gmx_fatal(FARGS, "Error reading header from pdo file\n");
446 ist >> Buffer0 >> Buffer1 >> Buffer2;
447 if (std::strcmp(Buffer1, "UMBRELLA"))
449 gmx_fatal(FARGS, "This does not appear to be a valid pdo file. Found %s, expected %s\n"
450 "(Found in first line: `%s')\n",
451 Buffer1, "UMBRELLA", line);
453 if (std::strcmp(Buffer2, "3.0"))
455 gmx_fatal(FARGS, "This does not appear to be a version 3.0 pdo file");
459 if (fgets(line, 2048, file) == NULL)
461 gmx_fatal(FARGS, "Error reading header from pdo file\n");
464 ist >> Buffer0 >> Buffer1 >> Buffer2 >> header->Dims[0] >> header->Dims[1] >> header->Dims[2];
465 /* printf("%d %d %d\n", header->Dims[0],header->Dims[1],header->Dims[2]); */
467 header->nDim = header->Dims[0] + header->Dims[1] + header->Dims[2];
468 if (header->nDim != 1)
470 gmx_fatal(FARGS, "Currently only supports one dimension");
474 if (fgets(line, 2048, file) == NULL)
476 gmx_fatal(FARGS, "Error reading header from pdo file\n");
479 ist >> Buffer0 >> Buffer1 >> header->nSkip;
482 if (fgets(line, 2048, file) == NULL)
484 gmx_fatal(FARGS, "Error reading header from pdo file\n");
487 ist >> Buffer0 >> Buffer1 >> Buffer2 >> header->Reference;
490 if (fgets(line, 2048, file) == NULL)
492 gmx_fatal(FARGS, "Error reading header from pdo file\n");
495 ist >> Buffer0 >> Buffer1 >> Buffer2 >> Buffer3 >> Buffer4 >> header->nPull;
499 printf("\tFound nPull=%d , nSkip=%d, ref=%s\n", header->nPull, header->nSkip,
503 for (i = 0; i < header->nPull; ++i)
505 if (fgets(line, 2048, file) == NULL)
507 gmx_fatal(FARGS, "Error reading header from pdo file\n");
510 ist >> Buffer0 >> Buffer1 >> Buffer2 >> header->PullName[i];
511 ist >> Buffer0 >> Buffer1 >> header->UmbPos[i][0];
512 ist >> Buffer0 >> Buffer1 >> header->UmbCons[i][0];
516 printf("\tpullgroup %d, pullname = %s, UmbPos = %g, UmbConst = %g\n",
517 i, header->PullName[i], header->UmbPos[i][0], header->UmbCons[i][0]);
521 if (fgets(line, 2048, file) == NULL)
523 gmx_fatal(FARGS, "Cannot read from file\n");
527 if (std::strcmp(Buffer3, "#####") != 0)
529 gmx_fatal(FARGS, "Expected '#####', found %s. Hick.\n", Buffer3);
534 static char *fgets3(FILE *fp, char ptr[], int *len)
539 if (fgets(ptr, *len-1, fp) == NULL)
544 while ((std::strchr(ptr, '\n') == NULL) && (!feof(fp)))
546 /* This line is longer than len characters, let's increase len! */
550 if (fgets(p-1, STRLEN, fp) == NULL)
555 slen = std::strlen(ptr);
556 if (ptr[slen-1] == '\n')
564 /*! \brief Read the data columns of and PDO file.
566 * TO DO: Get rid of the scanf function to avoid the clang warning.
567 * At the moment, this warning is avoided by hiding the format string
568 * the variable fmtlf.
570 void read_pdo_data(FILE * file, t_UmbrellaHeader * header,
571 int fileno, t_UmbrellaWindow * win,
572 t_UmbrellaOptions *opt,
573 gmx_bool bGetMinMax, real *mintmp, real *maxtmp)
575 int i, inttemp, bins, count, ntot;
576 real minval, maxval, minfound = 1e20, maxfound = -1e20;
577 double temp, time, time0 = 0, dt;
579 t_UmbrellaWindow * window = 0;
580 gmx_bool timeok, dt_ok = 1;
581 char *tmpbuf = 0, fmt[256], fmtign[256], fmtlf[5] = "%lf";
582 int len = STRLEN, dstep = 1;
583 const int blocklen = 4096;
593 /* Need to alocate memory and set up structure */
594 window->nPull = header->nPull;
597 snew(window->Histo, window->nPull);
598 snew(window->z, window->nPull);
599 snew(window->k, window->nPull);
600 snew(window->pos, window->nPull);
601 snew(window->N, window->nPull);
602 snew(window->Ntot, window->nPull);
603 snew(window->g, window->nPull);
604 snew(window->bsWeight, window->nPull);
606 window->bContrib = 0;
608 if (opt->bCalcTauInt)
610 snew(window->ztime, window->nPull);
616 snew(lennow, window->nPull);
618 for (i = 0; i < window->nPull; ++i)
621 window->bsWeight[i] = 1.;
622 snew(window->Histo[i], bins);
623 window->k[i] = header->UmbCons[i][0];
624 window->pos[i] = header->UmbPos[i][0];
628 if (opt->bCalcTauInt)
630 window->ztime[i] = 0;
634 /* Done with setup */
640 minval = maxval = bins = 0; /* Get rid of warnings */
645 while ( (ptr = fgets3(file, tmpbuf, &len)) != NULL)
649 if (ptr[0] == '#' || std::strlen(ptr) < 2)
654 /* Initiate format string */
656 std::strcat(fmtign, "%*s");
658 sscanf(ptr, fmtlf, &time); /* printf("Time %f\n",time); */
659 /* Round time to fs */
660 time = 1.0/1000*( static_cast<gmx_int64_t> (time*1000+0.5) );
662 /* get time step of pdo file */
672 dstep = static_cast<int>(opt->dt/dt+0.5);
680 window->dt = dt*dstep;
685 dt_ok = ((count-1)%dstep == 0);
686 timeok = (dt_ok && time >= opt->tmin && time <= opt->tmax);
688 printf(" time = %f, (tmin,tmax)=(%e,%e), dt_ok=%d timeok=%d\n",
689 time,opt->tmin, opt->tmax, dt_ok,timeok); */
693 for (i = 0; i < header->nPull; ++i)
695 std::strcpy(fmt, fmtign);
696 std::strcat(fmt, "%lf"); /* Creating a format stings such as "%*s...%*s%lf" */
697 std::strcat(fmtign, "%*s"); /* ignoring one more entry in the next loop */
698 if (sscanf(ptr, fmt, &temp))
700 temp += header->UmbPos[i][0];
714 if (opt->bCalcTauInt)
716 /* save time series for autocorrelation analysis */
717 ntot = window->Ntot[i];
718 if (ntot >= lennow[i])
720 lennow[i] += blocklen;
721 srenew(window->ztime[i], lennow[i]);
723 window->ztime[i][ntot] = temp;
727 temp /= (maxval-minval);
729 temp = std::floor(temp);
731 inttemp = static_cast<int> (temp);
738 else if (inttemp >= bins)
744 if (inttemp >= 0 && inttemp < bins)
746 window->Histo[i][inttemp] += 1.;
754 if (time > opt->tmax)
758 printf("time %f larger than tmax %f, stop reading pdo file\n", time, opt->tmax);
774 /*! \brief Set identical weights for all histograms
776 * Normally, the weight is given by the number data points in each
777 * histogram, together with the autocorrelation time. This can be overwritten
778 * by this routine (not recommended). Since we now support autocorrelations, it is better to set
779 * an appropriate autocorrelation times instead of using this function.
781 void enforceEqualWeights(t_UmbrellaWindow * window, int nWindows)
783 int i, k, j, NEnforced;
786 NEnforced = window[0].Ntot[0];
787 printf("\nFound -hist-eq. Enforcing equal weights for all histograms, \ni.e. doing a "
788 "non-weighted histogram analysis method. Ndata = %d\n", NEnforced);
789 /* enforce all histograms to have the same weight as the very first histogram */
791 for (j = 0; j < nWindows; ++j)
793 for (k = 0; k < window[j].nPull; ++k)
795 ratio = 1.0*NEnforced/window[j].Ntot[k];
796 for (i = 0; i < window[0].nBin; ++i)
798 window[j].Histo[k][i] *= ratio;
800 window[j].N[k] = static_cast<int>(ratio*window[j].N[k] + 0.5);
805 /*! \brief Simple linear interpolation between two given tabulated points
807 double tabulated_pot(double dist, t_UmbrellaOptions *opt)
810 double pl, pu, dz, dp;
812 jl = static_cast<int> (std::floor((dist-opt->tabMin)/opt->tabDz));
814 if (jl < 0 || ju >= opt->tabNbins)
816 gmx_fatal(FARGS, "Distance %f out of bounds of tabulated potential (jl=%d, ju=%d).\n"
817 "Provide an extended table.", dist, jl, ju);
821 dz = dist-opt->tabX[jl];
822 dp = (pu-pl)*dz/opt->tabDz;
828 * Check which bins substiantially contribute (accelerates WHAM)
830 * Don't worry, that routine does not mean we compute the PMF in limited precision.
831 * After rapid convergence (using only substiantal contributions), we always switch to
834 void setup_acc_wham(double *profile, t_UmbrellaWindow * window, int nWindows,
835 t_UmbrellaOptions *opt)
837 int i, j, k, nGrptot = 0, nContrib = 0, nTot = 0;
838 double U, min = opt->min, dz = opt->dz, temp, ztot_half, distance, ztot, contrib1, contrib2;
839 gmx_bool bAnyContrib;
840 static int bFirst = 1;
841 static double wham_contrib_lim;
845 for (i = 0; i < nWindows; ++i)
847 nGrptot += window[i].nPull;
849 wham_contrib_lim = opt->Tolerance/nGrptot;
852 ztot = opt->max-opt->min;
855 for (i = 0; i < nWindows; ++i)
857 if (!window[i].bContrib)
859 snew(window[i].bContrib, window[i].nPull);
861 for (j = 0; j < window[i].nPull; ++j)
863 if (!window[i].bContrib[j])
865 snew(window[i].bContrib[j], opt->bins);
868 for (k = 0; k < opt->bins; ++k)
870 temp = (1.0*k+0.5)*dz+min;
871 distance = temp - window[i].pos[j]; /* distance to umbrella center */
873 { /* in cyclic wham: */
874 if (distance > ztot_half) /* |distance| < ztot_half */
878 else if (distance < -ztot_half)
883 /* Note: there are two contributions to bin k in the wham equations:
884 i) N[j]*exp(- U/(BOLTZ*opt->Temperature) + window[i].z[j])
885 ii) exp(- U/(BOLTZ*opt->Temperature))
886 where U is the umbrella potential
887 If any of these number is larger wham_contrib_lim, I set contrib=TRUE
892 U = 0.5*window[i].k[j]*gmx::square(distance); /* harmonic potential assumed. */
896 U = tabulated_pot(distance, opt); /* Use tabulated potential */
898 contrib1 = profile[k]*std::exp(-U/(BOLTZ*opt->Temperature));
899 contrib2 = window[i].N[j]*std::exp(-U/(BOLTZ*opt->Temperature) + window[i].z[j]);
900 window[i].bContrib[j][k] = (contrib1 > wham_contrib_lim || contrib2 > wham_contrib_lim);
901 bAnyContrib = (bAnyContrib | window[i].bContrib[j][k]);
902 if (window[i].bContrib[j][k])
908 /* If this histo is far outside min and max all bContrib may be FALSE,
909 causing a floating point exception later on. To avoid that, switch
913 for (k = 0; k < opt->bins; ++k)
915 window[i].bContrib[j][k] = TRUE;
922 printf("Initialized rapid wham stuff (contrib tolerance %g)\n"
923 "Evaluating only %d of %d expressions.\n\n", wham_contrib_lim, nContrib, nTot);
928 printf("Updated rapid wham stuff. (evaluating only %d of %d contributions)\n",
934 //! Compute the PMF (one of the two main WHAM routines)
935 void calc_profile(double *profile, t_UmbrellaWindow * window, int nWindows,
936 t_UmbrellaOptions *opt, gmx_bool bExact)
938 double ztot_half, ztot, min = opt->min, dz = opt->dz;
940 ztot = opt->max-opt->min;
947 int nthreads = gmx_omp_get_max_threads();
948 int thread_id = gmx_omp_get_thread_num();
950 int i0 = thread_id*opt->bins/nthreads;
951 int i1 = std::min(opt->bins, ((thread_id+1)*opt->bins)/nthreads);
953 for (i = i0; i < i1; ++i)
956 double num, denom, invg, temp = 0, distance, U = 0;
958 for (j = 0; j < nWindows; ++j)
960 for (k = 0; k < window[j].nPull; ++k)
962 invg = 1.0/window[j].g[k] * window[j].bsWeight[k];
963 temp = (1.0*i+0.5)*dz+min;
964 num += invg*window[j].Histo[k][i];
966 if (!(bExact || window[j].bContrib[k][i]))
970 distance = temp - window[j].pos[k]; /* distance to umbrella center */
972 { /* in cyclic wham: */
973 if (distance > ztot_half) /* |distance| < ztot_half */
977 else if (distance < -ztot_half)
985 U = 0.5*window[j].k[k]*gmx::square(distance); /* harmonic potential assumed. */
989 U = tabulated_pot(distance, opt); /* Use tabulated potential */
991 denom += invg*window[j].N[k]*std::exp(-U/(BOLTZ*opt->Temperature) + window[j].z[k]);
994 profile[i] = num/denom;
997 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1001 //! Compute the free energy offsets z (one of the two main WHAM routines)
1002 double calc_z(double * profile, t_UmbrellaWindow * window, int nWindows,
1003 t_UmbrellaOptions *opt, gmx_bool bExact)
1005 double min = opt->min, dz = opt->dz, ztot_half, ztot;
1006 double maxglob = -1e20;
1008 ztot = opt->max-opt->min;
1011 #pragma omp parallel
1015 int nthreads = gmx_omp_get_max_threads();
1016 int thread_id = gmx_omp_get_thread_num();
1018 int i0 = thread_id*nWindows/nthreads;
1019 int i1 = std::min(nWindows, ((thread_id+1)*nWindows)/nthreads);
1020 double maxloc = -1e20;
1022 for (i = i0; i < i1; ++i)
1024 double total = 0, temp, distance, U = 0;
1027 for (j = 0; j < window[i].nPull; ++j)
1030 for (k = 0; k < window[i].nBin; ++k)
1032 if (!(bExact || window[i].bContrib[j][k]))
1036 temp = (1.0*k+0.5)*dz+min;
1037 distance = temp - window[i].pos[j]; /* distance to umbrella center */
1039 { /* in cyclic wham: */
1040 if (distance > ztot_half) /* |distance| < ztot_half */
1044 else if (distance < -ztot_half)
1052 U = 0.5*window[i].k[j]*gmx::square(distance); /* harmonic potential assumed. */
1056 U = tabulated_pot(distance, opt); /* Use tabulated potential */
1058 total += profile[k]*std::exp(-U/(BOLTZ*opt->Temperature));
1060 /* Avoid floating point exception if window is far outside min and max */
1063 total = -std::log(total);
1069 temp = std::abs(total - window[i].z[j]);
1074 window[i].z[j] = total;
1077 /* Now get maximum maxloc from the threads and put in maxglob */
1078 if (maxloc > maxglob)
1080 #pragma omp critical
1082 if (maxloc > maxglob)
1089 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1095 //! Make PMF symmetric around 0 (useful e.g. for membranes)
1096 void symmetrizeProfile(double* profile, t_UmbrellaOptions *opt)
1098 int i, j, bins = opt->bins;
1099 double *prof2, min = opt->min, max = opt->max, dz = opt->dz, zsym, deltaz, profsym;
1102 if (min > 0. || max < 0.)
1104 gmx_fatal(FARGS, "Cannot symmetrize profile around z=0 with min=%f and max=%f\n",
1105 opt->min, opt->max);
1110 for (i = 0; i < bins; i++)
1114 /* bin left of zsym */
1115 j = static_cast<int> (std::floor((zsym-min)/dz-0.5));
1116 if (j >= 0 && (j+1) < bins)
1118 /* interpolate profile linearly between bins j and j+1 */
1119 z1 = min+(j+0.5)*dz;
1121 profsym = profile[j] + (profile[j+1]-profile[j])/dz*deltaz;
1122 /* average between left and right */
1123 prof2[i] = 0.5*(profsym+profile[i]);
1127 prof2[i] = profile[i];
1131 std::memcpy(profile, prof2, bins*sizeof(double));
1135 //! Set energy unit (kJ/mol,kT,kCal/mol) and set it to zero at opt->zProf0
1136 void prof_normalization_and_unit(double * profile, t_UmbrellaOptions *opt)
1139 double unit_factor = 1., diff;
1143 /* Not log? Nothing to do! */
1149 /* Get profile in units of kT, kJ/mol, or kCal/mol */
1150 if (opt->unit == en_kT)
1154 else if (opt->unit == en_kJ)
1156 unit_factor = BOLTZ*opt->Temperature;
1158 else if (opt->unit == en_kCal)
1160 unit_factor = (BOLTZ/CAL2JOULE)*opt->Temperature;
1164 gmx_fatal(FARGS, "Sorry, I don't know this energy unit.");
1167 for (i = 0; i < bins; i++)
1169 if (profile[i] > 0.0)
1171 profile[i] = -std::log(profile[i])*unit_factor;
1175 /* shift to zero at z=opt->zProf0 */
1176 if (!opt->bProf0Set)
1182 /* Get bin with shortest distance to opt->zProf0
1183 (-0.5 from bin position and +0.5 from rounding cancel) */
1184 imin = static_cast<int>((opt->zProf0-opt->min)/opt->dz);
1189 else if (imin >= bins)
1193 diff = profile[imin];
1197 for (i = 0; i < bins; i++)
1203 //! Make an array of random integers (used for bootstrapping)
1204 void getRandomIntArray(int nPull, int blockLength, int* randomArray, gmx::DefaultRandomEngine * rng)
1206 gmx::UniformIntDistribution<int> dist(0, blockLength-1);
1208 int ipull, blockBase, nr, ipullRandom;
1210 if (blockLength == 0)
1212 blockLength = nPull;
1215 for (ipull = 0; ipull < nPull; ipull++)
1217 blockBase = (ipull/blockLength)*blockLength;
1219 { /* make sure nothing bad happens in the last block */
1220 nr = dist(*rng); // [0,blockLength-1]
1221 ipullRandom = blockBase + nr;
1223 while (ipullRandom >= nPull);
1224 if (ipullRandom < 0 || ipullRandom >= nPull)
1226 gmx_fatal(FARGS, "Ups, random iWin = %d, nPull = %d, nr = %d, "
1227 "blockLength = %d, blockBase = %d\n",
1228 ipullRandom, nPull, nr, blockLength, blockBase);
1230 randomArray[ipull] = ipullRandom;
1232 /*for (ipull=0; ipull<nPull; ipull++)
1233 printf("%d ",randomArray[ipull]); printf("\n"); */
1236 /*! \brief Set pull group information of a synthetic histogram
1238 * This is used when bootstapping new trajectories and thereby create new histogtrams,
1239 * but it is not required if we bootstrap complete histograms.
1241 void copy_pullgrp_to_synthwindow(t_UmbrellaWindow *synthWindow,
1242 t_UmbrellaWindow *thisWindow, int pullid)
1244 synthWindow->N [0] = thisWindow->N [pullid];
1245 synthWindow->Histo [0] = thisWindow->Histo [pullid];
1246 synthWindow->pos [0] = thisWindow->pos [pullid];
1247 synthWindow->z [0] = thisWindow->z [pullid];
1248 synthWindow->k [0] = thisWindow->k [pullid];
1249 synthWindow->bContrib[0] = thisWindow->bContrib [pullid];
1250 synthWindow->g [0] = thisWindow->g [pullid];
1251 synthWindow->bsWeight[0] = thisWindow->bsWeight [pullid];
1254 /*! \brief Calculate cumulative distribution function of of all histograms.
1256 * This allow to create random number sequences
1257 * which are distributed according to the histograms. Required to generate
1258 * the "synthetic" histograms for the Bootstrap method
1260 void calc_cumulatives(t_UmbrellaWindow *window, int nWindows,
1261 t_UmbrellaOptions *opt, const char *fnhist, const char *xlabel)
1265 char *fn = 0, *buf = 0;
1268 if (opt->bs_verbose)
1270 snew(fn, std::strlen(fnhist)+10);
1271 snew(buf, std::strlen(fnhist)+10);
1272 sprintf(fn, "%s_cumul.xvg", std::strncpy(buf, fnhist, std::strlen(fnhist)-4));
1273 fp = xvgropen(fn, "CDFs of umbrella windows", xlabel, "CDF", opt->oenv);
1277 for (i = 0; i < nWindows; i++)
1279 snew(window[i].cum, window[i].nPull);
1280 for (j = 0; j < window[i].nPull; j++)
1282 snew(window[i].cum[j], nbin+1);
1283 window[i].cum[j][0] = 0.;
1284 for (k = 1; k <= nbin; k++)
1286 window[i].cum[j][k] = window[i].cum[j][k-1]+window[i].Histo[j][k-1];
1289 /* normalize CDFs. Ensure cum[nbin]==1 */
1290 last = window[i].cum[j][nbin];
1291 for (k = 0; k <= nbin; k++)
1293 window[i].cum[j][k] /= last;
1298 printf("Cumulative distribution functions of all histograms created.\n");
1299 if (opt->bs_verbose)
1301 for (k = 0; k <= nbin; k++)
1303 fprintf(fp, "%g\t", opt->min+k*opt->dz);
1304 for (i = 0; i < nWindows; i++)
1306 for (j = 0; j < window[i].nPull; j++)
1308 fprintf(fp, "%g\t", window[i].cum[j][k]);
1313 printf("Wrote cumulative distribution functions to %s\n", fn);
1321 /*! \brief Return j such that xx[j] <= x < xx[j+1]
1323 * This is used to generate a random sequence distributed according to a histogram
1325 void searchCumulative(double xx[], int n, double x, int *j)
1347 else if (x == xx[n-1])
1357 //! Bootstrap new trajectories and thereby generate new (bootstrapped) histograms
1358 void create_synthetic_histo(t_UmbrellaWindow *synthWindow, t_UmbrellaWindow *thisWindow,
1359 int pullid, t_UmbrellaOptions *opt)
1361 int N, i, nbins, r_index, ibin;
1362 double r, tausteps = 0.0, a, ap, dt, x, invsqrt2, g, y, sig = 0., z, mu = 0.;
1365 N = thisWindow->N[pullid];
1366 dt = thisWindow->dt;
1367 nbins = thisWindow->nBin;
1369 /* tau = autocorrelation time */
1370 if (opt->tauBootStrap > 0.0)
1372 tausteps = opt->tauBootStrap/dt;
1374 else if (opt->bTauIntGiven || opt->bCalcTauInt)
1376 /* calc tausteps from g=1+2tausteps */
1377 g = thisWindow->g[pullid];
1383 "When generating hypothetical trajectories from given umbrella histograms,\n"
1384 "autocorrelation times (ACTs) are required. Otherwise the statistical error\n"
1385 "cannot be predicted. You have 3 options:\n"
1386 "1) Make gmx wham estimate the ACTs (options -ac and -acsig).\n"
1387 "2) Calculate the ACTs by yourself (e.g. with g_analyze) and provide them\n");
1389 " with option -iiact for all umbrella windows.\n"
1390 "3) If all ACTs are identical and know, you can define them with -bs-tau.\n"
1391 " Use option (3) only if you are sure what you're doing, you may severely\n"
1392 " underestimate the error if a too small ACT is given.\n");
1393 gmx_fatal(FARGS, errstr);
1396 synthWindow->N [0] = N;
1397 synthWindow->pos [0] = thisWindow->pos[pullid];
1398 synthWindow->z [0] = thisWindow->z[pullid];
1399 synthWindow->k [0] = thisWindow->k[pullid];
1400 synthWindow->bContrib[0] = thisWindow->bContrib[pullid];
1401 synthWindow->g [0] = thisWindow->g [pullid];
1402 synthWindow->bsWeight[0] = thisWindow->bsWeight[pullid];
1404 for (i = 0; i < nbins; i++)
1406 synthWindow->Histo[0][i] = 0.;
1409 if (opt->bsMethod == bsMethod_trajGauss)
1411 sig = thisWindow->sigma [pullid];
1412 mu = thisWindow->aver [pullid];
1415 /* Genrate autocorrelated Gaussian random variable with autocorrelation time tau
1417 If x and y are random numbers from N(0,1) (Gaussian with average 0 and sigma=1),
1419 z = a*x + sqrt(1-a^2)*y
1420 is also from N(0,1), and cov(z,x) = a. Thus, by gerenating a sequence
1421 x' = a*x + sqrt(1-a^2)*y, the sequnce x(t) is from N(0,1) and has an autocorrelation
1423 C(t) = exp(-t/tau) with tau=-1/ln(a)
1425 Then, use error function to turn the Gaussian random variable into a uniformly
1426 distributed one in [0,1]. Eventually, use cumulative distribution function of
1427 histogram to get random variables distributed according to histogram.
1428 Note: The ACT of the flat distribution and of the generated histogram is not
1429 100% exactly tau, but near tau (my test was 3.8 instead of 4).
1431 a = std::exp(-1.0/tausteps);
1432 ap = std::sqrt(1.0-a*a);
1433 invsqrt2 = 1.0/std::sqrt(2.0);
1435 /* init random sequence */
1436 x = opt->normalDistribution(opt->rng);
1438 if (opt->bsMethod == bsMethod_traj)
1440 /* bootstrap points from the umbrella histograms */
1441 for (i = 0; i < N; i++)
1443 y = opt->normalDistribution(opt->rng);
1445 /* get flat distribution in [0,1] using cumulative distribution function of Gauusian
1446 Note: CDF(Gaussian) = 0.5*{1+erf[x/sqrt(2)]}
1448 r = 0.5*(1+std::erf(x*invsqrt2));
1449 searchCumulative(thisWindow->cum[pullid], nbins+1, r, &r_index);
1450 synthWindow->Histo[0][r_index] += 1.;
1453 else if (opt->bsMethod == bsMethod_trajGauss)
1455 /* bootstrap points from a Gaussian with the same average and sigma
1456 as the respective umbrella histogram. The idea was, that -given
1457 limited sampling- the bootstrapped histograms are otherwise biased
1458 from the limited sampling of the US histos. However, bootstrapping from
1459 the Gaussian seems to yield a similar estimate. */
1463 y = opt->normalDistribution(opt->rng);
1466 ibin = static_cast<int> (std::floor((z-opt->min)/opt->dz));
1471 while ( (ibin += nbins) < 0)
1476 else if (ibin >= nbins)
1478 while ( (ibin -= nbins) >= nbins)
1485 if (ibin >= 0 && ibin < nbins)
1487 synthWindow->Histo[0][ibin] += 1.;
1494 gmx_fatal(FARGS, "Unknown bsMethod (id %d). That should not happen.\n", opt->bsMethod);
1498 /*! \brief Write all histograms to a file
1500 * If bs_index>=0, a number is added to the output file name to allow the ouput of all
1501 * sets of bootstrapped histograms.
1503 void print_histograms(const char *fnhist, t_UmbrellaWindow * window, int nWindows,
1504 int bs_index, t_UmbrellaOptions *opt, const char *xlabel)
1506 char *fn = 0, *buf = 0, title[256];
1512 snew(fn, std::strlen(fnhist)+10);
1513 snew(buf, std::strlen(fnhist)+1);
1514 sprintf(fn, "%s_bs%d.xvg", std::strncpy(buf, fnhist, std::strlen(fnhist)-4), bs_index);
1515 sprintf(title, "Umbrella histograms. Bootstrap #%d", bs_index);
1519 fn = gmx_strdup(fnhist);
1520 std::strcpy(title, "Umbrella histograms");
1523 fp = xvgropen(fn, title, xlabel, "count", opt->oenv);
1526 /* Write histograms */
1527 for (l = 0; l < bins; ++l)
1529 fprintf(fp, "%e\t", (l+0.5)*opt->dz+opt->min);
1530 for (i = 0; i < nWindows; ++i)
1532 for (j = 0; j < window[i].nPull; ++j)
1534 fprintf(fp, "%e\t", window[i].Histo[j][l]);
1541 printf("Wrote %s\n", fn);
1549 //! Used for qsort to sort random numbers
1550 int func_wham_is_larger(const void *a, const void *b)
1569 //! Make random weights for histograms for the Bayesian bootstrap of complete histograms)
1570 void setRandomBsWeights(t_UmbrellaWindow *synthwin, int nAllPull, t_UmbrellaOptions *opt)
1574 gmx::UniformRealDistribution<real> dist(0, nAllPull);
1578 /* generate ordered random numbers between 0 and nAllPull */
1579 for (i = 0; i < nAllPull-1; i++)
1581 r[i] = dist(opt->rng);
1583 qsort((void *)r, nAllPull-1, sizeof(double), &func_wham_is_larger);
1584 r[nAllPull-1] = 1.0*nAllPull;
1586 synthwin[0].bsWeight[0] = r[0];
1587 for (i = 1; i < nAllPull; i++)
1589 synthwin[i].bsWeight[0] = r[i]-r[i-1];
1592 /* avoid to have zero weight by adding a tiny value */
1593 for (i = 0; i < nAllPull; i++)
1595 if (synthwin[i].bsWeight[0] < 1e-5)
1597 synthwin[i].bsWeight[0] = 1e-5;
1604 //! The main bootstrapping routine
1605 void do_bootstrapping(const char *fnres, const char* fnprof, const char *fnhist,
1606 const char *xlabel, char* ylabel, double *profile,
1607 t_UmbrellaWindow * window, int nWindows, t_UmbrellaOptions *opt)
1609 t_UmbrellaWindow * synthWindow;
1610 double *bsProfile, *bsProfiles_av, *bsProfiles_av2, maxchange = 1e20, tmp, stddev;
1611 int i, j, *randomArray = 0, winid, pullid, ib;
1612 int iAllPull, nAllPull, *allPull_winId, *allPull_pullId;
1614 gmx_bool bExact = FALSE;
1616 /* init random generator */
1617 if (opt->bsSeed == 0)
1619 opt->bsSeed = static_cast<int>(gmx::makeRandomSeed());
1621 opt->rng.seed(opt->bsSeed);
1623 snew(bsProfile, opt->bins);
1624 snew(bsProfiles_av, opt->bins);
1625 snew(bsProfiles_av2, opt->bins);
1627 /* Create array of all pull groups. Note that different windows
1628 may have different nr of pull groups
1629 First: Get total nr of pull groups */
1631 for (i = 0; i < nWindows; i++)
1633 nAllPull += window[i].nPull;
1635 snew(allPull_winId, nAllPull);
1636 snew(allPull_pullId, nAllPull);
1638 /* Setup one array of all pull groups */
1639 for (i = 0; i < nWindows; i++)
1641 for (j = 0; j < window[i].nPull; j++)
1643 allPull_winId[iAllPull] = i;
1644 allPull_pullId[iAllPull] = j;
1649 /* setup stuff for synthetic windows */
1650 snew(synthWindow, nAllPull);
1651 for (i = 0; i < nAllPull; i++)
1653 synthWindow[i].nPull = 1;
1654 synthWindow[i].nBin = opt->bins;
1655 snew(synthWindow[i].Histo, 1);
1656 if (opt->bsMethod == bsMethod_traj || opt->bsMethod == bsMethod_trajGauss)
1658 snew(synthWindow[i].Histo[0], opt->bins);
1660 snew(synthWindow[i].N, 1);
1661 snew(synthWindow[i].pos, 1);
1662 snew(synthWindow[i].z, 1);
1663 snew(synthWindow[i].k, 1);
1664 snew(synthWindow[i].bContrib, 1);
1665 snew(synthWindow[i].g, 1);
1666 snew(synthWindow[i].bsWeight, 1);
1669 switch (opt->bsMethod)
1672 snew(randomArray, nAllPull);
1673 printf("\n\nWhen computing statistical errors by bootstrapping entire histograms:\n");
1674 please_cite(stdout, "Hub2006");
1676 case bsMethod_BayesianHist:
1677 /* just copy all histogams into synthWindow array */
1678 for (i = 0; i < nAllPull; i++)
1680 winid = allPull_winId [i];
1681 pullid = allPull_pullId[i];
1682 copy_pullgrp_to_synthwindow(synthWindow+i, window+winid, pullid);
1686 case bsMethod_trajGauss:
1687 calc_cumulatives(window, nWindows, opt, fnhist, xlabel);
1690 gmx_fatal(FARGS, "Unknown bootstrap method. That should not have happened.\n");
1693 /* do bootstrapping */
1694 fp = xvgropen(fnprof, "Bootstrap profiles", xlabel, ylabel, opt->oenv);
1695 for (ib = 0; ib < opt->nBootStrap; ib++)
1697 printf(" *******************************************\n"
1698 " ******** Start bootstrap nr %d ************\n"
1699 " *******************************************\n", ib+1);
1701 switch (opt->bsMethod)
1704 /* bootstrap complete histograms from given histograms */
1705 getRandomIntArray(nAllPull, opt->histBootStrapBlockLength, randomArray, &opt->rng);
1706 for (i = 0; i < nAllPull; i++)
1708 winid = allPull_winId [randomArray[i]];
1709 pullid = allPull_pullId[randomArray[i]];
1710 copy_pullgrp_to_synthwindow(synthWindow+i, window+winid, pullid);
1713 case bsMethod_BayesianHist:
1714 /* keep histos, but assign random weights ("Bayesian bootstrap") */
1715 setRandomBsWeights(synthWindow, nAllPull, opt);
1718 case bsMethod_trajGauss:
1719 /* create new histos from given histos, that is generate new hypothetical
1721 for (i = 0; i < nAllPull; i++)
1723 winid = allPull_winId[i];
1724 pullid = allPull_pullId[i];
1725 create_synthetic_histo(synthWindow+i, window+winid, pullid, opt);
1730 /* write histos in case of verbose output */
1731 if (opt->bs_verbose)
1733 print_histograms(fnhist, synthWindow, nAllPull, ib, opt, xlabel);
1740 std::memcpy(bsProfile, profile, opt->bins*sizeof(double)); /* use profile as guess */
1743 if ( (i%opt->stepUpdateContrib) == 0)
1745 setup_acc_wham(bsProfile, synthWindow, nAllPull, opt);
1747 if (maxchange < opt->Tolerance)
1751 if (((i%opt->stepchange) == 0 || i == 1) && i != 0)
1753 printf("\t%4d) Maximum change %e\n", i, maxchange);
1755 calc_profile(bsProfile, synthWindow, nAllPull, opt, bExact);
1758 while ( (maxchange = calc_z(bsProfile, synthWindow, nAllPull, opt, bExact)) > opt->Tolerance || !bExact);
1759 printf("\tConverged in %d iterations. Final maximum change %g\n", i, maxchange);
1763 prof_normalization_and_unit(bsProfile, opt);
1766 /* symmetrize profile around z=0 */
1769 symmetrizeProfile(bsProfile, opt);
1772 /* save stuff to get average and stddev */
1773 for (i = 0; i < opt->bins; i++)
1776 bsProfiles_av[i] += tmp;
1777 bsProfiles_av2[i] += tmp*tmp;
1778 fprintf(fp, "%e\t%e\n", (i+0.5)*opt->dz+opt->min, tmp);
1780 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
1784 /* write average and stddev */
1785 fp = xvgropen(fnres, "Average and stddev from bootstrapping", xlabel, ylabel, opt->oenv);
1786 if (output_env_get_print_xvgr_codes(opt->oenv))
1788 fprintf(fp, "@TYPE xydy\n");
1790 for (i = 0; i < opt->bins; i++)
1792 bsProfiles_av [i] /= opt->nBootStrap;
1793 bsProfiles_av2[i] /= opt->nBootStrap;
1794 tmp = bsProfiles_av2[i]-gmx::square(bsProfiles_av[i]);
1795 stddev = (tmp >= 0.) ? std::sqrt(tmp) : 0.; /* Catch rouding errors */
1796 fprintf(fp, "%e\t%e\t%e\n", (i+0.5)*opt->dz+opt->min, bsProfiles_av [i], stddev);
1799 printf("Wrote boot strap result to %s\n", fnres);
1802 //! Return type of input file based on file extension (xvg, pdo, or tpr)
1803 int whaminFileType(char *fn)
1806 len = std::strlen(fn);
1807 if (std::strcmp(fn+len-3, "tpr") == 0)
1811 else if (std::strcmp(fn+len-3, "xvg") == 0 || std::strcmp(fn+len-6, "xvg.gz") == 0)
1813 return whamin_pullxf;
1815 else if (std::strcmp(fn+len-3, "pdo") == 0 || std::strcmp(fn+len-6, "pdo.gz") == 0)
1821 gmx_fatal(FARGS, "Unknown file type of %s. Should be tpr, xvg, or pdo.\n", fn);
1823 return whamin_unknown;
1826 //! Read the files names in pdo-files.dat, pullf/x-files.dat, tpr-files.dat
1827 void read_wham_in(const char *fn, char ***filenamesRet, int *nfilesRet,
1828 t_UmbrellaOptions *opt)
1830 char **filename = 0, tmp[WHAM_MAXFILELEN+2];
1831 int nread, sizenow, i, block = 1;
1834 fp = gmx_ffopen(fn, "r");
1837 while (fgets(tmp, sizeof(tmp), fp) != NULL)
1839 if (std::strlen(tmp) >= WHAM_MAXFILELEN)
1841 gmx_fatal(FARGS, "Filename too long in %s. Only %d characters allowed.\n", fn, WHAM_MAXFILELEN);
1843 if (nread >= sizenow)
1846 srenew(filename, sizenow);
1847 for (i = sizenow-block; i < sizenow; i++)
1849 snew(filename[i], WHAM_MAXFILELEN);
1852 /* remove newline if there is one */
1853 if (tmp[std::strlen(tmp)-1] == '\n')
1855 tmp[std::strlen(tmp)-1] = '\0';
1857 std::strcpy(filename[nread], tmp);
1860 printf("Found file %s in %s\n", filename[nread], fn);
1864 *filenamesRet = filename;
1868 //! Open a file or a pipe to a gzipped file
1869 FILE *open_pdo_pipe(const char *fn, t_UmbrellaOptions *opt, gmx_bool *bPipeOpen)
1871 char Buffer[1024], gunzip[1024], *Path = 0;
1873 static gmx_bool bFirst = 1;
1875 /* gzipped pdo file? */
1876 if ((std::strcmp(fn+std::strlen(fn)-3, ".gz") == 0))
1878 /* search gunzip executable */
1879 if (!(Path = getenv("GMX_PATH_GZIP")))
1881 if (gmx_fexist("/bin/gunzip"))
1883 sprintf(gunzip, "%s", "/bin/gunzip");
1885 else if (gmx_fexist("/usr/bin/gunzip"))
1887 sprintf(gunzip, "%s", "/usr/bin/gunzip");
1891 gmx_fatal(FARGS, "Cannot find executable gunzip in /bin or /usr/bin.\n"
1892 "You may want to define the path to gunzip "
1893 "with the environment variable GMX_PATH_GZIP.", gunzip);
1898 sprintf(gunzip, "%s/gunzip", Path);
1899 if (!gmx_fexist(gunzip))
1901 gmx_fatal(FARGS, "Cannot find executable %s. Please define the path to gunzip"
1902 " in the environmental varialbe GMX_PATH_GZIP.", gunzip);
1907 printf("Using gunzip executable %s\n", gunzip);
1910 if (!gmx_fexist(fn))
1912 gmx_fatal(FARGS, "File %s does not exist.\n", fn);
1914 sprintf(Buffer, "%s -c < %s", gunzip, fn);
1917 printf("Executing command '%s'\n", Buffer);
1920 if ((pipe = popen(Buffer, "r")) == NULL)
1922 gmx_fatal(FARGS, "Unable to open pipe to `%s'\n", Buffer);
1925 gmx_fatal(FARGS, "Cannot open a compressed file on platform without pipe support");
1931 pipe = gmx_ffopen(fn, "r");
1938 //! Close file or pipe
1939 void pdo_close_file(FILE *fp)
1948 //! Reading all pdo files
1949 void read_pdo_files(char **fn, int nfiles, t_UmbrellaHeader* header,
1950 t_UmbrellaWindow *window, t_UmbrellaOptions *opt)
1953 real mintmp, maxtmp, done = 0.;
1956 /* char Buffer0[1000]; */
1960 gmx_fatal(FARGS, "No files found. Hick.");
1963 /* if min and max are not given, get min and max from the input files */
1966 printf("Automatic determination of boundaries from %d pdo files...\n", nfiles);
1969 for (i = 0; i < nfiles; ++i)
1971 file = open_pdo_pipe(fn[i], opt, &bPipeOpen);
1972 /*fgets(Buffer0,999,file);
1973 fprintf(stderr,"First line '%s'\n",Buffer0); */
1974 done = 100.0*(i+1)/nfiles;
1975 fprintf(stdout, "\rOpening %s ... [%2.0f%%]", fn[i], done); fflush(stdout);
1980 read_pdo_header(file, header, opt);
1981 /* here only determine min and max of this window */
1982 read_pdo_data(file, header, i, NULL, opt, TRUE, &mintmp, &maxtmp);
1983 if (maxtmp > opt->max)
1987 if (mintmp < opt->min)
1993 pdo_close_file(file);
2001 printf("\nDetermined boundaries to %f and %f\n\n", opt->min, opt->max);
2002 if (opt->bBoundsOnly)
2004 printf("Found option -boundsonly, now exiting.\n");
2008 /* store stepsize in profile */
2009 opt->dz = (opt->max-opt->min)/opt->bins;
2011 /* Having min and max, we read in all files */
2012 /* Loop over all files */
2013 for (i = 0; i < nfiles; ++i)
2015 done = 100.0*(i+1)/nfiles;
2016 fprintf(stdout, "\rOpening %s ... [%2.0f%%]", fn[i], done); fflush(stdout);
2021 file = open_pdo_pipe(fn[i], opt, &bPipeOpen);
2022 read_pdo_header(file, header, opt);
2023 /* load data into window */
2024 read_pdo_data(file, header, i, window, opt, FALSE, NULL, NULL);
2025 if ((window+i)->Ntot[0] == 0)
2027 fprintf(stderr, "\nWARNING, no data points read from file %s (check -b option)\n", fn[i]);
2031 pdo_close_file(file);
2039 for (i = 0; i < nfiles; ++i)
2046 //! translate 0/1 to N/Y to write pull dimensions
2047 #define int2YN(a) (((a) == 0) ? ("N") : ("Y"))
2049 //! Read pull groups from a tpr file (including position, force const, geometry, number of groups)
2050 void read_tpr_header(const char *fn, t_UmbrellaHeader* header, t_UmbrellaOptions *opt, t_coordselection *coordsel)
2055 static int first = 1;
2057 /* printf("Reading %s \n",fn); */
2058 read_tpx_state(fn, &ir, &state, NULL);
2062 gmx_fatal(FARGS, "This is not a tpr with COM pulling");
2064 if (ir.pull->ncoord == 0)
2066 gmx_fatal(FARGS, "No pull coordinates found in %s", fn);
2069 /* Read overall pull info */
2070 header->npullcrds = ir.pull->ncoord;
2071 header->bPrintCOM = ir.pull->bPrintCOM;
2072 header->bPrintRefValue = ir.pull->bPrintRefValue;
2073 header->bPrintComp = ir.pull->bPrintComp;
2075 /* Read pull coordinates */
2076 snew(header->pcrd, header->npullcrds);
2077 for (i = 0; i < ir.pull->ncoord; i++)
2079 header->pcrd[i].pull_type = ir.pull->coord[i].eType;
2080 header->pcrd[i].geometry = ir.pull->coord[i].eGeom;
2081 header->pcrd[i].ngroup = ir.pull->coord[i].ngroup;
2082 header->pcrd[i].k = ir.pull->coord[i].k;
2083 header->pcrd[i].init_dist = ir.pull->coord[i].init;
2085 copy_ivec(ir.pull->coord[i].dim, header->pcrd[i].dim);
2086 header->pcrd[i].ndim = header->pcrd[i].dim[XX] + header->pcrd[i].dim[YY] + header->pcrd[i].dim[ZZ];
2088 std::strcpy(header->pcrd[i].coord_unit,
2089 pull_coordinate_units(&ir.pull->coord[i]));
2091 if (ir.efep != efepNO && ir.pull->coord[i].k != ir.pull->coord[i].kB)
2093 gmx_fatal(FARGS, "Seems like you did free-energy perturbation, and you perturbed the force constant."
2094 " This is not supported.\n");
2096 if (coordsel && (coordsel->n != ir.pull->ncoord))
2098 gmx_fatal(FARGS, "Found %d pull coordinates in %s, but %d columns in the respective line\n"
2099 "coordinate selection file (option -is)\n", ir.pull->ncoord, fn, coordsel->n);
2103 /* Check pull coords for consistency */
2105 ivec thedim = { 0, 0, 0 };
2106 bool geometryIsSet = false;
2107 for (i = 0; i < ir.pull->ncoord; i++)
2109 if (coordsel == NULL || coordsel->bUse[i])
2111 if (header->pcrd[i].pull_type != epullUMBRELLA)
2113 gmx_fatal(FARGS, "%s: Pull coordinate %d is of type \"%s\", expected \"umbrella\". Only umbrella coodinates can enter WHAM.\n"
2114 "If you have umrella and non-umbrella coordinates, you can select the umbrella coordinates with gmx wham -is\n",
2115 fn, i+1, epull_names[header->pcrd[i].pull_type]);
2119 geom = header->pcrd[i].geometry;
2120 copy_ivec(header->pcrd[i].dim, thedim);
2121 geometryIsSet = true;
2123 if (geom != header->pcrd[i].geometry)
2125 gmx_fatal(FARGS, "%s: Your pull coordinates have different pull geometry (coordinate 1: %s, coordinate %d: %s)\n"
2126 "If you want to use only some pull coordinates in WHAM, please select them with option gmx wham -is\n",
2127 fn, epullg_names[geom], i+1, epullg_names[header->pcrd[i].geometry]);
2129 if (thedim[XX] != header->pcrd[i].dim[XX] || thedim[YY] != header->pcrd[i].dim[YY] || thedim[ZZ] != header->pcrd[i].dim[ZZ])
2131 gmx_fatal(FARGS, "%s: Your pull coordinates have different pull dimensions (coordinate 1: %s %s %s, coordinate %d: %s %s %s)\n"
2132 "If you want to use only some pull coordinates in WHAM, please select them with option gmx wham -is\n",
2133 fn, int2YN(thedim[XX]), int2YN(thedim[YY]), int2YN(thedim[ZZ]), i+1,
2134 int2YN(header->pcrd[i].dim[XX]), int2YN(header->pcrd[i].dim[YY]), int2YN(header->pcrd[i].dim[ZZ]));
2136 if (header->pcrd[i].geometry == epullgCYL)
2138 if (header->pcrd[i].dim[XX] || header->pcrd[i].dim[YY] || (!header->pcrd[i].dim[ZZ]))
2140 gmx_fatal(FARGS, "With pull geometry 'cylinder', expected pulling in Z direction only.\n"
2141 "However, found dimensions [%s %s %s]\n",
2142 int2YN(header->pcrd[i].dim[XX]), int2YN(header->pcrd[i].dim[YY]),
2143 int2YN(header->pcrd[i].dim[ZZ]));
2146 if (header->pcrd[i].k <= 0.0)
2148 gmx_fatal(FARGS, "%s: Pull coordinate %d has force constant of of %g in %s.\n"
2149 "That doesn't seem to be an Umbrella tpr.\n",
2150 fn, i+1, header->pcrd[i].k);
2155 if (opt->verbose || first)
2157 printf("\nFile %s, %d coordinates, with these options:\n", fn, header->npullcrds);
2159 for (i = 0; i < ir.pull->ncoord; i++)
2161 int lentmp = strlen(epullg_names[header->pcrd[i].geometry]);
2162 maxlen = (lentmp > maxlen) ? lentmp : maxlen;
2165 sprintf(fmt, "\tGeometry %%-%ds k = %%-8g position = %%-8g dimensions [%%s %%s %%s] (%%d dimensions). Used: %%s\n",
2167 for (i = 0; i < ir.pull->ncoord; i++)
2169 bool use = (coordsel == NULL || coordsel->bUse[i]);
2171 epullg_names[header->pcrd[i].geometry], header->pcrd[i].k, header->pcrd[i].init_dist,
2172 int2YN(header->pcrd[i].dim[XX]), int2YN(header->pcrd[i].dim[YY]), int2YN(header->pcrd[i].dim[ZZ]),
2173 header->pcrd[i].ndim, use ? "Yes" : "No");
2174 printf("\tPull group coordinates of %d groups expected in pullx files.\n", ir.pull->bPrintCOM ? header->pcrd[i].ngroup : 0);
2176 printf("\tReference value of the coordinate%s expected in pullx files.\n",
2177 header->bPrintRefValue ? "" : " not");
2179 if (!opt->verbose && first)
2181 printf("\tUse option -v to see this output for all input tpr files\n\n");
2187 //! 2-norm in a ndim-dimensional space
2188 double dist_ndim(double **dx, int ndim, int line)
2192 for (i = 0; i < ndim; i++)
2194 r2 += gmx::square(dx[i][line]);
2196 return std::sqrt(r2);
2199 //! Read pullx.xvg or pullf.xvg
2200 void read_pull_xf(const char *fn, t_UmbrellaHeader * header,
2201 t_UmbrellaWindow * window,
2202 t_UmbrellaOptions *opt,
2203 gmx_bool bGetMinMax, real *mintmp, real *maxtmp,
2204 t_coordselection *coordsel)
2206 double **y = 0, pos = 0., t, force, time0 = 0., dt;
2207 int ny, nt, bins, ibin, i, g, gUsed, dstep = 1;
2208 int nColExpect, ntot, column;
2209 real min, max, minfound = 1e20, maxfound = -1e20;
2210 gmx_bool dt_ok, timeok;
2211 const char *quantity;
2212 const int blocklen = 4096;
2214 static gmx_bool bFirst = TRUE;
2217 * Data columns in pull output:
2218 * - in force output pullf.xvg:
2219 * No reference columns, one column per pull coordinate
2221 * - in position output pullx.xvg:
2222 * * optionally, ndim columns for COMs of all groups (depending on on mdp options pull-print-com);
2223 * * The displacement, always one column. Note: with pull-print-components = yes, the dx/dy/dz would
2224 * be written separately into pullx file, but this is not supported and throws an error below;
2225 * * optionally, the position of the reference coordinate (depending on pull-print-ref-value)
2228 if (header->bPrintComp && opt->bPullx)
2230 gmx_fatal(FARGS, "gmx wham cannot read pullx files if the components of the coordinate was written\n"
2231 "(mdp option pull-print-components). Provide the pull force files instead (with option -if).\n");
2234 int *nColThisCrd, *nColCOMCrd, *nColRefCrd;
2235 snew(nColThisCrd, header->npullcrds);
2236 snew(nColCOMCrd, header->npullcrds);
2237 snew(nColRefCrd, header->npullcrds);
2239 if (opt->bPullx == FALSE)
2241 /* pullf reading: simply one column per coordinate */
2242 for (g = 0; g < header->npullcrds; g++)
2251 /* pullx reading. Note explanation above. */
2252 for (g = 0; g < header->npullcrds; g++)
2254 nColRefCrd[g] = (header->bPrintRefValue ? 1 : 0);
2255 nColCOMCrd[g] = (header->bPrintCOM ? header->pcrd[g].ndim*header->pcrd[g].ngroup : 0);
2256 nColThisCrd[g] = 1 + nColCOMCrd[g] + nColRefCrd[g];
2260 nColExpect = 1; /* time column */
2261 for (g = 0; g < header->npullcrds; g++)
2263 nColExpect += nColThisCrd[g];
2266 /* read pullf or pullx. Could be optimized if min and max are given. */
2267 nt = read_xvg(fn, &y, &ny);
2269 /* Check consistency */
2270 quantity = opt->bPullx ? "position" : "force";
2273 gmx_fatal(FARGS, "Empty pull %s file %s\n", quantity, fn);
2275 if (bFirst || opt->verbose)
2277 printf("\nReading pull %s file %s, expecting %d columns:\n", quantity, fn, nColExpect);
2278 for (i = 0; i < header->npullcrds; i++)
2280 printf("\tColumns for pull coordinate %d\n", i+1);
2281 printf("\t\treaction coordinate: %d\n"
2282 "\t\tcenter-of-mass of groups: %d\n"
2283 "\t\treference position column: %s\n",
2284 1, nColCOMCrd[i], (header->bPrintRefValue ? "Yes" : "No"));
2286 printf("\tFound %d times in %s\n", nt, fn);
2289 if (nColExpect != ny)
2291 gmx_fatal(FARGS, "Expected %d columns (including time column) in %s, but found %d."
2292 " Maybe you confused options -if and -ix ?", nColExpect, fn, ny);
2302 window->dt = y[0][1]-y[0][0];
2304 else if (opt->nBootStrap && opt->tauBootStrap != 0.0)
2306 fprintf(stderr, "\n *** WARNING, Could not determine time step in %s\n", fn);
2309 /* Need to alocate memory and set up structure for windows */
2312 /* Use only groups selected with option -is file */
2313 if (header->npullcrds != coordsel->n)
2315 gmx_fatal(FARGS, "tpr file contains %d pull groups, but expected %d from group selection file\n",
2316 header->npullcrds, coordsel->n);
2318 window->nPull = coordsel->nUse;
2322 window->nPull = header->npullcrds;
2325 window->nBin = bins;
2326 snew(window->Histo, window->nPull);
2327 snew(window->z, window->nPull);
2328 snew(window->k, window->nPull);
2329 snew(window->pos, window->nPull);
2330 snew(window->N, window->nPull);
2331 snew(window->Ntot, window->nPull);
2332 snew(window->g, window->nPull);
2333 snew(window->bsWeight, window->nPull);
2334 window->bContrib = 0;
2336 if (opt->bCalcTauInt)
2338 snew(window->ztime, window->nPull);
2342 window->ztime = NULL;
2344 snew(lennow, window->nPull);
2346 for (g = 0; g < window->nPull; ++g)
2349 window->bsWeight [g] = 1.;
2351 window->Ntot [g] = 0;
2353 snew(window->Histo[g], bins);
2355 if (opt->bCalcTauInt)
2357 window->ztime[g] = NULL;
2361 /* Copying umbrella center and force const is more involved since not
2362 all pull groups from header (tpr file) may be used in window variable */
2363 for (g = 0, gUsed = 0; g < header->npullcrds; ++g)
2365 if (coordsel && (coordsel->bUse[g] == FALSE))
2369 window->k [gUsed] = header->pcrd[g].k;
2370 window->pos[gUsed] = header->pcrd[g].init_dist;
2375 { /* only determine min and max */
2378 min = max = bins = 0; /* Get rid of warnings */
2382 for (i = 0; i < nt; i++)
2384 /* Do you want that time frame? */
2385 t = 1.0/1000*( static_cast<gmx_int64_t> ((y[0][i]*1000) + 0.5)); /* round time to fs */
2387 /* get time step of pdo file and get dstep from opt->dt */
2397 dstep = static_cast<int>(opt->dt/dt+0.5);
2405 window->dt = dt*dstep;
2409 dt_ok = (i%dstep == 0);
2410 timeok = (dt_ok && t >= opt->tmin && t <= opt->tmax);
2412 printf(" time = %f, (tmin,tmax)=(%e,%e), dt_ok=%d timeok=%d\n",
2413 t,opt->tmin, opt->tmax, dt_ok,timeok); */
2416 /* Note: if coordsel == NULL:
2417 * all groups in pullf/x file are stored in this window, and gUsed == g
2418 * if coordsel != NULL:
2419 * only groups with coordsel.bUse[g]==TRUE are stored. gUsed is not always equal g
2422 for (g = 0; g < header->npullcrds; ++g)
2424 /* was this group selected for application in WHAM? */
2425 if (coordsel && (coordsel->bUse[g] == FALSE))
2433 /* y has 1 time column y[0] and one column per force y[1],...,y[nCrds] */
2435 pos = -force/header->pcrd[g].k + header->pcrd[g].init_dist;
2439 /* Pick the correct column index.
2440 Note that there is always exactly one displacement column.
2443 for (int j = 0; j < g; j++)
2445 column += nColThisCrd[j];
2450 /* printf("crd %d dpos %f poseq %f pos %f \n",g,dpos,poseq,pos); */
2464 if (gUsed >= window->nPull)
2466 gmx_fatal(FARGS, "gUsed too large (%d, nPull=%d). This error should have been caught before.\n",
2467 gUsed, window->nPull);
2470 if (opt->bCalcTauInt && !bGetMinMax)
2472 /* save time series for autocorrelation analysis */
2473 ntot = window->Ntot[gUsed];
2474 /* printf("i %d, ntot %d, lennow[g] = %d\n",i,ntot,lennow[g]); */
2475 if (ntot >= lennow[gUsed])
2477 lennow[gUsed] += blocklen;
2478 srenew(window->ztime[gUsed], lennow[gUsed]);
2480 window->ztime[gUsed][ntot] = pos;
2483 ibin = static_cast<int> (std::floor((pos-min)/(max-min)*bins));
2488 while ( (ibin += bins) < 0)
2493 else if (ibin >= bins)
2495 while ( (ibin -= bins) >= bins)
2501 if (ibin >= 0 && ibin < bins)
2503 window->Histo[gUsed][ibin] += 1.;
2506 window->Ntot[gUsed]++;
2510 else if (t > opt->tmax)
2514 printf("time %f larger than tmax %f, stop reading this pullx/pullf file\n", t, opt->tmax);
2526 for (i = 0; i < ny; i++)
2532 //! read pullf-files.dat or pullx-files.dat and tpr-files.dat
2533 void read_tpr_pullxf_files(char **fnTprs, char **fnPull, int nfiles,
2534 t_UmbrellaHeader* header,
2535 t_UmbrellaWindow *window, t_UmbrellaOptions *opt)
2538 real mintmp, maxtmp;
2540 printf("Reading %d tpr and pullf files\n", nfiles);
2542 /* min and max not given? */
2545 printf("Automatic determination of boundaries...\n");
2548 for (i = 0; i < nfiles; i++)
2550 if (whaminFileType(fnTprs[i]) != whamin_tpr)
2552 gmx_fatal(FARGS, "Expected the %d'th file in input file to be a tpr file\n", i);
2554 read_tpr_header(fnTprs[i], header, opt, (opt->nCoordsel > 0) ? &opt->coordsel[i] : NULL);
2555 if (whaminFileType(fnPull[i]) != whamin_pullxf)
2557 gmx_fatal(FARGS, "Expected the %d'th file in input file to be a xvg (pullx/pullf) file\n", i);
2559 read_pull_xf(fnPull[i], header, NULL, opt, TRUE, &mintmp, &maxtmp,
2560 (opt->nCoordsel > 0) ? &opt->coordsel[i] : NULL);
2561 if (maxtmp > opt->max)
2565 if (mintmp < opt->min)
2570 printf("\nDetermined boundaries to %f and %f\n\n", opt->min, opt->max);
2571 if (opt->bBoundsOnly)
2573 printf("Found option -boundsonly, now exiting.\n");
2577 /* store stepsize in profile */
2578 opt->dz = (opt->max-opt->min)/opt->bins;
2580 bool foundData = false;
2581 for (i = 0; i < nfiles; i++)
2583 if (whaminFileType(fnTprs[i]) != whamin_tpr)
2585 gmx_fatal(FARGS, "Expected the %d'th file in input file to be a tpr file\n", i);
2587 read_tpr_header(fnTprs[i], header, opt, (opt->nCoordsel > 0) ? &opt->coordsel[i] : NULL);
2588 if (whaminFileType(fnPull[i]) != whamin_pullxf)
2590 gmx_fatal(FARGS, "Expected the %d'th file in input file to be a xvg (pullx/pullf) file\n", i);
2592 read_pull_xf(fnPull[i], header, window+i, opt, FALSE, NULL, NULL,
2593 (opt->nCoordsel > 0) ? &opt->coordsel[i] : NULL);
2594 if (window[i].Ntot[0] == 0)
2596 fprintf(stderr, "\nWARNING, no data points read from file %s (check -b option)\n", fnPull[i]);
2605 gmx_fatal(FARGS, "No data points were found in pullf/pullx files. Maybe you need to specify the -b option?\n");
2608 for (i = 0; i < nfiles; i++)
2617 /*! \brief Read integrated autocorrelation times from input file (option -iiact)
2619 * Note: Here we consider tau[int] := int_0^inf ACF(t) as the integrated autocorrelation times.
2620 * The factor `g := 1 + 2*tau[int]` subsequently enters the uncertainty.
2622 void readIntegratedAutocorrelationTimes(t_UmbrellaWindow *window, int nwins, const char* fn)
2624 int nlines, ny, i, ig;
2627 printf("Readging Integrated autocorrelation times from %s ...\n", fn);
2628 nlines = read_xvg(fn, &iact, &ny);
2629 if (nlines != nwins)
2631 gmx_fatal(FARGS, "Found %d lines with integrated autocorrelation times in %s.\nExpected %d",
2634 for (i = 0; i < nlines; i++)
2636 if (window[i].nPull != ny)
2638 gmx_fatal(FARGS, "You are providing autocorrelation times with option -iiact and the\n"
2639 "number of pull groups is different in different simulations. That is not\n"
2640 "supported yet. Sorry.\n");
2642 for (ig = 0; ig < window[i].nPull; ig++)
2644 /* compare Kumar et al, J Comp Chem 13, 1011-1021 (1992) */
2645 window[i].g[ig] = 1+2*iact[ig][i]/window[i].dt;
2647 if (iact[ig][i] <= 0.0)
2649 fprintf(stderr, "\nWARNING, IACT = %f (window %d, group %d)\n", iact[ig][i], i, ig);
2656 /*! \brief Smooth autocorreltion times along the reaction coordinate.
2659 * if the ACT is subject to high uncertainty in case if limited sampling. Note
2660 * that -in case of limited sampling- the ACT may be severely underestimated.
2661 * Note: the g=1+2tau are overwritten.
2662 * If opt->bAllowReduceIact==FALSE, the ACTs are never reduced, only increased
2665 void smoothIact(t_UmbrellaWindow *window, int nwins, t_UmbrellaOptions *opt)
2668 double pos, dpos2, siglim, siglim2, gaufact, invtwosig2, w, weight, tausm;
2670 /* only evaluate within +- 3sigma of the Gausian */
2671 siglim = 3.0*opt->sigSmoothIact;
2672 siglim2 = gmx::square(siglim);
2673 /* pre-factor of Gaussian */
2674 gaufact = 1.0/(std::sqrt(2*M_PI)*opt->sigSmoothIact);
2675 invtwosig2 = 0.5/gmx::square(opt->sigSmoothIact);
2677 for (i = 0; i < nwins; i++)
2679 snew(window[i].tausmooth, window[i].nPull);
2680 for (ig = 0; ig < window[i].nPull; ig++)
2684 pos = window[i].pos[ig];
2685 for (j = 0; j < nwins; j++)
2687 for (jg = 0; jg < window[j].nPull; jg++)
2689 dpos2 = gmx::square(window[j].pos[jg]-pos);
2690 if (dpos2 < siglim2)
2692 w = gaufact*std::exp(-dpos2*invtwosig2);
2694 tausm += w*window[j].tau[jg];
2695 /*printf("Weight %g dpos2=%g pos=%g gaufact=%g invtwosig2=%g\n",
2696 w,dpos2,pos,gaufact,invtwosig2); */
2701 if (opt->bAllowReduceIact || tausm > window[i].tau[ig])
2703 window[i].tausmooth[ig] = tausm;
2707 window[i].tausmooth[ig] = window[i].tau[ig];
2709 window[i].g[ig] = 1+2*tausm/window[i].dt;
2714 //! Stop integrating autoccorelation function when ACF drops under this value
2715 #define WHAM_AC_ZERO_LIMIT 0.05
2717 /*! \brief Try to compute the autocorrelation time for each umbrealla window
2719 void calcIntegratedAutocorrelationTimes(t_UmbrellaWindow *window, int nwins,
2720 t_UmbrellaOptions *opt, const char *fn, const char *xlabel)
2722 int i, ig, ncorr, ntot, j, k, *count, restart;
2723 real *corr, c0, dt, tmp;
2724 real *ztime, av, tausteps;
2725 FILE *fp, *fpcorr = 0;
2729 fpcorr = xvgropen("hist_autocorr.xvg", "Autocorrelation functions of umbrella windows",
2730 "time [ps]", "autocorrelation function", opt->oenv);
2734 for (i = 0; i < nwins; i++)
2736 fprintf(stdout, "\rEstimating integrated autocorrelation times ... [%2.0f%%] ...", 100.*(i+1)/nwins);
2738 ntot = window[i].Ntot[0];
2740 /* using half the maximum time as length of autocorrelation function */
2744 gmx_fatal(FARGS, "Tryig to estimtate autocorrelation time from only %d"
2745 " points. Provide more pull data!", ntot);
2748 /* snew(corrSq,ncorr); */
2751 snew(window[i].tau, window[i].nPull);
2752 restart = static_cast<int>(opt->acTrestart/dt+0.5);
2758 for (ig = 0; ig < window[i].nPull; ig++)
2760 if (ntot != window[i].Ntot[ig])
2762 gmx_fatal(FARGS, "Encountered different nr of frames in different pull groups.\n"
2763 "That should not happen. (%d and %d)\n", ntot, window[i].Ntot[ig]);
2765 ztime = window[i].ztime[ig];
2767 /* calc autocorrelation function C(t) = < [z(tau)-<z>]*[z(tau+t)-<z>]> */
2768 for (j = 0, av = 0; (j < ntot); j++)
2773 for (k = 0; (k < ncorr); k++)
2778 for (j = 0; (j < ntot); j += restart)
2780 for (k = 0; (k < ncorr) && (j+k < ntot); k++)
2782 tmp = (ztime[j]-av)*(ztime[j+k]-av);
2784 /* corrSq[k] += tmp*tmp; */
2788 /* divide by nr of frames for each time displacement */
2789 for (k = 0; (k < ncorr); k++)
2791 /* count probably = (ncorr-k+(restart-1))/restart; */
2792 corr[k] = corr[k]/count[k];
2793 /* variance of autocorrelation function */
2794 /* corrSq[k]=corrSq[k]/count[k]; */
2796 /* normalize such that corr[0] == 0 */
2798 for (k = 0; (k < ncorr); k++)
2801 /* corrSq[k]*=c0*c0; */
2804 /* write ACFs in verbose mode */
2807 for (k = 0; (k < ncorr); k++)
2809 fprintf(fpcorr, "%g %g\n", k*dt, corr[k]);
2811 fprintf(fpcorr, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
2814 /* esimate integrated correlation time, fitting is too unstable */
2815 tausteps = 0.5*corr[0];
2816 /* consider corr below WHAM_AC_ZERO_LIMIT as noise */
2817 for (j = 1; (j < ncorr) && (corr[j] > WHAM_AC_ZERO_LIMIT); j++)
2819 tausteps += corr[j];
2822 /* g = 1+2*tau, see. Ferrenberg/Swendsen, PRL 63:1195 (1989) or
2823 Kumar et al, eq. 28 ff. */
2824 window[i].tau[ig] = tausteps*dt;
2825 window[i].g[ig] = 1+2*tausteps;
2826 /* printf("win %d, group %d, estimated correlation time = %g ps\n",i,ig,window[i].tau[ig]); */
2837 /* plot IACT along reaction coordinate */
2838 fp = xvgropen(fn, "Integrated autocorrelation times", xlabel, "IACT [ps]", opt->oenv);
2839 if (output_env_get_print_xvgr_codes(opt->oenv))
2841 fprintf(fp, "@ s0 symbol 1\n@ s0 symbol size 0.5\n@ s0 line linestyle 0\n");
2842 fprintf(fp, "# WIN tau(gr1) tau(gr2) ...\n");
2843 for (i = 0; i < nwins; i++)
2845 fprintf(fp, "# %3d ", i);
2846 for (ig = 0; ig < window[i].nPull; ig++)
2848 fprintf(fp, " %11g", window[i].tau[ig]);
2853 for (i = 0; i < nwins; i++)
2855 for (ig = 0; ig < window[i].nPull; ig++)
2857 fprintf(fp, "%8g %8g\n", window[i].pos[ig], window[i].tau[ig]);
2860 if (opt->sigSmoothIact > 0.0)
2862 printf("Smoothing autocorrelation times along reaction coordinate with Gaussian of sig = %g\n",
2863 opt->sigSmoothIact);
2864 /* smooth IACT along reaction coordinate and overwrite g=1+2tau */
2865 smoothIact(window, nwins, opt);
2866 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(opt->oenv) ? "&" : "");
2867 if (output_env_get_print_xvgr_codes(opt->oenv))
2869 fprintf(fp, "@ s1 symbol 1\n@ s1 symbol size 0.5\n@ s1 line linestyle 0\n");
2870 fprintf(fp, "@ s1 symbol color 2\n");
2872 for (i = 0; i < nwins; i++)
2874 for (ig = 0; ig < window[i].nPull; ig++)
2876 fprintf(fp, "%8g %8g\n", window[i].pos[ig], window[i].tausmooth[ig]);
2881 printf("Wrote %s\n", fn);
2885 * compute average and sigma of each umbrella histogram
2887 void averageSigma(t_UmbrellaWindow *window, int nwins)
2890 real av, sum2, sig, diff, *ztime, nSamplesIndep;
2892 for (i = 0; i < nwins; i++)
2894 snew(window[i].aver, window[i].nPull);
2895 snew(window[i].sigma, window[i].nPull);
2897 ntot = window[i].Ntot[0];
2898 for (ig = 0; ig < window[i].nPull; ig++)
2900 ztime = window[i].ztime[ig];
2901 for (k = 0, av = 0.; k < ntot; k++)
2906 for (k = 0, sum2 = 0.; k < ntot; k++)
2911 sig = std::sqrt(sum2/ntot);
2912 window[i].aver[ig] = av;
2914 /* Note: This estimate for sigma is biased from the limited sampling.
2915 Correct sigma by n/(n-1) where n = number of independent
2916 samples. Only possible if IACT is known.
2920 nSamplesIndep = window[i].N[ig]/(window[i].tau[ig]/window[i].dt);
2921 window[i].sigma[ig] = sig * nSamplesIndep/(nSamplesIndep-1);
2925 window[i].sigma[ig] = sig;
2927 printf("win %d, aver = %f sig = %f\n", i, av, window[i].sigma[ig]);
2934 * Use histograms to compute average force on pull group.
2936 void computeAverageForce(t_UmbrellaWindow *window, int nWindows, t_UmbrellaOptions *opt)
2938 int i, j, bins = opt->bins, k;
2939 double dz, min = opt->min, max = opt->max, displAv, temp, distance, ztot, ztot_half, w, weight;
2942 dz = (max-min)/bins;
2943 ztot = opt->max-min;
2946 /* Compute average displacement from histograms */
2947 for (j = 0; j < nWindows; ++j)
2949 snew(window[j].forceAv, window[j].nPull);
2950 for (k = 0; k < window[j].nPull; ++k)
2954 for (i = 0; i < opt->bins; ++i)
2956 temp = (1.0*i+0.5)*dz+min;
2957 distance = temp - window[j].pos[k];
2959 { /* in cyclic wham: */
2960 if (distance > ztot_half) /* |distance| < ztot_half */
2964 else if (distance < -ztot_half)
2969 w = window[j].Histo[k][i]/window[j].g[k];
2970 displAv += w*distance;
2972 /* Are we near min or max? We are getting wrong forces from the histgrams since
2973 the histograms are zero outside [min,max). Therefore, assume that the position
2974 on the other side of the histomgram center is equally likely. */
2977 posmirrored = window[j].pos[k]-distance;
2978 if (posmirrored >= max || posmirrored < min)
2980 displAv += -w*distance;
2987 /* average force from average displacement */
2988 window[j].forceAv[k] = displAv*window[j].k[k];
2989 /* sigma from average square displacement */
2990 /* window[j].sigma [k] = sqrt(displAv2); */
2991 /* printf("Win %d, sigma = %f\n",j,sqrt(displAv2)); */
2997 * Check if the complete reaction coordinate is covered by the histograms
2999 void checkReactionCoordinateCovered(t_UmbrellaWindow *window, int nwins,
3000 t_UmbrellaOptions *opt)
3002 int i, ig, j, bins = opt->bins, bBoundary;
3003 real avcount = 0, z, relcount, *count;
3004 snew(count, opt->bins);
3006 for (j = 0; j < opt->bins; ++j)
3008 for (i = 0; i < nwins; i++)
3010 for (ig = 0; ig < window[i].nPull; ig++)
3012 count[j] += window[i].Histo[ig][j];
3015 avcount += 1.0*count[j];
3018 for (j = 0; j < bins; ++j)
3020 relcount = count[j]/avcount;
3021 z = (j+0.5)*opt->dz+opt->min;
3022 bBoundary = ( j<bins/20 || (bins-j)>bins/20 );
3023 /* check for bins with no data */
3026 fprintf(stderr, "\nWARNING, no data point in bin %d (z=%g) !\n"
3027 "You may not get a reasonable profile. Check your histograms!\n", j, z);
3029 /* and check for poor sampling */
3030 else if (relcount < 0.005 && !bBoundary)
3032 fprintf(stderr, "Warning, poor sampling bin %d (z=%g). Check your histograms!\n", j, z);
3038 /*! \brief Compute initial potential by integrating the average force
3040 * This speeds up the convergence by roughly a factor of 2
3042 void guessPotByIntegration(t_UmbrellaWindow *window, int nWindows, t_UmbrellaOptions *opt, const char *xlabel)
3044 int i, j, ig, bins = opt->bins, nHist, winmin, groupmin;
3045 double dz, min = opt->min, *pot, pos, hispos, dist, diff, fAv, distmin, *f;
3048 dz = (opt->max-min)/bins;
3050 printf("Getting initial potential by integration.\n");
3052 /* Compute average displacement from histograms */
3053 computeAverageForce(window, nWindows, opt);
3055 /* Get force for each bin from all histograms in this bin, or, alternatively,
3056 if no histograms are inside this bin, from the closest histogram */
3059 for (j = 0; j < opt->bins; ++j)
3061 pos = (1.0*j+0.5)*dz+min;
3065 groupmin = winmin = 0;
3066 for (i = 0; i < nWindows; i++)
3068 for (ig = 0; ig < window[i].nPull; ig++)
3070 hispos = window[i].pos[ig];
3071 dist = std::abs(hispos-pos);
3072 /* average force within bin */
3076 fAv += window[i].forceAv[ig];
3078 /* at the same time, remember closest histogram */
3087 /* if no histogram found in this bin, use closest histogram */
3094 fAv = window[winmin].forceAv[groupmin];
3098 for (j = 1; j < opt->bins; ++j)
3100 pot[j] = pot[j-1] - 0.5*dz*(f[j-1]+f[j]);
3103 /* cyclic wham: linearly correct possible offset */
3106 diff = (pot[bins-1]-pot[0])/(bins-1);
3107 for (j = 1; j < opt->bins; ++j)
3114 fp = xvgropen("pmfintegrated.xvg", "PMF from force integration", xlabel, "PMF (kJ/mol)", opt->oenv);
3115 for (j = 0; j < opt->bins; ++j)
3117 fprintf(fp, "%g %g\n", (j+0.5)*dz+opt->min, pot[j]);
3120 printf("verbose mode: wrote %s with PMF from interated forces\n", "pmfintegrated.xvg");
3123 /* get initial z=exp(-F[i]/kT) from integrated potential, where F[i] denote the free
3124 energy offsets which are usually determined by wham
3125 First: turn pot into probabilities:
3127 for (j = 0; j < opt->bins; ++j)
3129 pot[j] = std::exp(-pot[j]/(BOLTZ*opt->Temperature));
3131 calc_z(pot, window, nWindows, opt, TRUE);
3137 //! Count number of words in an line
3138 static int wordcount(char *ptr)
3143 if (std::strlen(ptr) == 0)
3147 /* fprintf(stderr,"ptr='%s'\n",ptr); */
3149 for (i = 0; (ptr[i] != '\0'); i++)
3151 is[cur] = isspace(ptr[i]);
3152 if ((i > 0) && (is[cur] && !is[1-cur]))
3161 /*! \brief Read input file for pull group selection (option -is)
3163 * TO DO: ptr=fgets(...) is never freed (small memory leak)
3165 void readPullCoordSelection(t_UmbrellaOptions *opt, char **fnTpr, int nTpr)
3168 int i, iline, n, len = STRLEN, temp;
3169 char *ptr = 0, *tmpbuf = 0;
3170 char fmt[1024], fmtign[1024];
3171 int block = 1, sizenow;
3173 fp = gmx_ffopen(opt->fnCoordSel, "r");
3174 opt->coordsel = NULL;
3179 while ( (ptr = fgets3(fp, tmpbuf, &len)) != NULL)
3184 if (iline >= sizenow)
3187 srenew(opt->coordsel, sizenow);
3189 opt->coordsel[iline].n = n;
3190 opt->coordsel[iline].nUse = 0;
3191 snew(opt->coordsel[iline].bUse, n);
3194 for (i = 0; i < n; i++)
3196 std::strcpy(fmt, fmtign);
3197 std::strcat(fmt, "%d");
3198 if (sscanf(ptr, fmt, &temp))
3200 opt->coordsel[iline].bUse[i] = (temp > 0);
3201 if (opt->coordsel[iline].bUse[i])
3203 opt->coordsel[iline].nUse++;
3206 std::strcat(fmtign, "%*s");
3210 opt->nCoordsel = iline;
3211 if (nTpr != opt->nCoordsel)
3213 gmx_fatal(FARGS, "Found %d tpr files but %d lines in %s\n", nTpr, opt->nCoordsel,
3217 printf("\nUse only these pull coordinates:\n");
3218 for (iline = 0; iline < nTpr; iline++)
3220 printf("%s (%d of %d coordinates):", fnTpr[iline], opt->coordsel[iline].nUse, opt->coordsel[iline].n);
3221 for (i = 0; i < opt->coordsel[iline].n; i++)
3223 if (opt->coordsel[iline].bUse[i])
3236 #define WHAMBOOLXOR(a, b) ( ((!(a)) && (b)) || ((a) && (!(b))))
3238 //! Number of elements in fnm (used for command line parsing)
3239 #define NFILE asize(fnm)
3241 //! The main gmx wham routine
3242 int gmx_wham(int argc, char *argv[])
3244 const char *desc[] = {
3245 "[THISMODULE] is an analysis program that implements the Weighted",
3246 "Histogram Analysis Method (WHAM). It is intended to analyze",
3247 "output files generated by umbrella sampling simulations to ",
3248 "compute a potential of mean force (PMF).[PAR]",
3250 "[THISMODULE] is currently not fully up to date. It only supports pull setups",
3251 "where the first pull coordinate(s) is/are umbrella pull coordinates",
3252 "and, if multiple coordinates need to be analyzed, all used the same",
3253 "geometry and dimensions. In most cases this is not an issue.[PAR]",
3254 "At present, three input modes are supported.",
3256 "* With option [TT]-it[tt], the user provides a file which contains the",
3257 " file names of the umbrella simulation run-input files ([REF].tpr[ref] files),",
3258 " AND, with option [TT]-ix[tt], a file which contains file names of",
3259 " the pullx [TT]mdrun[tt] output files. The [REF].tpr[ref] and pullx files must",
3260 " be in corresponding order, i.e. the first [REF].tpr[ref] created the",
3261 " first pullx, etc.",
3262 "* Same as the previous input mode, except that the the user",
3263 " provides the pull force output file names ([TT]pullf.xvg[tt]) with option [TT]-if[tt].",
3264 " From the pull force the position in the umbrella potential is",
3265 " computed. This does not work with tabulated umbrella potentials.",
3266 "* With option [TT]-ip[tt], the user provides file names of (gzipped) [REF].pdo[ref] files, i.e.",
3267 " the GROMACS 3.3 umbrella output files. If you have some unusual"
3268 " reaction coordinate you may also generate your own [REF].pdo[ref] files and",
3269 " feed them with the [TT]-ip[tt] option into to [THISMODULE]. The [REF].pdo[ref] file header",
3270 " must be similar to the following::",
3273 " # Component selection: 0 0 1",
3275 " # Ref. Group 'TestAtom'",
3276 " # Nr. of pull groups 2",
3277 " # Group 1 'GR1' Umb. Pos. 5.0 Umb. Cons. 1000.0",
3278 " # Group 2 'GR2' Umb. Pos. 2.0 Umb. Cons. 500.0",
3281 " The number of pull groups, umbrella positions, force constants, and names ",
3282 " may (of course) differ. Following the header, a time column and ",
3283 " a data column for each pull group follows (i.e. the displacement",
3284 " with respect to the umbrella center). Up to four pull groups are possible ",
3285 " per [REF].pdo[ref] file at present.[PAR]",
3286 "By default, all pull coordinates found in all pullx/pullf files are used in WHAM. If only ",
3287 "some of the pull coordinates should be used, a pull coordinate selection file (option [TT]-is[tt]) can ",
3288 "be provided. The selection file must contain one line for each tpr file in tpr-files.dat.",
3289 "Each of these lines must contain one digit (0 or 1) for each pull coordinate in the tpr file. ",
3290 "Here, 1 indicates that the pull coordinate is used in WHAM, and 0 means it is omitted. Example:",
3291 "If you have three tpr files, each containing 4 pull coordinates, but only pull coordinates 1 and 2 should be ",
3292 "used, coordsel.dat looks like this::",
3298 "By default, the output files are::",
3300 " [TT]-o[tt] PMF output file",
3301 " [TT]-hist[tt] Histograms output file",
3303 "Always check whether the histograms sufficiently overlap.[PAR]",
3304 "The umbrella potential is assumed to be harmonic and the force constants are ",
3305 "read from the [REF].tpr[ref] or [REF].pdo[ref] files. If a non-harmonic umbrella force was applied ",
3306 "a tabulated potential can be provided with [TT]-tab[tt].",
3311 "* [TT]-bins[tt] Number of bins used in analysis",
3312 "* [TT]-temp[tt] Temperature in the simulations",
3313 "* [TT]-tol[tt] Stop iteration if profile (probability) changed less than tolerance",
3314 "* [TT]-auto[tt] Automatic determination of boundaries",
3315 "* [TT]-min,-max[tt] Boundaries of the profile",
3317 "The data points that are used to compute the profile",
3318 "can be restricted with options [TT]-b[tt], [TT]-e[tt], and [TT]-dt[tt]. ",
3319 "Adjust [TT]-b[tt] to ensure sufficient equilibration in each ",
3320 "umbrella window.[PAR]",
3321 "With [TT]-log[tt] (default) the profile is written in energy units, otherwise ",
3322 "(with [TT]-nolog[tt]) as probability. The unit can be specified with [TT]-unit[tt]. ",
3323 "With energy output, the energy in the first bin is defined to be zero. ",
3324 "If you want the free energy at a different ",
3325 "position to be zero, set [TT]-zprof0[tt] (useful with bootstrapping, see below).[PAR]",
3326 "For cyclic or periodic reaction coordinates (dihedral angle, channel PMF",
3327 "without osmotic gradient), the option [TT]-cycl[tt] is useful.",
3328 "[THISMODULE] will make use of the",
3329 "periodicity of the system and generate a periodic PMF. The first and the last bin of the",
3330 "reaction coordinate will assumed be be neighbors.[PAR]",
3331 "Option [TT]-sym[tt] symmetrizes the profile around z=0 before output, ",
3332 "which may be useful for, e.g. membranes.",
3337 "If available, the number of OpenMP threads used by gmx wham can be controlled by setting",
3338 "the [TT]OMP_NUM_THREADS[tt] environment variable.",
3343 "With [TT]-ac[tt], [THISMODULE] estimates the integrated autocorrelation ",
3344 "time (IACT) [GRK]tau[grk] for each umbrella window and weights the respective ",
3345 "window with 1/[1+2*[GRK]tau[grk]/dt]. The IACTs are written ",
3346 "to the file defined with [TT]-oiact[tt]. In verbose mode, all ",
3347 "autocorrelation functions (ACFs) are written to [TT]hist_autocorr.xvg[tt]. ",
3348 "Because the IACTs can be severely underestimated in case of limited ",
3349 "sampling, option [TT]-acsig[tt] allows one to smooth the IACTs along the ",
3350 "reaction coordinate with a Gaussian ([GRK]sigma[grk] provided with [TT]-acsig[tt], ",
3351 "see output in [TT]iact.xvg[tt]). Note that the IACTs are estimated by simple ",
3352 "integration of the ACFs while the ACFs are larger 0.05.",
3353 "If you prefer to compute the IACTs by a more sophisticated (but possibly ",
3354 "less robust) method such as fitting to a double exponential, you can ",
3355 "compute the IACTs with [gmx-analyze] and provide them to [THISMODULE] with the file ",
3356 "[TT]iact-in.dat[tt] (option [TT]-iiact[tt]), which should contain one line per ",
3357 "input file ([REF].pdo[ref] or pullx/f file) and one column per pull coordinate in the respective file.",
3362 "Statistical errors may be estimated with bootstrap analysis. Use it with care, ",
3363 "otherwise the statistical error may be substantially underestimated. ",
3364 "More background and examples for the bootstrap technique can be found in ",
3365 "Hub, de Groot and Van der Spoel, JCTC (2010) 6: 3713-3720.",
3366 "[TT]-nBootstrap[tt] defines the number of bootstraps (use, e.g., 100). ",
3367 "Four bootstrapping methods are supported and ",
3368 "selected with [TT]-bs-method[tt].",
3370 "* [TT]b-hist[tt] Default: complete histograms are considered as independent ",
3371 " data points, and the bootstrap is carried out by assigning random weights to the ",
3372 " histograms (\"Bayesian bootstrap\"). Note that each point along the reaction coordinate",
3373 " must be covered by multiple independent histograms (e.g. 10 histograms), otherwise the ",
3374 " statistical error is underestimated.",
3375 "* [TT]hist[tt] Complete histograms are considered as independent data points. ",
3376 " For each bootstrap, N histograms are randomly chosen from the N given histograms ",
3377 " (allowing duplication, i.e. sampling with replacement).",
3378 " To avoid gaps without data along the reaction coordinate blocks of histograms ",
3379 " ([TT]-histbs-block[tt]) may be defined. In that case, the given histograms are ",
3380 " divided into blocks and only histograms within each block are mixed. Note that ",
3381 " the histograms within each block must be representative for all possible histograms, ",
3382 " otherwise the statistical error is underestimated.",
3383 "* [TT]traj[tt] The given histograms are used to generate new random trajectories,",
3384 " such that the generated data points are distributed according the given histograms ",
3385 " and properly autocorrelated. The autocorrelation time (ACT) for each window must be ",
3386 " known, so use [TT]-ac[tt] or provide the ACT with [TT]-iiact[tt]. If the ACT of all ",
3387 " windows are identical (and known), you can also provide them with [TT]-bs-tau[tt]. ",
3388 " Note that this method may severely underestimate the error in case of limited sampling, ",
3389 " that is if individual histograms do not represent the complete phase space at ",
3390 " the respective positions.",
3391 "* [TT]traj-gauss[tt] The same as method [TT]traj[tt], but the trajectories are ",
3392 " not bootstrapped from the umbrella histograms but from Gaussians with the average ",
3393 " and width of the umbrella histograms. That method yields similar error estimates ",
3394 " like method [TT]traj[tt].",
3396 "Bootstrapping output:",
3398 "* [TT]-bsres[tt] Average profile and standard deviations",
3399 "* [TT]-bsprof[tt] All bootstrapping profiles",
3401 "With [TT]-vbs[tt] (verbose bootstrapping), the histograms of each bootstrap are written, ",
3402 "and, with bootstrap method [TT]traj[tt], the cumulative distribution functions of ",
3406 const char *en_unit[] = {NULL, "kJ", "kCal", "kT", NULL};
3407 const char *en_unit_label[] = {"", "E (kJ mol\\S-1\\N)", "E (kcal mol\\S-1\\N)", "E (kT)", NULL};
3408 const char *en_bsMethod[] = { NULL, "b-hist", "hist", "traj", "traj-gauss", NULL };
3409 static t_UmbrellaOptions opt;
3412 { "-min", FALSE, etREAL, {&opt.min},
3413 "Minimum coordinate in profile"},
3414 { "-max", FALSE, etREAL, {&opt.max},
3415 "Maximum coordinate in profile"},
3416 { "-auto", FALSE, etBOOL, {&opt.bAuto},
3417 "Determine min and max automatically"},
3418 { "-bins", FALSE, etINT, {&opt.bins},
3419 "Number of bins in profile"},
3420 { "-temp", FALSE, etREAL, {&opt.Temperature},
3422 { "-tol", FALSE, etREAL, {&opt.Tolerance},
3424 { "-v", FALSE, etBOOL, {&opt.verbose},
3426 { "-b", FALSE, etREAL, {&opt.tmin},
3427 "First time to analyse (ps)"},
3428 { "-e", FALSE, etREAL, {&opt.tmax},
3429 "Last time to analyse (ps)"},
3430 { "-dt", FALSE, etREAL, {&opt.dt},
3431 "Analyse only every dt ps"},
3432 { "-histonly", FALSE, etBOOL, {&opt.bHistOnly},
3433 "Write histograms and exit"},
3434 { "-boundsonly", FALSE, etBOOL, {&opt.bBoundsOnly},
3435 "Determine min and max and exit (with [TT]-auto[tt])"},
3436 { "-log", FALSE, etBOOL, {&opt.bLog},
3437 "Calculate the log of the profile before printing"},
3438 { "-unit", FALSE, etENUM, {en_unit},
3439 "Energy unit in case of log output" },
3440 { "-zprof0", FALSE, etREAL, {&opt.zProf0},
3441 "Define profile to 0.0 at this position (with [TT]-log[tt])"},
3442 { "-cycl", FALSE, etBOOL, {&opt.bCycl},
3443 "Create cyclic/periodic profile. Assumes min and max are the same point."},
3444 { "-sym", FALSE, etBOOL, {&opt.bSym},
3445 "Symmetrize profile around z=0"},
3446 { "-hist-eq", FALSE, etBOOL, {&opt.bHistEq},
3447 "HIDDENEnforce equal weight for all histograms. (Non-Weighed-HAM)"},
3448 { "-ac", FALSE, etBOOL, {&opt.bCalcTauInt},
3449 "Calculate integrated autocorrelation times and use in wham"},
3450 { "-acsig", FALSE, etREAL, {&opt.sigSmoothIact},
3451 "Smooth autocorrelation times along reaction coordinate with Gaussian of this [GRK]sigma[grk]"},
3452 { "-ac-trestart", FALSE, etREAL, {&opt.acTrestart},
3453 "When computing autocorrelation functions, restart computing every .. (ps)"},
3454 { "-acred", FALSE, etBOOL, {&opt.bAllowReduceIact},
3455 "HIDDENWhen smoothing the ACTs, allows one to reduce ACTs. Otherwise, only increase ACTs "
3456 "during smoothing"},
3457 { "-nBootstrap", FALSE, etINT, {&opt.nBootStrap},
3458 "nr of bootstraps to estimate statistical uncertainty (e.g., 200)" },
3459 { "-bs-method", FALSE, etENUM, {en_bsMethod},
3460 "Bootstrap method" },
3461 { "-bs-tau", FALSE, etREAL, {&opt.tauBootStrap},
3462 "Autocorrelation time (ACT) assumed for all histograms. Use option [TT]-ac[tt] if ACT is unknown."},
3463 { "-bs-seed", FALSE, etINT, {&opt.bsSeed},
3464 "Seed for bootstrapping. (-1 = use time)"},
3465 { "-histbs-block", FALSE, etINT, {&opt.histBootStrapBlockLength},
3466 "When mixing histograms only mix within blocks of [TT]-histbs-block[tt]."},
3467 { "-vbs", FALSE, etBOOL, {&opt.bs_verbose},
3468 "Verbose bootstrapping. Print the CDFs and a histogram file for each bootstrap."},
3469 { "-stepout", FALSE, etINT, {&opt.stepchange},
3470 "HIDDENWrite maximum change every ... (set to 1 with [TT]-v[tt])"},
3471 { "-updateContr", FALSE, etINT, {&opt.stepUpdateContrib},
3472 "HIDDENUpdate table with significan contributions to WHAM every ... iterations"},
3476 { efDAT, "-ix", "pullx-files", ffOPTRD}, /* wham input: pullf.xvg's and tprs */
3477 { efDAT, "-if", "pullf-files", ffOPTRD}, /* wham input: pullf.xvg's and tprs */
3478 { efDAT, "-it", "tpr-files", ffOPTRD}, /* wham input: tprs */
3479 { efDAT, "-ip", "pdo-files", ffOPTRD}, /* wham input: pdo files (gmx3 style) */
3480 { efDAT, "-is", "coordsel", ffOPTRD}, /* input: select pull coords to use */
3481 { efXVG, "-o", "profile", ffWRITE }, /* output file for profile */
3482 { efXVG, "-hist", "histo", ffWRITE}, /* output file for histograms */
3483 { efXVG, "-oiact", "iact", ffOPTWR}, /* writing integrated autocorrelation times */
3484 { efDAT, "-iiact", "iact-in", ffOPTRD}, /* reading integrated autocorrelation times */
3485 { efXVG, "-bsres", "bsResult", ffOPTWR}, /* average and errors of bootstrap analysis */
3486 { efXVG, "-bsprof", "bsProfs", ffOPTWR}, /* output file for bootstrap profiles */
3487 { efDAT, "-tab", "umb-pot", ffOPTRD}, /* Tabulated umbrella potential (if not harmonic) */
3490 int i, j, l, nfiles, nwins, nfiles2;
3491 t_UmbrellaHeader header;
3492 t_UmbrellaWindow * window = NULL;
3493 double *profile, maxchange = 1e20;
3494 gmx_bool bMinSet, bMaxSet, bAutoSet, bExact = FALSE;
3495 char **fninTpr, **fninPull, **fninPdo;
3497 FILE *histout, *profout;
3498 char xlabel[STRLEN], ylabel[256], title[256];
3501 opt.verbose = FALSE;
3502 opt.bHistOnly = FALSE;
3511 opt.coordsel = NULL;
3513 /* bootstrapping stuff */
3515 opt.bsMethod = bsMethod_hist;
3516 opt.tauBootStrap = 0.0;
3518 opt.histBootStrapBlockLength = 8;
3519 opt.bs_verbose = FALSE;
3524 opt.Temperature = 298;
3525 opt.Tolerance = 1e-6;
3526 opt.bBoundsOnly = FALSE;
3528 opt.bCalcTauInt = FALSE;
3529 opt.sigSmoothIact = 0.0;
3530 opt.bAllowReduceIact = TRUE;
3531 opt.bInitPotByIntegration = TRUE;
3532 opt.acTrestart = 1.0;
3533 opt.stepchange = 100;
3534 opt.stepUpdateContrib = 100;
3536 if (!parse_common_args(&argc, argv, 0,
3537 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &opt.oenv))
3542 opt.unit = nenum(en_unit);
3543 opt.bsMethod = nenum(en_bsMethod);
3545 opt.bProf0Set = opt2parg_bSet("-zprof0", asize(pa), pa);
3547 opt.bTab = opt2bSet("-tab", NFILE, fnm);
3548 opt.bPdo = opt2bSet("-ip", NFILE, fnm);
3549 opt.bTpr = opt2bSet("-it", NFILE, fnm);
3550 opt.bPullx = opt2bSet("-ix", NFILE, fnm);
3551 opt.bPullf = opt2bSet("-if", NFILE, fnm);
3552 opt.bTauIntGiven = opt2bSet("-iiact", NFILE, fnm);
3553 if (opt.bTab && opt.bPullf)
3555 gmx_fatal(FARGS, "Force input does not work with tabulated potentials. "
3556 "Provide pullx.xvg or pdo files!");
3559 if (!opt.bPdo && !WHAMBOOLXOR(opt.bPullx, opt.bPullf))
3561 gmx_fatal(FARGS, "Give either pullx (-ix) OR pullf (-if) data. Not both.");
3563 if (!opt.bPdo && !(opt.bTpr || opt.bPullf || opt.bPullx))
3565 gmx_fatal(FARGS, "gmx wham supports three input modes, pullx, pullf, or pdo file input."
3566 "\n\n Check gmx wham -h !");
3569 opt.fnPdo = opt2fn("-ip", NFILE, fnm);
3570 opt.fnTpr = opt2fn("-it", NFILE, fnm);
3571 opt.fnPullf = opt2fn("-if", NFILE, fnm);
3572 opt.fnPullx = opt2fn("-ix", NFILE, fnm);
3573 opt.fnCoordSel = opt2fn_null("-is", NFILE, fnm);
3575 bMinSet = opt2parg_bSet("-min", asize(pa), pa);
3576 bMaxSet = opt2parg_bSet("-max", asize(pa), pa);
3577 bAutoSet = opt2parg_bSet("-auto", asize(pa), pa);
3578 if ( (bMinSet || bMaxSet) && bAutoSet)
3580 gmx_fatal(FARGS, "With -auto, do not give -min or -max\n");
3583 if ( (bMinSet && !bMaxSet) || (!bMinSet && bMaxSet))
3585 gmx_fatal(FARGS, "When giving -min, you must give -max (and vice versa), too\n");
3588 if (bMinSet && opt.bAuto)
3590 printf("Note: min and max given, switching off -auto.\n");
3594 if (opt.bTauIntGiven && opt.bCalcTauInt)
3596 gmx_fatal(FARGS, "Either read (option -iiact) or calculate (option -ac) the\n"
3597 "the autocorrelation times. Not both.");
3600 if (opt.tauBootStrap > 0.0 && opt2parg_bSet("-ac", asize(pa), pa))
3602 gmx_fatal(FARGS, "Either compute autocorrelation times (ACTs) (option -ac) or "
3603 "provide it with -bs-tau for bootstrapping. Not Both.\n");
3605 if (opt.tauBootStrap > 0.0 && opt2bSet("-iiact", NFILE, fnm))
3607 gmx_fatal(FARGS, "Either provide autocorrelation times (ACTs) with file iact-in.dat "
3608 "(option -iiact) or define all ACTs with -bs-tau for bootstrapping\n. Not Both.");
3611 /* Reading gmx4/gmx5 pull output and tpr files */
3612 if (opt.bTpr || opt.bPullf || opt.bPullx)
3614 read_wham_in(opt.fnTpr, &fninTpr, &nfiles, &opt);
3616 fnPull = opt.bPullf ? opt.fnPullf : opt.fnPullx;
3617 read_wham_in(fnPull, &fninPull, &nfiles2, &opt);
3618 printf("Found %d tpr and %d pull %s files in %s and %s, respectively\n",
3619 nfiles, nfiles2, opt.bPullf ? "force" : "position", opt.fnTpr, fnPull);
3620 if (nfiles != nfiles2)
3622 gmx_fatal(FARGS, "Found %d file names in %s, but %d in %s\n", nfiles,
3623 opt.fnTpr, nfiles2, fnPull);
3626 /* Read file that selects the pull group to be used */
3627 if (opt.fnCoordSel != NULL)
3629 readPullCoordSelection(&opt, fninTpr, nfiles);
3632 window = initUmbrellaWindows(nfiles);
3633 read_tpr_pullxf_files(fninTpr, fninPull, nfiles, &header, window, &opt);
3636 { /* reading pdo files */
3637 if (opt.fnCoordSel != NULL)
3639 gmx_fatal(FARGS, "Reading a -is file is not supported with PDO input files.\n"
3640 "Use awk or a similar tool to pick the required pull groups from your PDO files\n");
3642 read_wham_in(opt.fnPdo, &fninPdo, &nfiles, &opt);
3643 printf("Found %d pdo files in %s\n", nfiles, opt.fnPdo);
3644 window = initUmbrellaWindows(nfiles);
3645 read_pdo_files(fninPdo, nfiles, &header, window, &opt);
3648 /* It is currently assumed that all pull coordinates have the same geometry, so they also have the same coordinate units.
3649 We can therefore get the units for the xlabel from the first coordinate. */
3650 sprintf(xlabel, "\\xx\\f{} (%s)", header.pcrd[0].coord_unit);
3654 /* enforce equal weight for all histograms? */
3657 enforceEqualWeights(window, nwins);
3660 /* write histograms */
3661 histout = xvgropen(opt2fn("-hist", NFILE, fnm), "Umbrella histograms",
3662 xlabel, "count", opt.oenv);
3663 for (l = 0; l < opt.bins; ++l)
3665 fprintf(histout, "%e\t", (l+0.5)/opt.bins*(opt.max-opt.min)+opt.min);
3666 for (i = 0; i < nwins; ++i)
3668 for (j = 0; j < window[i].nPull; ++j)
3670 fprintf(histout, "%e\t", window[i].Histo[j][l]);
3673 fprintf(histout, "\n");
3676 printf("Wrote %s\n", opt2fn("-hist", NFILE, fnm));
3679 printf("Wrote histograms to %s, now exiting.\n", opt2fn("-hist", NFILE, fnm));
3683 /* Using tabulated umbrella potential */
3686 setup_tab(opt2fn("-tab", NFILE, fnm), &opt);
3689 /* Integrated autocorrelation times provided ? */
3690 if (opt.bTauIntGiven)
3692 readIntegratedAutocorrelationTimes(window, nwins, opt2fn("-iiact", NFILE, fnm));
3695 /* Compute integrated autocorrelation times */
3696 if (opt.bCalcTauInt)
3698 calcIntegratedAutocorrelationTimes(window, nwins, &opt, opt2fn("-oiact", NFILE, fnm), xlabel);
3701 /* calc average and sigma for each histogram
3702 (maybe required for bootstrapping. If not, this is fast anyhow) */
3703 if (opt.nBootStrap && opt.bsMethod == bsMethod_trajGauss)
3705 averageSigma(window, nwins);
3708 /* Get initial potential by simple integration */
3709 if (opt.bInitPotByIntegration)
3711 guessPotByIntegration(window, nwins, &opt, xlabel);
3714 /* Check if complete reaction coordinate is covered */
3715 checkReactionCoordinateCovered(window, nwins, &opt);
3717 /* Calculate profile */
3718 snew(profile, opt.bins);
3726 if ( (i%opt.stepUpdateContrib) == 0)
3728 setup_acc_wham(profile, window, nwins, &opt);
3730 if (maxchange < opt.Tolerance)
3733 /* if (opt.verbose) */
3734 printf("Switched to exact iteration in iteration %d\n", i);
3736 calc_profile(profile, window, nwins, &opt, bExact);
3737 if (((i%opt.stepchange) == 0 || i == 1) && i != 0)
3739 printf("\t%4d) Maximum change %e\n", i, maxchange);
3743 while ( (maxchange = calc_z(profile, window, nwins, &opt, bExact)) > opt.Tolerance || !bExact);
3744 printf("Converged in %d iterations. Final maximum change %g\n", i, maxchange);
3746 /* calc error from Kumar's formula */
3747 /* Unclear how the error propagates along reaction coordinate, therefore
3749 /* calc_error_kumar(profile,window, nwins,&opt); */
3751 /* Write profile in energy units? */
3754 prof_normalization_and_unit(profile, &opt);
3755 std::strcpy(ylabel, en_unit_label[opt.unit]);
3756 std::strcpy(title, "Umbrella potential");
3760 std::strcpy(ylabel, "Density of states");
3761 std::strcpy(title, "Density of states");
3764 /* symmetrize profile around z=0? */
3767 symmetrizeProfile(profile, &opt);
3770 /* write profile or density of states */
3771 profout = xvgropen(opt2fn("-o", NFILE, fnm), title, xlabel, ylabel, opt.oenv);
3772 for (i = 0; i < opt.bins; ++i)
3774 fprintf(profout, "%e\t%e\n", (i+0.5)/opt.bins*(opt.max-opt.min)+opt.min, profile[i]);
3777 printf("Wrote %s\n", opt2fn("-o", NFILE, fnm));
3779 /* Bootstrap Method */
3782 do_bootstrapping(opt2fn("-bsres", NFILE, fnm), opt2fn("-bsprof", NFILE, fnm),
3783 opt2fn("-hist", NFILE, fnm),
3784 xlabel, ylabel, profile, window, nwins, &opt);
3788 freeUmbrellaWindows(window, nfiles);
3790 printf("\nIn case you use results from gmx wham for a publication, please cite:\n");
3791 please_cite(stdout, "Hub2010");