Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALInference 4.1.9.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALInferenceKDE.h File Reference

Prototypes

LALInferenceKDELALInferenceNewKDE (REAL8 *pts, INT4 npts, INT4 dim, INT4 *mask)
 Allocate, fill, and tune a Gaussian kernel density estimate from an array of samples. More...
 
LALInferenceKDELALInferenceNewKDEfromMat (gsl_matrix *data, INT4 *mask)
 Allocate, fill, and tune a Gaussian kernel density estimate from a matrix of points. More...
 
LALInferenceKDELALInferenceInitKDE (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...
 
REAL8LALInferenceDrawKDESample (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...
 

Function Documentation

◆ LALInferenceNewKDE()

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.

Parameters
[in]ptsArray containing the points to estimate the distribution of.
[in]nptsThe total number of points contained in pts.
[in]dimThe number of dimensions of the points contained in pts.
[in]maskAn optional npts long 0/1 array to mask points in pts.
Returns
A LALInferenceKDE structure containing the KDE of the distribution.

Definition at line 63 of file LALInferenceKDE.c.

◆ LALInferenceNewKDEfromMat()

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.

Parameters
[in]dataGSL matrix with rows of samples from a target distribution.
[in]maskAn optional npts long 0/1 array to mask points in pts.
Returns
A LALInferenceKDE structure containing the KDE of the distribution.
See also
LALInferenceNewKDE()

Definition at line 109 of file LALInferenceKDE.c.

◆ LALInferenceInitKDE()

LALInferenceKDE * LALInferenceInitKDE ( INT4  npts,
INT4  dim 
)

Construct an empty KDE structure.

Create an empty LALInferenceKDE structure, allocated to handle the given size and dimension.

Parameters
[in]nptsNumber of samples that will be used to estimate the distribution.
[in]dimNumber of dimensions that the probability density function will be estimated in.
Returns
An allocated, empty LALInferenceKDE structure.
See also
LALInferenceKDE, LALInferenceSetKDEBandwidth()

Definition at line 155 of file LALInferenceKDE.c.

◆ LALInferenceDestroyKDE()

void LALInferenceDestroyKDE ( LALInferenceKDE kde)

Free an allocated KDE structure.

Frees all memory allocated for a given KDE structure.

Parameters
[in]kdeThe KDE structure to be freed.
See also
LALInferenceKDE, LALInferenceInitKDE

Definition at line 189 of file LALInferenceKDE.c.

◆ LALInferenceSetKDEBandwidth()

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.

Parameters
[in]kdeThe kernel density estimate to estimate the bandwidth of.

Definition at line 247 of file LALInferenceKDE.c.

◆ LALInferenceKDEEvaluatePoint()

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.

Parameters
[in]kdeThe kernel density estimate to evaluate.
[in]pointAn array containing the point to evaluate the PDF at.
Returns
The value of the estimated probability density function at point.

Definition at line 300 of file LALInferenceKDE.c.

◆ LALInferenceDrawKDESample()

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.

Parameters
[in]kdeThe kernel density estimate to draw point from.
[in]rngGSL random number generator to use.
Returns
The sample drawn from the distribution.

Definition at line 429 of file LALInferenceKDE.c.

◆ LALInferenceCholeskyDecompose()

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.

Parameters
matThe matrix to decompose (in place).
Returns
Status of call to gsl_linalg_cholesky_decomp().

Definition at line 217 of file LALInferenceKDE.c.

◆ LALInferenceComputeCovariance()

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.

Parameters
[out]covGSL matrix containing the covariance matrix of data.
[in]dataGSL matrix data set to compute the covariance matrix of.

Definition at line 554 of file LALInferenceKDE.c.

◆ LALInferenceComputeMean()

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.

Parameters
[out]meanGSL vector containing the mean.
[in]dataGSL matrix data set to compute the mean of.

Definition at line 528 of file LALInferenceKDE.c.

◆ LALInferenceMatrixDet()

REAL8 LALInferenceMatrixDet ( gsl_matrix *  mat)

Calculate the determinant of a matrix.

Calculates the determinant of a GSL matrix.

Parameters
[in]matThe matrix to calculate the determinant of.
Returns
The determinant of mat.

Definition at line 500 of file LALInferenceKDE.c.

◆ log_add_exps()

REAL8 log_add_exps ( REAL8 vals,
INT4  size 
)

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.

Parameters
[in]valsArray of values to be exponentiated, summed, and logged.
[in]sizeNumber of elements in vals.
Returns
The log of the sum of elements in vals.

Definition at line 601 of file LALInferenceKDE.c.