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
Go to the documentation of this file.
1/*
2 * LALInferenceKDE.h: Bayesian Followup, kernel density estimator.
3 *
4 * Copyright (C) 2013 Ben Farr
5 *
6 *
7 * This program is free software; you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation; either version 2 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with with program; see the file COPYING. If not, write to the
19 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
20 * MA 02110-1301 USA
21 */
22#ifndef LALInferenceKDE_h
23#define LALInferenceKDE_h
24
25#include <lal/LALStdlib.h>
26#include <lal/LALConstants.h>
27#include <lal/LALDatatypes.h>
28#include <lal/LALInference.h>
29
30struct tagkmeans;
31
32/**
33 * Structure containing the Guassian kernel density of a set of samples.
34 */
35typedef struct
36tagKDE
37{
38 gsl_matrix *data; /**< Data to estimate the underlying distribution of */
39 INT4 dim; /**< Dimension of points in \a data. */
40 INT4 npts; /**< Number of points in \a data. */
41 REAL8 bandwidth; /**< Bandwidth of kernels. */
42 REAL8 log_norm_factor; /**< Normalization factor of the KDE. */
43 gsl_vector * mean; /**< The mean of \a data */
44 gsl_matrix * cov; /**< The covariance matrix of \a data */
45 gsl_matrix * cholesky_decomp_cov; /**< The Cholesky decomposition of the
46 covariance matrix, containing both terms
47 as returned by gsl_linalg_cholesky_decomp(). */
48 gsl_matrix * cholesky_decomp_cov_lower; /**< Just the lower portion of \a cholesky_decomp_cov. */
49
50 LALInferenceParamVaryType * lower_bound_types; /**< Array of param boundary types */
51 LALInferenceParamVaryType * upper_bound_types; /**< Array of param boundary types */
52 REAL8 * lower_bounds; /**< Lower param bounds */
53 REAL8 * upper_bounds; /**< Upper param bounds */
55
56/* Allocate, fill, and tune a Gaussian kernel density estimate given an array of points. */
57LALInferenceKDE *LALInferenceNewKDE(REAL8 *pts, INT4 npts, INT4 dim, INT4 *mask);
58
59/* Allocate, fill, and tune a Gaussian kernel density estimate given a matrix of points. */
60LALInferenceKDE *LALInferenceNewKDEfromMat(gsl_matrix *data, INT4 *mask);
61
62/* Allocate, fill, and tune a Gaussian kernel density estimate given a matrix of points. */
64
65/* Free an allocated KDE structure. */
67
68/* Calculate the bandwidth and normalization factor for a KDE. */
70
71/* Evaluate the (log) PDF from a KDE at a single point. */
73
74/* Draw a sample from a kernel density estimate. */
76
77/* Compute the Cholesky decomposition of a matrix. */
78INT4 LALInferenceCholeskyDecompose(gsl_matrix *mat);
79
80/* Calculate the covariance matrix of a data set. */
81void LALInferenceComputeCovariance(gsl_matrix *cov, gsl_matrix *pts);
82
83/* Calculate the mean of a data set. */
84void LALInferenceComputeMean(gsl_vector *mean, gsl_matrix *points);
85
86/* Calculate the determinant of a matrix. */
87REAL8 LALInferenceMatrixDet(gsl_matrix *mat);
88
89/* Determine the log of the sum of an array of exponentials. */
91
92#endif
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)
double REAL8
int32_t INT4
LALInferenceParamVaryType
An enumerated type for denoting the topology of a parameter.
Definition: LALInference.h:127
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.