gmm
Functions
-
int main(int argc, char *argv[])
gmm [ option ] [ infile ]
-l int
length of vector \((1 \le L)\)
-m int
order of vector \((0 \le L - 1)\)
-k int
number of mixtures \((1 \le K)\)
-i int
number of iterations \((1 \le N)\)
-d double
convergence threshold \((0 \le \epsilon)\)
-w double
floor value of weight \((0 \le F_w \le 1/K)\)
-v double
floor value of variance \((0 \le F_v)\)
-M double
MAP smoothing parameter \((0 \le \alpha \le 1)\)
-U str
double-type initial GMM parameters
-S str
double-type total log-likelihood
-f
use full covariance
-V
show average log likelihood at each iteration
-B int+
block size of covariance matrix
infile str
double-type training data sequencea
stdout
double-type GMM parameters
The following examples show four types of covariance.
gmm -l 10 < data.d > diag.gmm gmm -l 10 -f < data.d > full.gmm gmm -l 10 -B 5 5 < data.d > block-wise-diag.gmm gmm -l 10 -f -B 5 5 < data.d > block-diag.gmm
If -M option is specified, the MAP estimates of the GMM paramaeters are obtained using universal background model.
gmm -k 8 < data1.d > ubm.gmm gmm -k 8 -U ubm.gmm -M 0.1 < data2.d > map.gmm
- Parameters:
argc – [in] Number of arguments.
argv – [in] Argument vector.
- Returns:
0 on success, 1 on failure.
-
class GaussianMixtureModeling
Estimate model parameters of GMM.
The input is the \(M\)-th order input vectors:
\[ \begin{array}{cccc} \boldsymbol{x}_0, & \boldsymbol{x}_1, & \ldots, & \boldsymbol{x}_{T-1}, \end{array} \]where \(T\) is the number of vectors. The output is the set of GMM parameters \(\boldsymbol{\lambda}\):\[ \begin{array}{ccccccc} w_0, & \boldsymbol{\mu}_0, & \boldsymbol{\varSigma}_0, & \ldots, & w_{K-1}, & \boldsymbol{\mu}_{K-1}, & \boldsymbol{\varSigma}_{K-1} \end{array} \]where \(K\) is the number of mixture components, \(w_k\) is the \(k\)-th mixture weight, \(\mu_k\) is the \(k\)-th mean vector, and \(\boldsymbol{\varSigma}_k\) is the \(k\)-th covariance matrix. The mixture weights satisfy\[ \sum_{k=0}^{K-1} w_k = 1. \]The covariance matrix can be full, diagonal, block-diagonal, or block-wise diagonal.The \(\boldsymbol{\lambda}\) is iteratively updated by the following update formulae:
\[\begin{split}\begin{eqnarray} \hat{w}_k &=& \frac{1}{T} \sum_{t=0}^{T-1} \gamma_{k,t}, \\ \boldsymbol{\mu}_k &=& \frac{\displaystyle\sum_{t=0}^{T-1} \gamma_{k,t} \boldsymbol{x}_t} {\displaystyle\sum_{t=0}^{T-1} \gamma_{k,t}}, \\ \boldsymbol{\varSigma}_k &=& \frac{\displaystyle\sum_{t=0}^{T-1} \gamma_{k,t} \boldsymbol{x}_t \boldsymbol{x}_t^{\mathsf{T}}} {\displaystyle\sum_{t=0}^{T-1} \gamma_{k,t}} - \boldsymbol{\mu}_t \boldsymbol{\mu}_t^{\mathsf{T}}, \end{eqnarray}\end{split}\]where\[ \gamma_{k,t} = \frac{w_k \mathcal{N}( \boldsymbol{x}_t \,|\, \boldsymbol{\mu}_k, \boldsymbol{\varSigma}_k)} {\displaystyle\sum_{k=0}^{K-1} w_k \mathcal{N}( \boldsymbol{x}_t \,|\, \boldsymbol{\mu}_k, \boldsymbol{\varSigma}_k)} \]is the posterior probability of being in the \(k\)-th component at time \(t\).If the universal background model (UBM) \(\boldsymbol{\lambda}'\) is given, the \(\boldsymbol{\lambda}\) is estimated by the maximum a posteriori method. The joint prior density is the product of Dirichlet and normal-Wishart densities. The update formulae are obtained as
\[\begin{split}\begin{eqnarray} \hat{w}_k &=& \frac{\xi_k + \displaystyle\sum_{t=0}^{T-1} \gamma_{k,t}} {\displaystyle\sum_{k=0}^{K-1} \xi_k + T}, \\ \boldsymbol{\mu}_k &=& \frac{\xi_k \boldsymbol{\mu}'_k + \displaystyle\sum_{t=0}^{T-1} \gamma_{k,t} \boldsymbol{x}_t} {\xi_k + \displaystyle\sum_{t=0}^{T-1} \gamma_{k,t}}, \\ \boldsymbol{\varSigma}_k &=& \frac{\xi_k \boldsymbol{\varSigma}'_k + \xi_k (\boldsymbol{\mu}'_k - \boldsymbol{\mu}_k) (\boldsymbol{\mu}'_k - \boldsymbol{\mu}_k)^{\mathsf{T}} + \displaystyle\sum_{t=0}^{T-1} \gamma_{k,t} (\boldsymbol{x}_t - \boldsymbol{\mu}_k) (\boldsymbol{x}_t - \boldsymbol{\mu}_k)^{\mathsf{T}}} {\xi_k + \displaystyle\sum_{t=0}^{T-1} \gamma_{k,t}}. \end{eqnarray}\end{split}\]where\[ \xi_k = \alpha w'_k. \]and \(\alpha\) controlls the importance of the UBM.Public Types
Public Functions
-
GaussianMixtureModeling(int num_order, int num_mixture, int num_iteration, double convergence_threshold, CovarianceType covariance_type, std::vector<int> block_size, double weight_floor, double variance_floor, InitializationType initialization_type, int log_interval, double smoothing_parameter = 0.0, const std::vector<double> &ubm_weights = {}, const std::vector<std::vector<double>> &ubm_mean_vectors = {}, const std::vector<SymmetricMatrix> &ubm_covariance_matrices = {})
- Parameters:
num_order – [in] Order of vector, \(M\).
num_mixture – [in] Number of mixtures, \(K\).
num_iteration – [in] Number of iterations, \(N\).
convergence_threshold – [in] Convergence threshold.
covariance_type – [in] Type of covariance.
block_size – [in] Block size of covariance.
weight_floor – [in] Floor value of weight.
variance_floor – [in] Floor value of variance.
initialization_type – [in] Type of initialization.
log_interval – [in] Show log-likelihood every this step.
smoothing_parameter – [in] MAP hyperparameter (optional), \(\alpha\).
ubm_weights – [in] Weights of UBM-GMM (optional).
ubm_mean_vectors – [in] Means of UBM-GMM (optional).
ubm_covariance_matrices – [in] Covariances of UBM-GMM (optional).
-
inline int GetNumOrder() const
- Returns:
Order of vector.
-
inline int GetNumMixture() const
- Returns:
Number of mixture components.
-
inline int GetNumIteration() const
- Returns:
Number of iterations.
-
inline double GetConvergenceThreshold() const
- Returns:
Convergence threshold.
-
inline CovarianceType GetCovarianceType() const
- Returns:
Type of covariance.
-
inline double GetWeightFloor() const
- Returns:
Floor value of weight.
-
inline double GetVarianceFloor() const
- Returns:
Floor value of variance.
-
inline InitializationType GetInitializationType() const
- Returns:
Type of initialization.
-
inline double GetSmoothingParameter() const
- Returns:
MAP smoothing parameter.
-
inline bool IsDiagonal() const
- Returns:
True if covariance is pure diagonal.
-
inline bool IsValid() const
- Returns:
True if this object is valid.
-
bool Run(const std::vector<std::vector<double>> &input_vectors, std::vector<double> *weights, std::vector<std::vector<double>> *mean_vectors, std::vector<SymmetricMatrix> *covariance_matrices, double *total_log_likelihood) const
- Parameters:
input_vectors – [in] \(M\)-th order input vectors. The shape is \([T, M+1]\).
weights – [inout] \(K\) mixture weights.
mean_vectors – [inout] \(K\) mean vectors. The shape is \([K, M+1]\).
covariance_matrices – [inout] \(K\) covariance matrices. The shape is \([K, M+1, M+1]\).
total_log_likelihood – [out] Total log-likelihood.
- Returns:
True on success, false on failure.
Public Static Functions
-
static bool CalculateLogProbability(int num_order, int num_mixture, bool is_diagonal, bool check_size, const std::vector<double> &input_vector, const std::vector<double> &weights, const std::vector<std::vector<double>> &mean_vectors, const std::vector<SymmetricMatrix> &covariance_matrices, std::vector<double> *components_of_log_probability, double *log_probability, GaussianMixtureModeling::Buffer *buffer)
Calculate log-probablity of data.
- Parameters:
num_order – [in] Order of input vector.
num_mixture – [in] Number of mixture components.
is_diagonal – [in] If true, diagonal covariance is assumed.
check_size – [in] If true, check sanity of input GMM parameters.
input_vector – [in] \(M\)-th order input vector.
weights – [in] \(K\) mixture weights.
mean_vectors – [in] \(K\) mean vectors.
covariance_matrices – [in] \(K\) covariance matrices.
components_of_log_probability – [out] Components of log-probability.
log_probability – [out] Log-probability of input vector.
buffer – [out] Buffer.
- Returns:
True on success, false on failure.
-
class Buffer
-
GaussianMixtureModeling(int num_order, int num_mixture, int num_iteration, double convergence_threshold, CovarianceType covariance_type, std::vector<int> block_size, double weight_floor, double variance_floor, InitializationType initialization_type, int log_interval, double smoothing_parameter = 0.0, const std::vector<double> &ubm_weights = {}, const std::vector<std::vector<double>> &ubm_mean_vectors = {}, const std::vector<SymmetricMatrix> &ubm_covariance_matrices = {})