MixtureofGaussianHMM Class Reference

A wrapper class for HMM functionals in Mixture of Gaussion case. More...

Collaboration diagram for MixtureofGaussianHMM:
[legend]

List of all members.

Public Member Functions

void ComputeLogLikelihood (const ArrayList< Matrix > &list_data_seq, ArrayList< double > *list_likelihood) const
 Compute the log-likelihood of a list of sequences.
double ComputeLogLikelihood (const Matrix &data_seq) const
 Compute the log-likelihood of a sequence.
void ComputeViterbiStateSequence (const Matrix &data_seq, Vector *state_seq) const
 Compute the most probable sequence (Viterbi).
void DecodeInit (const Matrix &data_seq, Matrix *state_prob_mat, Matrix *forward_prob_mat, Matrix *backward_prob_mat, Vector *scale_vec) const
 A decode version that initialized the output matrices.
void DecodeOverwrite (const Matrix &data_seq, Matrix *state_prob_mat, Matrix *forward_prob_mat, Matrix *backward_prob_mat, Vector *scale_vec) const
 Decode a sequence into probabilities of each state at each time step using scaled forward-backward algorithm.
void EstimateModel (int numstate, int numcluster, const Matrix &data_seq, const Vector &state_seq)
void EstimateModel (int numcluster, const Matrix &data_seq, const Vector &state_seq)
 Estimate the matrices by a data sequence and a state sequence Must be already initialized.
void GenerateSequence (int L, Matrix *data_seq, Vector *state_seq) const
 Generate a random data sequence of a given length.
void Init ()
 Initializes empty object.
void Init (const Matrix &transmission, const ArrayList< MixtureGauss > &list_mixture_gauss)
 Initializes from computed transmission and Mixture of Gaussian parameters.
void InitFromFile (const char *profile)
 Initializes by loading from a file.
const ArrayList< MixtureGauss > & list_mixture_gauss () const
void LoadProfile (const char *profile)
 Load from file, used when already initialized.
void SaveProfile (const char *profile) const
 Save matrices to file.
void setModel (const Matrix &transmission, const ArrayList< MixtureGauss > &list_mixture_gauss)
 Setter used when already initialized.
void TrainBaumWelch (const ArrayList< Matrix > &list_data_seq, int max_iteration, double tolerance)
 Train the model with a list of sequences, must be already initialized using Baum-Welch EM algorithm.
void TrainViterbi (const ArrayList< Matrix > &list_data_seq, int max_iteration, double tolerance)
 Train the model with a list of sequences, must be already initialized using Viterbi algorithm to determine the state sequence of each sequence.
const Matrix & transmission () const
 Getters.

Static Public Member Functions

static void BackwardProcedure (int L, const Matrix &trans, const Matrix &emis_prob, const Vector &scales, Matrix *bs)
static void CalculateEmissionProb (const Matrix &seq, const ArrayList< MixtureGauss > &mixs, Matrix *emis_prob)
static double Decode (int L, const Matrix &trans, const Matrix &emis_prob, Matrix *pstates, Matrix *fs, Matrix *bs, Vector *scales)
static double Decode (const Matrix &trans, const Matrix &emis_prob, Matrix *pstates, Matrix *fs, Matrix *bs, Vector *scales)
static void EstimateInit (int numStates, int NumClusters, const Matrix &seq, const Vector &states, Matrix *trans, ArrayList< MixtureGauss > *mixs)
static void EstimateInit (int NumClusters, const Matrix &seq, const Vector &states, Matrix *trans, ArrayList< MixtureGauss > *mixs)
 Estimate transition and emission distribution from sequence and states.
static void ForwardProcedure (int L, const Matrix &trans, const Matrix &emis_prob, Vector *scales, Matrix *fs)
 Calculate posteriori probabilities of states at each steps Scaled Forward - Backward procedure trans: Transition probabilities, size M x M emis_prob: Emission probabilities along the sequence, size M x L (L is the sequence length) pstates: size M x L fs: scaled forward probabities, size M x L bs: scaled backward probabities, size M x L scales: scale factors, length L RETURN: log probabilities of sequence.
static void GenerateInit (int L, const Matrix &trans, const ArrayList< MixtureGauss > &mixs, Matrix *seq, Vector *states)
 Generating a sequence and states using transition and emission probabilities.
