Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Typedefs | Enumerations | Functions
statistics.h File Reference
#include <cstdio>
#include <tuple>
#include "gromacs/utility/real.h"
+ Include dependency graph for statistics.h:
+ This graph shows which files directly or indirectly include this file:

Description

Declares simple statistics toolbox.

Authors
David van der Spoel david.nosp@m..van.nosp@m.dersp.nosp@m.oel@.nosp@m.icm.u.nosp@m.u.se

Typedefs

typedef struct gmx_stats * gmx_stats_t
 Abstract container type.
 

Enumerations

enum  {
  elsqWEIGHT_NONE, elsqWEIGHT_X, elsqWEIGHT_Y, elsqWEIGHT_XY,
  elsqWEIGHT_NR
}
 Enum for statistical weights.
 

Functions

gmx_stats_t gmx_stats_init ()
 Initiate a data structure. More...
 
void gmx_stats_free (gmx_stats_t stats)
 Destroy a data structure. More...
 
void gmx_stats_add_point (gmx_stats_t stats, double x, double y, double dx, double dy)
 Add a point to the data set. More...
 
void gmx_stats_get_ab (gmx_stats_t stats, int weight, real *a, real *b, real *da, real *db, real *chi2, real *Rfit)
 Fit the data to y = ax + b, possibly weighted, if uncertainties have been input. da and db may be NULL. More...
 
real gmx_stats_get_average (gmx_stats_t stats)
 Computes and returns the average value. More...
 
std::tuple< real, real, realgmx_stats_get_ase (gmx_stats_t stats)
 Pointers may be null, in which case no assignment will be done. More...
 
void lsq_y_ax_b (int n, real x[], real y[], real *a, real *b, real *r, real *chi2)
 Fit a straight line y=ax+b thru the n data points x, y. More...
 
void lsq_y_ax_b_xdouble (int n, double x[], real y[], real *a, real *b, real *r, real *chi2)
 Fit a straight line y=ax+b thru the n data points x, y. More...
 
void lsq_y_ax_b_error (int n, real x[], real y[], real dy[], real *a, real *b, real *da, real *db, real *r, real *chi2)
 Fit a straight line y=ax+b thru the n data points x, y. More...
 

Function Documentation

void gmx_stats_add_point ( gmx_stats_t  stats,
double  x,
double  y,
double  dx,
double  dy 
)

Add a point to the data set.

Parameters
[in]statsThe data structure
[in]xThe x value
[in]yThe y value
[in]dxThe error in the x value
[in]dyThe error in the y value
void gmx_stats_free ( gmx_stats_t  stats)

Destroy a data structure.

Parameters
statsThe data structure
void gmx_stats_get_ab ( gmx_stats_t  stats,
int  weight,
real a,
real b,
real da,
real db,
real chi2,
real Rfit 
)

Fit the data to y = ax + b, possibly weighted, if uncertainties have been input. da and db may be NULL.

Parameters
[in]statsThe data structure
[in]weighttype of weighting
[out]aslope
[out]bintercept
[out]dasigma in a
[out]dbsigma in b
[out]chi2normalized quality of fit
[out]Rfitcorrelation coefficient
std::tuple<real, real, real> gmx_stats_get_ase ( gmx_stats_t  stats)

Pointers may be null, in which case no assignment will be done.

Parameters
[in]statsThe data structure
Returns
Tuple of (average value, its standard deviation, its standard error)
Exceptions
InconsistentInputErrorif given no points to analyze
real gmx_stats_get_average ( gmx_stats_t  stats)

Computes and returns the average value.

Parameters
[in]statsThe data structure
Returns
Average value
Exceptions
InconsistentInputErrorif given no points to average
gmx_stats_t gmx_stats_init ( )

Initiate a data structure.

Returns
the data structure
void lsq_y_ax_b ( int  n,
real  x[],
real  y[],
real a,
real b,
real r,
real chi2 
)

Fit a straight line y=ax+b thru the n data points x, y.

Parameters
[in]nnumber of points
[in]xdata points x
[in]ydata point y
[out]aslope
[out]bintercept
[out]rcorrelation coefficient
[out]chi2quality of fit
Exceptions
InconsistentInputErrorif given no points to fit
void lsq_y_ax_b_error ( int  n,
real  x[],
real  y[],
real  dy[],
real a,
real b,
real da,
real db,
real r,
real chi2 
)

Fit a straight line y=ax+b thru the n data points x, y.

Parameters
[in]nnumber of points
[in]xdata points x
[in]ydata point y
[in]dyuncertainty in data point y
[out]aslope
[out]bintercept
[out]daerror in slope
[out]dberror in intercept
[out]rcorrelation coefficient
[out]chi2quality of fit
Exceptions
InconsistentInputErrorif given no points to fit
void lsq_y_ax_b_xdouble ( int  n,
double  x[],
real  y[],
real a,
real b,
real r,
real chi2 
)

Fit a straight line y=ax+b thru the n data points x, y.

Parameters
[in]nnumber of points
[in]xdata points x
[in]ydata point y
[out]aslope
[out]bintercept
[out]rcorrelation coefficient
[out]chi2quality of fit
Exceptions
InconsistentInputErrorif given no points to fit Suits cases where x is already always computed in double precision even in a mixed-precision build configuration.