New entry
[gromacs/rigid-bodies.git] / include / gmx_statistics.h
blob5d28891b98faa401154d8f801dd8a0eb4480e17c
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 #include "typedefs.h"
44 typedef struct gmx_stats *gmx_stats_t;
46 /* Error codes returned by the routines */
47 enum { estatsOK, estatsNO_POINTS, estatsNO_MEMORY, estatsERROR,
48 estatsINVALID_INPUT, estatsNOT_IMPLEMENTED, estatsNR };
50 enum { elsqWEIGHT_NONE, elsqWEIGHT_X, elsqWEIGHT_Y,
51 elsqWEIGHT_XY, elsqWEIGHT_NR };
53 enum { ehistoX, ehistoY, ehistoNR };
55 gmx_stats_t gmx_stats_init();
57 int gmx_stats_done(gmx_stats_t stats);
59 /* Remove outliers from a straight line, where level in units of
60 sigma. Level needs to be larger than one obviously. */
61 int gmx_stats_remove_outliers(gmx_stats_t stats,double level);
63 int gmx_stats_add_point(gmx_stats_t stats,double x,double y,
64 double dx,double dy);
66 /* The arrays dx and dy may be NULL if no uncertainties are available,
67 in that case zero uncertainties will be assumed. */
68 int gmx_stats_add_points(gmx_stats_t stats,int n,real *x,real *y,
69 real *dx,real *dy);
71 /* Return the data points one by one. Return estatsOK while there are
72 more points, and returns estatsNOPOINTS when the last point has
73 been returned. Should be used in a while loop. Variables for either
74 pointer may be NULL, in which case the routine can be used as an
75 expensive point counter. */
76 int gmx_stats_get_point(gmx_stats_t stats,real *x,real *y,
77 real *dx,real *dy);
79 /* Fit the data to y = ax + b, possibly weighted, if uncertainties
80 have been input. Returns slope in *a and intercept in b, *return
81 sigmas in *da and *db respectively. Returns normalized *quality of
82 fit in *chi2 and correlation of fit with data in Rfit. chi2, Rfit,
83 da and db may be NULL. */
84 int gmx_stats_get_ab(gmx_stats_t stats,int weight,
85 real *a,real *b,
86 real *da,real *db,real *chi2,real *Rfit);
88 /* Fit the data to y = ax, possibly weighted, if uncertainties have
89 been input. Returns slope in *a, sigma in a in *da, and normalized
90 quality of fit in *chi2 and correlation of fit with data in
91 Rfit. chi2, Rfit and da may be NULL. */
92 int gmx_stats_get_a(gmx_stats_t stats,int weight,
93 real *a,real *da,real *chi2,real *Rfit);
95 /* Return the correlation coefficient between the data (x and y) as
96 input to the structure. */
97 int gmx_stats_get_corr_coeff(gmx_stats_t stats,real *R);
99 /* Returns the root mean square deviation between x and y values. */
100 int gmx_stats_get_rmsd(gmx_stats_t gstats,real *rmsd);
102 int gmx_stats_get_npoints(gmx_stats_t stats,int *N);
104 int gmx_stats_get_average(gmx_stats_t stats,real *aver);
106 int gmx_stats_get_sigma(gmx_stats_t stats,real *sigma);
108 int gmx_stats_get_error(gmx_stats_t stats,real *error);
110 /* Get all three of the above. Pointers may be null, in which case no
111 assignment will be done. */
112 int gmx_stats_get_ase(gmx_stats_t gstats,real *aver,real *sigma,real *error);
114 /* Dump the x, y, dx, dy data to a text file */
115 int gmx_stats_dump_xy(gmx_stats_t gstats,FILE *fp);
117 /* Make a histogram of the data present. Uses either bindwith to
118 determine the number of bins, or nbins to determine the binwidth,
119 therefore one of these should be zero, but not the other. If *nbins = 0
120 the number of bins will be returned in this variable. ehisto should be one of
121 ehistoX or ehistoY. If
122 normalized not equal to zero, the integral of the histogram will be
123 normalized to one. The output is in two arrays, *x and *y, to which
124 you should pass a pointer. Memory for the arrays will be allocated
125 as needed. Function returns one of the estats codes. */
126 int gmx_stats_make_histogram(gmx_stats_t gstats,real binwidth,int *nbins,
127 int ehisto,
128 int normalized,real **x,real **y);
130 /* Return message belonging to error code */
131 const char *gmx_stats_message(int estats);
133 /****************************************************
134 * Some statistics utilities for convenience: useful when a complete data
135 * set is available already from another source, e.g. an xvg file.
136 ****************************************************/
137 int lsq_y_ax(int n, real x[], real y[], real *a);
138 /* Fit a straight line y=ax thru the n data points x, y, return the
139 slope in *a. Return value can be estatsOK, or something else. */
141 int lsq_y_ax_b(int n, real x[], real y[], real *a, real *b,real *r,
142 real *chi2);
143 /* Fit a straight line y=ax+b thru the n data points x,y.
144 * Returns the "fit quality" sigma = sqrt(chi^2/(n-2)).
145 * The correlation coefficient is returned in r.
148 int lsq_y_ax_b_xdouble(int n, double x[], real y[],
149 real *a, real *b,real *r,real *chi2);
150 /* As lsq_y_ax_b, but with x in double precision.
153 int lsq_y_ax_b_error(int n, real x[], real y[], real dy[],
154 real *a, real *b, real *da, real *db,
155 real *r,real *chi2);
156 /* Fit a straight line y=ax+b thru the n data points x,y, with sigma dy
157 * Returns the "fit quality" sigma = sqrt(chi^2/(n-2)).
158 * The correlation coefficient is returned in r.
161 #ifdef __cplusplus
163 #endif
165 #endif