static success_t LoadProfile (const char *profile, Matrix *trans, ArrayList< MixtureGauss > *mixs)
static success_t SaveProfile (const char *profile, const Matrix &trans, const ArrayList< MixtureGauss > &mixs)
static void Train (const ArrayList< Matrix > &seqs, Matrix *guessTR, ArrayList< MixtureGauss > *guessMG, int max_iter, double tol)
 Baum-Welch and Viterbi estimation of transition and emission distribution (Gaussian).
static void TrainViterbi (const ArrayList< Matrix > &seqs, Matrix *guessTR, ArrayList< MixtureGauss > *guessMG, int max_iter, double tol)
static double ViterbiInit (int L, const Matrix &trans, const Matrix &emis_prob, Vector *states)
static double ViterbiInit (const Matrix &trans, const Matrix &emis_prob, Vector *states)
 Calculate the most probable states for a sequence Viterbi algorithm trans: Transition probabilities, size M x M emis_prob: Emission probabilities, size M x L states: Unitialized, will have length L RETURN: log probability of the most probable sequence.

Detailed Description

A wrapper class for HMM functionals in Mixture of Gaussion case.

This class maintains transition probabilities and Mixture of Gaussian parameters (mean and covariance) for each state, more details below.

Definition at line 53 of file mixgaussHMM.h.


Member Function Documentation

void MixtureofGaussianHMM::ComputeLogLikelihood ( const ArrayList< Matrix > &  list_data_seq,
ArrayList< double > *  list_likelihood 
) const

Compute the log-likelihood of a list of sequences.

Definition at line 133 of file mixgaussHMM.cc.

References ForwardProcedure(), ArrayList< TElem >::Init(), GenVector< T >::Init(), ArrayList< TElem >::PushBackCopy(), and ArrayList< TElem >::size().

double MixtureofGaussianHMM::ComputeLogLikelihood ( const Matrix &  data_seq  )  const

Compute the log-likelihood of a sequence.

Definition at line 118 of file mixgaussHMM.cc.

References ForwardProcedure(), and GenVector< T >::Init().

void MixtureofGaussianHMM::ComputeViterbiStateSequence ( const Matrix &  data_seq,
Vector state_seq 
) const

Compute the most probable sequence (Viterbi).

Definition at line 153 of file mixgaussHMM.cc.

References ViterbiInit().

void MixtureofGaussianHMM::DecodeInit ( const Matrix &  data_seq,
Matrix *  state_prob_mat,
Matrix *  forward_prob_mat,
Matrix *  backward_prob_mat,
Vector scale_vec 
) const

A decode version that initialized the output matrices.

Definition at line 102 of file mixgaussHMM.cc.

References GenVector< T >::Init().

void MixtureofGaussianHMM::DecodeOverwrite ( const Matrix &  data_seq,
Matrix *  state_prob_mat,
Matrix *  forward_prob_mat,
Matrix *  backward_prob_mat,
Vector scale_vec 
) const

Decode a sequence into probabilities of each state at each time step using scaled forward-backward algorithm.

Also return forward, backward probabilities and scale factors

Definition at line 90 of file mixgaussHMM.cc.

void MixtureofGaussianHMM::EstimateInit ( int  NumClusters,
const Matrix &  seq,
const Vector states,
Matrix *  trans,
ArrayList< MixtureGauss > *  mixs 
) [static]

Estimate transition and emission distribution from sequence and states.

Definition at line 320 of file mixgaussHMM.cc.

References GenVector< T >::length().

Referenced by EstimateModel().

void MixtureofGaussianHMM::EstimateModel ( int  numcluster,
const Matrix &  data_seq,
const Vector state_seq 
)

Estimate the matrices by a data sequence and a state sequence Must be already initialized.

Definition at line 74 of file mixgaussHMM.cc.

References EstimateInit().

void MixtureofGaussianHMM::ForwardProcedure ( int  L,
const Matrix &  trans,
const Matrix &  emis_prob,
Vector scales,
Matrix *  fs 
) [static]

Calculate posteriori probabilities of states at each steps Scaled Forward - Backward procedure trans: Transition probabilities, size M x M emis_prob: Emission probabilities along the sequence, size M x L (L is the sequence length) pstates: size M x L fs: scaled forward probabities, size M x L bs: scaled backward probabities, size M x L scales: scale factors, length L RETURN: log probabilities of sequence.

Definition at line 329 of file mixgaussHMM.cc.

