Prototypes | |
LALInferenceKDE * | LALInferenceNewKDE (REAL8 *pts, INT4 npts, INT4 dim, INT4 *mask) |
Allocate, fill, and tune a Gaussian kernel density estimate from an array of samples. More... | |
LALInferenceKDE * | LALInferenceNewKDEfromMat (gsl_matrix *data, INT4 *mask) |
Allocate, fill, and tune a Gaussian kernel density estimate from a matrix of points. More... | |
LALInferenceKDE * | LALInferenceInitKDE (INT4 npts, INT4 dim) |
Construct an empty KDE structure. More... | |
void | LALInferenceDestroyKDE (LALInferenceKDE *kde) |
Free an allocated KDE structure. More... | |
void | LALInferenceSetKDEBandwidth (LALInferenceKDE *kde) |
Calculate the bandwidth and normalization factor for a KDE. More... | |
REAL8 | LALInferenceKDEEvaluatePoint (LALInferenceKDE *kde, REAL8 *point) |
Evaluate the (log) PDF from a KDE at a single point. More... | |
REAL8 * | LALInferenceDrawKDESample (LALInferenceKDE *kde, gsl_rng *rng) |
Draw a sample from a kernel density estimate. More... | |
INT4 | LALInferenceCholeskyDecompose (gsl_matrix *mat) |
Compute the Cholesky decomposition of a matrix. More... | |
void | LALInferenceComputeCovariance (gsl_matrix *cov, gsl_matrix *pts) |
Calculate the covariance matrix of a data set. More... | |
void | LALInferenceComputeMean (gsl_vector *mean, gsl_matrix *points) |
Calculate the mean of a data set. More... | |
REAL8 | LALInferenceMatrixDet (gsl_matrix *mat) |
Calculate the determinant of a matrix. More... | |
REAL8 | log_add_exps (REAL8 *vals, INT4 size) |
Determine the log of the sum of an array of exponentials. More... | |
Go to the source code of this file.
Data Structures | |
struct | LALInferenceKDE |
Structure containing the Guassian kernel density of a set of samples. More... | |
LALInferenceKDE * LALInferenceNewKDE | ( | REAL8 * | pts, |
INT4 | npts, | ||
INT4 | dim, | ||
INT4 * | mask | ||
) |
Allocate, fill, and tune a Gaussian kernel density estimate from an array of samples.
Given a set of points, the distribution from which the points were sampled from is estimated by a kernel density estimate (KDE), using a Gaussian kernal. This routine initializes a new KDE and sets the bandwidth according to Scott's rule. From this it is possible to both evaluate the estimate of the probability density function at an arbitrary location, and generate more samples from the estimated distribution. A mask can be provided to select a subset of the provided data, otherwise all data is used.
[in] | pts | Array containing the points to estimate the distribution of. |
[in] | npts | The total number of points contained in pts. |
[in] | dim | The number of dimensions of the points contained in pts. |
[in] | mask | An optional npts long 0/1 array to mask points in pts. |
Definition at line 63 of file LALInferenceKDE.c.
LALInferenceKDE * LALInferenceNewKDEfromMat | ( | gsl_matrix * | data, |
INT4 * | mask | ||
) |
Allocate, fill, and tune a Gaussian kernel density estimate from a matrix of points.
Estimate the underlying distribution of samples in a GSL matrix using a Gaussian KDE.
[in] | data | GSL matrix with rows of samples from a target distribution. |
[in] | mask | An optional npts long 0/1 array to mask points in pts. |
Definition at line 109 of file LALInferenceKDE.c.
LALInferenceKDE * LALInferenceInitKDE | ( | INT4 | npts, |
INT4 | dim | ||
) |
Construct an empty KDE structure.
Create an empty LALInferenceKDE structure, allocated to handle the given size and dimension.
[in] | npts | Number of samples that will be used to estimate the distribution. |
[in] | dim | Number of dimensions that the probability density function will be estimated in. |
Definition at line 155 of file LALInferenceKDE.c.
void LALInferenceDestroyKDE | ( | LALInferenceKDE * | kde | ) |
Free an allocated KDE structure.
Frees all memory allocated for a given KDE structure.
[in] | kde | The KDE structure to be freed. |
Definition at line 189 of file LALInferenceKDE.c.
void LALInferenceSetKDEBandwidth | ( | LALInferenceKDE * | kde | ) |
Calculate the bandwidth and normalization factor for a KDE.
Use Scott's rule to determine the bandwidth, and corresponding normalization factor, for a KDE.
[in] | kde | The kernel density estimate to estimate the bandwidth of. |
Definition at line 247 of file LALInferenceKDE.c.
REAL8 LALInferenceKDEEvaluatePoint | ( | LALInferenceKDE * | kde, |
REAL8 * | point | ||
) |
Evaluate the (log) PDF from a KDE at a single point.
Calculate the (log) value of the probability density function estimate from a kernel density estimate at a single point.
[in] | kde | The kernel density estimate to evaluate. |
[in] | point | An array containing the point to evaluate the PDF at. |
Definition at line 300 of file LALInferenceKDE.c.
REAL8 * LALInferenceDrawKDESample | ( | LALInferenceKDE * | kde, |
gsl_rng * | rng | ||
) |
Draw a sample from a kernel density estimate.
Draw a sample from a distribution, as estimated by a kernel density estimator.
[in] | kde | The kernel density estimate to draw point from. |
[in] | rng | GSL random number generator to use. |
Definition at line 429 of file LALInferenceKDE.c.
INT4 LALInferenceCholeskyDecompose | ( | gsl_matrix * | mat | ) |
Compute the Cholesky decomposition of a matrix.
A wrapper for gsl_linalg_cholesky_decomp() that avoids halting if the matrix is found not to be positive-definite. This often happens when decomposing a covariance matrix that is poorly estimated due to low sample size.
mat | The matrix to decompose (in place). |
Definition at line 217 of file LALInferenceKDE.c.
void LALInferenceComputeCovariance | ( | gsl_matrix * | cov, |
gsl_matrix * | data | ||
) |
Calculate the covariance matrix of a data set.
Calculate the covariance matrix of a data set contained in the rows of a GSL matrix.
[out] | cov | GSL matrix containing the covariance matrix of data. |
[in] | data | GSL matrix data set to compute the covariance matrix of. |
Definition at line 554 of file LALInferenceKDE.c.
void LALInferenceComputeMean | ( | gsl_vector * | mean, |
gsl_matrix * | data | ||
) |
Calculate the mean of a data set.
Calculate the mean vector of a data set contained in the rows of a GSL matrix.
[out] | mean | GSL vector containing the mean. |
[in] | data | GSL matrix data set to compute the mean of. |
Definition at line 528 of file LALInferenceKDE.c.
REAL8 LALInferenceMatrixDet | ( | gsl_matrix * | mat | ) |
Calculate the determinant of a matrix.
Calculates the determinant of a GSL matrix.
[in] | mat | The matrix to calculate the determinant of. |
Definition at line 500 of file LALInferenceKDE.c.
Determine the log of the sum of an array of exponentials.
Utility for calculating the log of the sum of an array of exponentials. Useful for avoiding overflows.
[in] | vals | Array of values to be exponentiated, summed, and logged. |
[in] | size | Number of elements in vals. |
Definition at line 601 of file LALInferenceKDE.c.