22#ifndef LALInferenceClusteredKDE_h
23#define LALInferenceClusteredKDE_h
25#include <lal/LALStdlib.h>
26#include <lal/LALConstants.h>
27#include <lal/LALDatatypes.h>
29#include <lal/LALInferenceKDE.h>
REAL8 LALInferenceKmeansBIC(LALInferenceKmeans *kmeans)
Calculate the BIC using the KDEs of the cluster distributions.
REAL8 LALInferenceKmeansMaxLogL(LALInferenceKmeans *kmeans)
Calculate max likelihood of a kmeans assuming spherical Gaussian clusters.
void LALInferenceKmeansRecursiveSplit(LALInferenceKmeans *kmeans, INT4 ntrials, gsl_rng *rng)
Recursively split a k=1 kmeans.
REAL8 euclidean_dist_squared(gsl_vector *x, gsl_vector *y)
Calculate the squared Euclidean distance bewteen two points.
void LALInferenceKmeansImposeBounds(LALInferenceKmeans *kmeans, LALInferenceVariables *params, LALInferenceVariables *priorArgs, INT4 cyclic_reflective)
Impose boundaries on individual KDEs.
LALInferenceKmeans * LALInferenceCreateKmeans(INT4 k, gsl_matrix *data, gsl_rng *rng)
Generate a new kmeans struct from a set of data.
void LALInferenceKmeansForgyInitialize(LALInferenceKmeans *kmeans)
Randomly select points from the data as initial centroids for a kmeans run.
LALInferenceKmeans * LALInferenceRecursiveKmeans(gsl_matrix *data, INT4 ntrials, gsl_rng *rng)
Run kmeans, recursively splitting centroids, accepting if the BIC increases.
void accumulate_vectors(gsl_matrix **mat_ptr, gsl_vector *vec)
Add a vector to the end of a matrix, allocating if necessary.
gsl_matrix * mask_data(gsl_matrix *data, INT4 *mask)
Generate a new matrix by masking an existing one.
void LALInferenceKmeansBuildKDE(LALInferenceKmeans *kmeans)
Build the kernel density estimate from a kmeans clustering.
LALInferenceKmeans * LALInferenceXmeans(gsl_matrix *data, INT4 ntrials, gsl_rng *rng)
Cluster data xmeans-style, splitting centroids according to the BIC.
REAL8 LALInferenceKmeansPDF(LALInferenceKmeans *kmeans, REAL8 *pt)
Evaluate the estimated (log) PDF from a clustered-KDE at a point.
LALInferenceKmeans * LALInferenceIncrementalKmeans(gsl_matrix *data, INT4 ntrials, gsl_rng *rng)
Kmeans cluster data, increasing k until the BIC stops increasing.
gsl_matrix * LALInferenceWhitenSamples(gsl_matrix *samples)
Transform a data set to obtain a 0-mean and identity covariance matrix.
LALInferenceKmeans * LALInferenceOptimizedKmeans(gsl_matrix *data, INT4 ntrials, gsl_rng *rng)
Kmeans cluster data, maximizing BIC assuming BIC(k) is concave-down.
void LALInferenceKmeansRun(LALInferenceKmeans *kmeans)
Run the kmeans algorithm until cluster assignments don't change.
REAL8 LALInferenceWhitenedKmeansPDF(LALInferenceKmeans *kmeans, REAL8 *pt)
Estimate (log) PDF from a clustered-KDE at an already whitened point.
void euclidean_centroid(gsl_vector *centroid, gsl_matrix *data, INT4 *mask)
Find the centroid of a masked data set.
void append_vec_to_mat(gsl_matrix **mat_ptr, gsl_vector *vec)
Add a vector to the end of a matrix.
void LALInferenceKmeansAssignment(LALInferenceKmeans *kmeans)
The assignment step of the kmeans algorithm.
void LALInferenceKmeansDestroy(LALInferenceKmeans *kmeans)
Destroy a kmeans instance.
void LALInferenceKmeansUpdate(LALInferenceKmeans *kmeans)
The update step of the kmeans algorithm.
LALInferenceKmeans * LALInferenceKmeansExtractCluster(LALInferenceKmeans *kmeans, INT4 cluster_id)
Extract a single cluster from an existing kmeans as a new 1-means.
void LALInferenceKmeansConstructMask(LALInferenceKmeans *kmeans, INT4 *mask, INT4 cluster_id)
Construct a mask to select only the data assigned to a single cluster.
void LALInferenceKmeansSeededInitialize(LALInferenceKmeans *kmeans)
A 'k-means++'-like seeded initialization of centroids for a kmeans run.
LALInferenceKmeans * LALInferenceKmeansRunBestOf(INT4 k, gsl_matrix *samples, INT4 ntrials, gsl_rng *rng)
Run a kmeans several times and return the best.
REAL8 * LALInferenceKmeansDraw(LALInferenceKmeans *kmeans)
Draw a sample from a kmeans-KDE estimate of a distribution.
void LALInferenceKmeansRandomPartitionInitialize(LALInferenceKmeans *kmeans)
Use the resulting centroids from a random cluster assignment as the initial centroids of a kmeans run...
Structure containing the Guassian kernel density of a set of samples.
Structure for performing the kmeans clustering algorithm on a set of samples.
INT4 npts
Number of points being clustered.
INT4 * sizes
Cluster sizes.
REAL8 error
Error of current clustering.
gsl_vector * std
Std deviations of unwhitened data.
INT4 * assignments
Cluster assignments.
INT4 k
Number of clusters.
gsl_matrix * recursive_centroids
Matrix used to accumulate tested centroids.
gsl_matrix * data
Data to be clustered (typically whitened)
REAL8 * weights
Fraction of data points in each cluster.
gsl_rng * rng
Random number generator.
gsl_matrix * cov
Covariance matrix of data.
gsl_matrix * centroids
Array with rows containing cluster centroids.
LALInferenceKDE ** KDEs
Array of KDEs, one for each cluster.
INT4 has_changed
Flag indicating a change to cluster assignmens.
INT4 dim
Dimension of data.
gsl_vector * mean
Mean of unwhitened data (for tranforming points)
INT4 * mask
Mask used to select data from individual clusters.
void(* centroid)(gsl_vector *centroid, gsl_matrix *data, INT4 *mask)
Find centroid.
The LALInferenceVariables structure to contain a set of parameters Implemented as a linked list of LA...