Merge branch 'master' of git@git.gromacs.org:gromacs
[gromacs/rigid-bodies.git] / include / gmx_statistics.h
blob4640c78f9f1d0eb9e6d44e773829367b8181b2a8
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
31 * And Hey:
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 #ifndef _GMX_STATS_H
36 #define _GMX_STATS_H
38 #ifdef __cplusplus
39 extern "C" {
40 #endif
42 typedef struct gmx_stats *gmx_stats_t;
44 /* Error codes returned by the routines */
45 enum { estatsOK, estatsNO_POINTS, estatsNO_MEMORY, estatsERROR,
46 estatsINVALID_INPUT, estatsNOT_IMPLEMENTED, estatsNR };
48 enum { elsqWEIGHT_NONE, elsqWEIGHT_X, elsqWEIGHT_Y,
49 elsqWEIGHT_XY, elsqWEIGHT_NR };
51 extern gmx_stats_t gmx_stats_init();
53 extern int gmx_stats_done(gmx_stats_t stats);
55 /* Remove outliers from a straight line, where level in units of
56 sigma. Level needs to be larger than one obviously. */
57 extern int gmx_stats_remove_outliers(gmx_stats_t stats,double level);
59 extern int gmx_stats_add_point(gmx_stats_t stats,double x,double y,
60 double dx,double dy);
62 /* The arrays dx and dy may be NULL if no uncertainties are available,
63 in that case zero uncertainties will be assumed. */
64 extern int gmx_stats_add_points(gmx_stats_t stats,int n,real *x,real *y,
65 real *dx,real *dy);
67 /* Return the data points one by one. Return estatsOK while there are
68 more points, and returns estatsNOPOINTS when the last point has
69 been returned. Should be used in a while loop. Variables for either
70 pointer may be NULL, in which case the routine can be used as an
71 expensive point counter. */
72 extern int gmx_stats_get_point(gmx_stats_t stats,real *x,real *y,
73 real *dx,real *dy);
75 /* Fit the data to y = ax + b, possibly weighted, if uncertainties
76 have been input. Returns slope in *a and intercept in b, *return
77 sigmas in *da and *db respectively. Returns normalized *quality of
78 fit in *chi2 and correlation of fit with data in Rfit. chi2, Rfit,
79 da and db may be NULL. */
80 extern int gmx_stats_get_ab(gmx_stats_t stats,int weight,
81 real *a,real *b,
82 real *da,real *db,real *chi2,real *Rfit);
84 /* Fit the data to y = ax, possibly weighted, if uncertainties have
85 been input. Returns slope in *a, sigma in a in *da, and normalized
86 quality of fit in *chi2 and correlation of fit with data in
87 Rfit. chi2, Rfit and da may be NULL. */
88 extern int gmx_stats_get_a(gmx_stats_t stats,int weight,
89 real *a,real *da,real *chi2,real *Rfit);
91 /* Return the correlation coefficient between the data (x and y) as
92 input to the structure. */
93 extern int gmx_stats_get_corr_coeff(gmx_stats_t stats,real *R);
95 /* Returns the root mean square deviation between x and y values. */
96 extern int gmx_stats_get_rmsd(gmx_stats_t gstats,real *rmsd);
98 extern int gmx_stats_get_npoints(gmx_stats_t stats,int *N);
100 extern int gmx_stats_get_average(gmx_stats_t stats,real *aver);
102 extern int gmx_stats_get_sigma(gmx_stats_t stats,real *sigma);
104 extern int gmx_stats_get_error(gmx_stats_t stats,real *error);
106 /* Get all three of the above. Pointers may be null, in which case no
107 assignment will be done. */
108 extern int gmx_stats_get_ase(gmx_stats_t gstats,real *aver,real *sigma,real *error);
110 /* Dump the x, y, dx, dy data to a text file */
111 extern int gmx_stats_dump_xy(gmx_stats_t gstats,FILE *fp);
113 /* Make a histogram of the data present. Uses either bindwith to
114 determine the number of bins, or nbins to determine the binwidth,
115 therefore one of these should be zero, but not the other. If
116 normalized not equal to zero, the integral of the histogram will be
117 normalized to one. The output is in two arrays, *x and *y, to which
118 you should pass a pointer. Memory for the arrays will be allocated
119 as needed. Function returns one of the estats codes. */
120 extern int gmx_stats_make_histogram(gmx_stats_t gstats,real binwidth,int nbins,
121 int normalized,real **x,real **y);
123 /* Return message belonging to error code */
124 extern const char *gmx_stats_message(int estats);
126 /****************************************************
127 * Some statistics utilities for convenience: useful when a complete data
128 * set is available already from another source, e.g. an xvg file.
129 ****************************************************/
130 extern int lsq_y_ax(int n, real x[], real y[], real *a);
131 /* Fit a straight line y=ax thru the n data points x, y, return the
132 slope in *a. Return value can be estatsOK, or something else. */
134 extern int lsq_y_ax_b(int n, real x[], real y[], real *a, real *b,real *r,
135 real *chi2);
136 /* Fit a straight line y=ax+b thru the n data points x,y.
137 * Returns the "fit quality" sigma = sqrt(chi^2/(n-2)).
138 * The correlation coefficient is returned in r.
141 extern int lsq_y_ax_b_xdouble(int n, double x[], real y[],
142 real *a, real *b,real *r,real *chi2);
143 /* As lsq_y_ax_b, but with x in double precision.
146 extern int lsq_y_ax_b_error(int n, real x[], real y[], real dy[],
147 real *a, real *b, real *da, real *db,
148 real *r,real *chi2);
149 /* Fit a straight line y=ax+b thru the n data points x,y, with sigma dy
150 * Returns the "fit quality" sigma = sqrt(chi^2/(n-2)).
151 * The correlation coefficient is returned in r.
154 #ifdef __cplusplus
156 #endif
158 #endif