1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*- */
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Green Red Orange Magenta Azure Cyan Skyblue
44 #include "gmx_statistics.h"
46 static int gmx_dnint(double x)
51 typedef struct gmx_stats {
52 double aa,a,b,sigma_aa,sigma_a,sigma_b,aver,sigma_aver,error;
53 double rmsd,Rdata,Rfit,Rfitaa,chi2,chi2aa;
59 gmx_stats_t gmx_stats_init()
65 return (gmx_stats_t) stats;
68 int gmx_stats_get_npoints(gmx_stats_t gstats, int *N)
70 gmx_stats *stats = (gmx_stats *) gstats;
77 int gmx_stats_done(gmx_stats_t gstats)
79 gmx_stats *stats = (gmx_stats *) gstats;
93 int gmx_stats_add_point(gmx_stats_t gstats,double x,double y,
96 gmx_stats *stats = (gmx_stats *) gstats;
99 if (stats->np+1 >= stats->nalloc)
101 if (stats->nalloc == 0)
102 stats->nalloc = 1024;
105 srenew(stats->x,stats->nalloc);
106 srenew(stats->y,stats->nalloc);
107 srenew(stats->dx,stats->nalloc);
108 srenew(stats->dy,stats->nalloc);
109 for(i=stats->np; (i<stats->nalloc); i++)
117 stats->x[stats->np] = x;
118 stats->y[stats->np] = y;
119 stats->dx[stats->np] = dx;
120 stats->dy[stats->np] = dy;
127 int gmx_stats_get_point(gmx_stats_t gstats,real *x,real *y,
130 gmx_stats *stats = (gmx_stats *) gstats;
132 if (stats->np_c < stats->np)
134 if (NULL != x) *x = stats->x[stats->np_c];
135 if (NULL != y) *y = stats->y[stats->np_c];
136 if (NULL != dx) *dx = stats->dx[stats->np_c];
137 if (NULL != dy) *dy = stats->dy[stats->np_c];
144 return estatsNO_POINTS;
147 int gmx_stats_add_points(gmx_stats_t gstats,int n,real *x,real *y,
154 if ((ok = gmx_stats_add_point(gstats,x[i],y[i],
155 (NULL != dx) ? dx[i] : 0,
156 (NULL != dy) ? dy[i] : 0)) != estatsOK)
164 static int gmx_stats_compute(gmx_stats *stats,int weight)
166 double yy,yx,xx,sx,sy,dy,chi2,chi2aa,d2;
167 double ssxx,ssyy,ssxy;
168 double w,wtot,yx_nw,sy_nw,sx_nw,yy_nw,xx_nw,dx2,dy2;
172 if (stats->computed == 0)
176 return estatsNO_POINTS;
188 d2 += dsqr(stats->x[i]-stats->y[i]);
189 if ((stats->dy[i]) && (weight == elsqWEIGHT_Y))
191 w = 1/dsqr(stats->dy[i]);
200 xx += w*dsqr(stats->x[i]);
201 xx_nw += dsqr(stats->x[i]);
203 yy += w*dsqr(stats->y[i]);
204 yy_nw += dsqr(stats->y[i]);
206 yx += w*stats->y[i]*stats->x[i];
207 yx_nw += stats->y[i]*stats->x[i];
210 sx_nw += stats->x[i];
213 sy_nw += stats->y[i];
216 /* Compute average, sigma and error */
217 stats->aver = sy_nw/N;
218 stats->sigma_aver = sqrt(yy_nw/N - dsqr(sy_nw/N));
219 stats->error = stats->sigma_aver/sqrt(N);
221 /* Compute RMSD between x and y */
222 stats->rmsd = sqrt(d2/N);
224 /* Correlation coefficient for data */
230 ssxx = N*(xx_nw - dsqr(sx_nw));
231 ssyy = N*(yy_nw - dsqr(sy_nw));
232 ssxy = N*(yx_nw - (sx_nw*sy_nw));
233 stats->Rdata = sqrt(dsqr(ssxy)/(ssxx*ssyy));
235 /* Compute straight line through datapoints, either with intercept
236 zero (result in aa) or with intercept variable (results in a
244 stats->a = (yx-sx*sy)/(xx-sx*sx);
245 stats->b = (sy)-(stats->a)*(sx);
247 /* Compute chi2, deviation from a line y = ax+b. Also compute
248 chi2aa which returns the deviation from a line y = ax. */
253 if (stats->dy[i] > 0)
261 chi2aa += dsqr((stats->y[i]-(stats->aa*stats->x[i]))/dy);
262 chi2 += dsqr((stats->y[i]-(stats->a*stats->x[i]+stats->b))/dy);
266 stats->chi2 = sqrt(chi2/(N-2));
267 stats->chi2aa = sqrt(chi2aa/(N-2));
269 /* Look up equations! */
272 stats->sigma_a = sqrt(stats->chi2/((N-2)*dx2));
273 stats->sigma_b = stats->sigma_a*sqrt(xx);
274 stats->Rfit = fabs(ssxy)/sqrt(ssxx*ssyy);
275 /*stats->a*sqrt(dx2/dy2);*/
276 stats->Rfitaa = stats->aa*sqrt(dx2/dy2);
294 int gmx_stats_get_ab(gmx_stats_t gstats,int weight,
295 real *a,real *b,real *da,real *db,
296 real *chi2,real *Rfit)
298 gmx_stats *stats = (gmx_stats *) gstats;
301 if ((ok = gmx_stats_compute(stats,weight)) != estatsOK)
313 *da = stats->sigma_a;
317 *db = stats->sigma_b;
331 int gmx_stats_get_a(gmx_stats_t gstats,int weight,real *a,real *da,
332 real *chi2,real *Rfit)
334 gmx_stats *stats = (gmx_stats *) gstats;
337 if ((ok = gmx_stats_compute(stats,weight)) != estatsOK)
339 if (NULL != a) *a = stats->aa;
340 if (NULL != da) *da = stats->sigma_aa;
341 if (NULL != chi2) *chi2 = stats->chi2aa;
342 if (NULL != Rfit) *Rfit = stats->Rfitaa;
347 int gmx_stats_get_average(gmx_stats_t gstats,real *aver)
349 gmx_stats *stats = (gmx_stats *) gstats;
352 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
362 int gmx_stats_get_ase(gmx_stats_t gstats,real *aver,real *sigma,real *error)
364 gmx_stats *stats = (gmx_stats *) gstats;
367 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
378 *sigma = stats->sigma_aver;
382 *error = stats->error;
388 int gmx_stats_get_sigma(gmx_stats_t gstats,real *sigma)
390 gmx_stats *stats = (gmx_stats *) gstats;
393 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
396 *sigma = stats->sigma_aver;
401 int gmx_stats_get_error(gmx_stats_t gstats,real *error)
403 gmx_stats *stats = (gmx_stats *) gstats;
406 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
409 *error = stats->error;
414 int gmx_stats_get_corr_coeff(gmx_stats_t gstats,real *R)
416 gmx_stats *stats = (gmx_stats *) gstats;
419 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
427 int gmx_stats_get_rmsd(gmx_stats_t gstats,real *rmsd)
429 gmx_stats *stats = (gmx_stats *) gstats;
432 if ((ok = gmx_stats_compute(stats,elsqWEIGHT_NONE)) != estatsOK)
442 int gmx_stats_dump_xy(gmx_stats_t gstats,FILE *fp)
444 gmx_stats *stats = (gmx_stats *) gstats;
447 for(i=0; (i<stats->np); i++)
449 fprintf(fp,"%12g %12g %12g %12g\n",stats->x[i],stats->y[i],
450 stats->dx[i],stats->dy[i]);
456 int gmx_stats_remove_outliers(gmx_stats_t gstats,double level)
458 gmx_stats *stats = (gmx_stats *) gstats;
459 int i,iter=1,done=0,ok;
462 while ((stats->np >= 10) && !done)
464 if ((ok = gmx_stats_get_rmsd(gstats,&rmsd)) != estatsOK)
469 for(i=0; (i<stats->np); )
471 r = fabs(stats->x[i]-stats->y[i]);
474 fprintf(stderr,"Removing outlier, iter = %d, rmsd = %g, x = %g, y = %g\n",
475 iter,rmsd,stats->x[i],stats->y[i]);
478 stats->x[i] = stats->x[stats->np-1];
479 stats->y[i] = stats->y[stats->np-1];
480 stats->dx[i] = stats->dx[stats->np-1];
481 stats->dy[i] = stats->dy[stats->np-1];
497 int gmx_stats_make_histogram(gmx_stats_t gstats,real binwidth,int *nb,
498 int ehisto,int normalized,real **x,real **y)
500 gmx_stats *stats = (gmx_stats *) gstats;
501 int i,ok,index=0,nbins=*nb,*nindex;
502 double minx,maxx,maxy,miny,delta,dd,minh;
504 if (((binwidth <= 0) && (nbins <= 0)) ||
505 ((binwidth > 0) && (nbins > 0)))
507 return estatsINVALID_INPUT;
511 return estatsNO_POINTS;
513 minx = maxx = stats->x[0];
514 miny = maxy = stats->y[0];
515 for(i=1; (i<stats->np); i++)
517 miny = (stats->y[i] < miny) ? stats->y[i] : miny;
518 maxy = (stats->y[i] > maxy) ? stats->y[i] : maxy;
519 minx = (stats->x[i] < minx) ? stats->x[i] : minx;
520 maxx = (stats->x[i] > maxx) ? stats->x[i] : maxx;
522 if (ehisto == ehistoX)
527 else if (ehisto == ehistoY)
533 return estatsINVALID_INPUT;
537 binwidth = (delta)/nbins;
541 nbins = gmx_dnint((delta)/binwidth + 0.5);
545 for(i=0; (i<nbins); i++)
547 (*x)[i] = minh + binwidth*(i+0.5);
555 dd = 1.0/(binwidth*stats->np);
559 for(i=0; (i<stats->np); i++)
561 if (ehisto == ehistoY)
562 index = (stats->y[i]-miny)/binwidth;
563 else if (ehisto == ehistoX)
564 index = (stats->x[i]-minx)/binwidth;
578 for(i=0; (i<nbins); i++)
580 (*y)[i] /= nindex[i];
587 static const char *stats_error[estatsNR] =
589 "All well in STATS land",
592 "Invalid histogram input",
594 "Not implemented yet"
597 const char *gmx_stats_message(int estats)
599 if ((estats >= 0) && (estats < estatsNR))
601 return stats_error[estats];
605 return stats_error[estatsERROR];
609 /* Old convenience functions, should be merged with the core
611 int lsq_y_ax(int n, real x[], real y[], real *a)
613 gmx_stats_t lsq = gmx_stats_init();
617 gmx_stats_add_points(lsq,n,x,y,0,0);
618 if ((ok = gmx_stats_get_a(lsq,elsqWEIGHT_NONE,a,&da,&chi2,&Rfit)) != estatsOK)
627 for (i=0; i<n; i++) {
636 static int low_lsq_y_ax_b(int n, real *xr, double *xd, real yr[],
637 real *a, real *b,real *r,real *chi2)
642 lsq = gmx_stats_init();
645 if ((ok = gmx_stats_add_point(lsq,(NULL != xd) ? xd[i] : xr[i],yr[i],0,0))
651 if ((ok = gmx_stats_get_ab(lsq,elsqWEIGHT_NONE,a,b,NULL,NULL,chi2,r)) != estatsOK)
658 double x,y,yx,xx,yy,sx,sy,chi2;
661 for (i=0; i<n; i++) {
675 *a = (n*yx-sy*sx)/(n*xx-sx*sx);
677 *r = sqrt((xx-sx*sx)/(yy-sy*sy));
682 chi2 += dsqr(yr[i] - ((*a)*xd[i] + (*b)));
685 chi2 += dsqr(yr[i] - ((*a)*xr[i] + (*b)));
689 return sqrt(chi2/(n-2));
695 int lsq_y_ax_b(int n, real x[], real y[], real *a, real *b,real *r,real *chi2)
697 return low_lsq_y_ax_b(n,x,NULL,y,a,b,r,chi2);
700 int lsq_y_ax_b_xdouble(int n, double x[], real y[], real *a, real *b,
703 return low_lsq_y_ax_b(n,NULL,x,y,a,b,r,chi2);
706 int lsq_y_ax_b_error(int n, real x[], real y[], real dy[],
707 real *a, real *b, real *da, real *db,
713 lsq = gmx_stats_init();
716 if ((ok = gmx_stats_add_point(lsq,x[i],y[i],0,dy[i])) != estatsOK)
721 if ((ok = gmx_stats_get_ab(lsq,elsqWEIGHT_Y,a,b,da,db,chi2,r)) != estatsOK)
725 if ((ok = gmx_stats_done(lsq)) != estatsOK)
733 double sxy,sxx,syy,sx,sy,w,s_2,dx2,dy2,mins;
735 sxy=sxx=syy=sx=sy=w=0.0;
738 mins = min(mins,dy[i]);
740 gmx_fatal(FARGS,"Zero or negative weigths in linear regression analysis");
742 for (i=0; i<n; i++) {
743 s_2 = dsqr(1.0/dy[i]);
744 sxx += s_2*dsqr(x[i]);
745 sxy += s_2*y[i]*x[i];
746 syy += s_2*dsqr(y[i]);
763 *chi2+=dsqr((y[i]-((*a)*x[i]+(*b)))/dy[i]);
766 *da = sqrt(*chi2/((n-2)*dx2));
768 *r = *a*sqrt(dx2/dy2);
771 fprintf(debug,"sx = %g, sy = %g, sxy = %g, sxx = %g, w = %g\n"
772 "chi2 = %g, dx2 = %g\n",
773 sx,sy,sxy,sxx,w,*chi2,dx2);
776 *chi2 = sqrt(*chi2/(n-2));