Referenced by ComputeLogLikelihood().

void MixtureofGaussianHMM::GenerateInit ( int  L,
const Matrix &  trans,
const ArrayList< MixtureGauss > &  mixs,
Matrix *  seq,
Vector states 
) [static]

Generating a sequence and states using transition and emission probabilities.

L: sequence length trans: Matrix M x M (M states) means: list of mean vectors length N (emission vector of length N) covs: list of square root of covariance matrices size N x N seq: generated sequence, uninitialized matrix, will have size N x L states: generated states, uninitialized vector, will have length L

Definition at line 216 of file mixgaussHMM.cc.

References GenVector< T >::Init(), and ArrayList< TElem >::size().

Referenced by GenerateSequence().

void MixtureofGaussianHMM::GenerateSequence ( int  L,
Matrix *  data_seq,
Vector state_seq 
) const

Generate a random data sequence of a given length.

Definition at line 70 of file mixgaussHMM.cc.

References GenerateInit().

void MixtureofGaussianHMM::Init (  )  [inline]

Initializes empty object.

Definition at line 82 of file mixgaussHMM.h.

References ArrayList< TElem >::Init().

void MixtureofGaussianHMM::Init ( const Matrix &  transmission,
const ArrayList< MixtureGauss > &  list_mixture_gauss 
)

Initializes from computed transmission and Mixture of Gaussian parameters.

void MixtureofGaussianHMM::InitFromFile ( const char *  profile  ) 

Initializes by loading from a file.

Definition at line 55 of file mixgaussHMM.cc.

References LoadProfile().

Referenced by LoadProfile().

void MixtureofGaussianHMM::LoadProfile ( const char *  profile  ) 

Load from file, used when already initialized.

Definition at line 60 of file mixgaussHMM.cc.

References InitFromFile().

Referenced by InitFromFile().

void MixtureofGaussianHMM::SaveProfile ( const char *  profile  )  const

Save matrices to file.

Definition at line 66 of file mixgaussHMM.cc.

void MixtureofGaussianHMM::setModel ( const Matrix &  transmission,
const ArrayList< MixtureGauss > &  list_mixture_gauss 
)

Setter used when already initialized.

Definition at line 45 of file mixgaussHMM.cc.

References ArrayList< TElem >::InitCopy(), and ArrayList< TElem >::size().

void MixtureofGaussianHMM::Train ( const ArrayList< Matrix > &  seqs,
Matrix *  guessTR,
ArrayList< MixtureGauss > *  guessMG,
int  max_iter,
double  tol 
) [static]

Baum-Welch and Viterbi estimation of transition and emission distribution (Gaussian).

Definition at line 368 of file mixgaussHMM.cc.

References GenVector< T >::get(), ArrayList< TElem >::Init(), GenVector< T >::Init(), GenMatrix< T >::Init(), ArrayList< TElem >::PushBackCopy(), GenVector< T >::SetZero(), and ArrayList< TElem >::size().

Referenced by TrainBaumWelch().

void MixtureofGaussianHMM::TrainBaumWelch ( const ArrayList< Matrix > &  list_data_seq,
int  max_iteration,
double  tolerance 
)

Train the model with a list of sequences, must be already initialized using Baum-Welch EM algorithm.

Definition at line 161 of file mixgaussHMM.cc.

References Train().

void MixtureofGaussianHMM::TrainViterbi ( const ArrayList< Matrix > &  list_data_seq,
int  max_iteration,
double  tolerance 
)

Train the model with a list of sequences, must be already initialized using Viterbi algorithm to determine the state sequence of each sequence.

Definition at line 165 of file mixgaussHMM.cc.

const Matrix& MixtureofGaussianHMM::transmission (  )  const [inline]

Getters.

Definition at line 68 of file mixgaussHMM.h.

double MixtureofGaussianHMM::ViterbiInit ( const Matrix &  trans,
const Matrix &  emis_prob,
Vector states 
) [static]

Calculate the most probable states for a sequence Viterbi algorithm trans: Transition probabilities, size M x M emis_prob: Emission probabilities, size M x L states: Unitialized, will have length L RETURN: log probability of the most probable sequence.

Definition at line 349 of file mixgaussHMM.cc.

Referenced by ComputeViterbiStateSequence().


The documentation for this class was generated from the following files:
Generated on Mon Jan 24 12:04:40 2011 for FASTlib by  doxygen 1.6.3