22#ifndef LALInferenceKDE_h
23#define LALInferenceKDE_h
25#include <lal/LALStdlib.h>
26#include <lal/LALConstants.h>
27#include <lal/LALDatatypes.h>
28#include <lal/LALInference.h>
INT4 LALInferenceCholeskyDecompose(gsl_matrix *mat)
Compute the Cholesky decomposition of a matrix.
void LALInferenceSetKDEBandwidth(LALInferenceKDE *kde)
Calculate the bandwidth and normalization factor for a KDE.
LALInferenceKDE * LALInferenceNewKDE(REAL8 *pts, INT4 npts, INT4 dim, INT4 *mask)
Allocate, fill, and tune a Gaussian kernel density estimate from an array of samples.
REAL8 * LALInferenceDrawKDESample(LALInferenceKDE *kde, gsl_rng *rng)
Draw a sample from a kernel density estimate.
REAL8 LALInferenceKDEEvaluatePoint(LALInferenceKDE *kde, REAL8 *point)
Evaluate the (log) PDF from a KDE at a single point.
void LALInferenceDestroyKDE(LALInferenceKDE *kde)
Free an allocated KDE structure.
REAL8 log_add_exps(REAL8 *vals, INT4 size)
Determine the log of the sum of an array of exponentials.
void LALInferenceComputeMean(gsl_vector *mean, gsl_matrix *points)
Calculate the mean of a data set.
LALInferenceKDE * LALInferenceInitKDE(INT4 npts, INT4 dim)
Construct an empty KDE structure.
REAL8 LALInferenceMatrixDet(gsl_matrix *mat)
Calculate the determinant of a matrix.
void LALInferenceComputeCovariance(gsl_matrix *cov, gsl_matrix *pts)
Calculate the covariance matrix of a data set.
LALInferenceKDE * LALInferenceNewKDEfromMat(gsl_matrix *data, INT4 *mask)
Allocate, fill, and tune a Gaussian kernel density estimate from a matrix of points.
static REAL8 mean(REAL8 *array, int N)
LALInferenceParamVaryType
An enumerated type for denoting the topology of a parameter.
Structure containing the Guassian kernel density of a set of samples.
REAL8 log_norm_factor
Normalization factor of the KDE.
REAL8 * upper_bounds
Upper param bounds.
LALInferenceParamVaryType * lower_bound_types
Array of param boundary types.
INT4 npts
Number of points in data.
gsl_matrix * cov
The covariance matrix of data.
LALInferenceParamVaryType * upper_bound_types
Array of param boundary types.
gsl_matrix * cholesky_decomp_cov
The Cholesky decomposition of the covariance matrix, containing both terms as returned by gsl_linalg_...
gsl_matrix * data
Data to estimate the underlying distribution of.
REAL8 bandwidth
Bandwidth of kernels.
INT4 dim
Dimension of points in data.
REAL8 * lower_bounds
Lower param bounds.
gsl_vector * mean
The mean of data.
gsl_matrix * cholesky_decomp_cov_lower
Just the lower portion of cholesky_decomp_cov.