HMMRegression Class Reference

A class that represents a Hidden Markov Model Regression (HMMR). More...

Inheritance diagram for HMMRegression:

Public Member Functions

 HMMRegression (const size_t states, const distribution::RegressionDistribution emissions, const double tolerance=1e-5)
 Create the Hidden Markov Model Regression with the given number of hidden states and the given default regression emission. More...

 
 HMMRegression (const arma::vec &initial, const arma::mat &transition, const std::vector< distribution::RegressionDistribution > &emission, const double tolerance=1e-5)
 Create the Hidden Markov Model Regression with the given initial probability vector, the given transition matrix, and the given regression emission distributions. More...

 
double Estimate (const arma::mat &predictors, const arma::vec &responses, arma::mat &stateProb, arma::mat &forwardProb, arma::mat &backwardProb, arma::vec &scales) const
 Estimate the probabilities of each hidden state at each time step for each given data observation, using the Forward-Backward algorithm. More...

 
double Estimate (const arma::mat &predictors, const arma::vec &responses, arma::mat &stateProb) const
 Estimate the probabilities of each hidden state at each time step of each given data observation, using the Forward-Backward algorithm. More...

 
void Filter (const arma::mat &predictors, const arma::vec &responses, arma::vec &filterSeq, size_t ahead=0) const
 HMMR filtering. More...

 
double LogLikelihood (const arma::mat &predictors, const arma::vec &responses) const
 Compute the log-likelihood of the given predictors and responses. More...

 
double Predict (const arma::mat &predictors, const arma::vec &responses, arma::Row< size_t > &stateSeq) const
 Compute the most probable hidden state sequence for the given predictors and responses, using the Viterbi algorithm, returning the log-likelihood of the most likely state sequence. More...

 
void Smooth (const arma::mat &predictors, const arma::vec &responses, arma::vec &smoothSeq) const
 HMM smoothing. More...

 
void Train (const std::vector< arma::mat > &predictors, const std::vector< arma::vec > &responses)
 Train the model using the Baum-Welch algorithm, with only the given predictors and responses. More...

 
void Train (const std::vector< arma::mat > &predictors, const std::vector< arma::vec > &responses, const std::vector< arma::Row< size_t > > &stateSeq)
 Train the model using the given labeled observations; the transition and regression emissions are directly estimated. More...

 
- Public Member Functions inherited from HMM< distribution::RegressionDistribution >
 HMM (const size_t states=0, const distribution::RegressionDistribution emissions=distribution::RegressionDistribution(), const double tolerance=1e-5)
 Create the Hidden Markov Model with the given number of hidden states and the given default distribution for emissions. More...

 
 HMM (const arma::vec &initial, const arma::mat &transition, const std::vector< distribution::RegressionDistribution > &emission, const double tolerance=1e-5)
 Create the Hidden Markov Model with the given initial probability vector, the given transition matrix, and the given emission distributions. More...

 
size_t Dimensionality () const
 Get the dimensionality of observations. More...

 
size_t & Dimensionality ()
 Set the dimensionality of observations. More...

 
const std::vector< distribution::RegressionDistribution > & Emission () const
 Return the emission distributions. More...

 
std::vector< distribution::RegressionDistribution > & Emission ()
 Return a modifiable emission probability matrix reference. More...

 
double EmissionLogLikelihood (const arma::vec &emissionLogProb, double &logLikelihood, arma::vec &forwardLogProb) const
 Compute the log-likelihood of the given emission probability up to time t, storing the result in logLikelihood. More...

 
double EmissionLogScaleFactor (const arma::vec &emissionLogProb, arma::vec &forwardLogProb) const
 Compute the log of the scaling factor of the given emission probability at time t. More...

 
double Estimate (const arma::mat &dataSeq, arma::mat &stateProb, arma::mat &forwardProb, arma::mat &backwardProb, arma::vec &scales) const
 Estimate the probabilities of each hidden state at each time step for each given data observation, using the Forward-Backward algorithm. More...

 
double Estimate (const arma::mat &dataSeq, arma::mat &stateProb) const
 Estimate the probabilities of each hidden state at each time step of each given data observation, using the Forward-Backward algorithm. More...

 
void Filter (const arma::mat &dataSeq, arma::mat &filterSeq, size_t ahead=0) const
 HMM filtering. More...

 
void Generate (const size_t length, arma::mat &dataSequence, arma::Row< size_t > &stateSequence, const size_t startState=0) const
 Generate a random data sequence of the given length. More...

 
const arma::vec & Initial () const
 Return the vector of initial state probabilities. More...

 
arma::vec & Initial ()
 Modify the vector of initial state probabilities. More...

 
void load (Archive &ar, const uint32_t version)
 Load the object. More...

 
double LogEstimate (const arma::mat &dataSeq, arma::mat &stateLogProb, arma::mat &forwardLogProb, arma::mat &backwardLogProb, arma::vec &logScales) const
 Estimate the probabilities of each hidden state at each time step for each given data observation, using the Forward-Backward algorithm. More...

 
double LogLikelihood (const arma::mat &dataSeq) const
 Compute the log-likelihood of the given data sequence. More...

 
double LogLikelihood (const arma::vec &data, double &logLikelihood, arma::vec &forwardLogProb) const
 Compute the log-likelihood of the given data up to time t, storing the result in logLikelihood. More...

 
double LogScaleFactor (const arma::vec &data, arma::vec &forwardLogProb) const
 Compute the log of the scaling factor of the given data at time t. More...

 
double Predict (const arma::mat &dataSeq, arma::Row< size_t > &stateSeq) const
 Compute the most probable hidden state sequence for the given data sequence, using the Viterbi algorithm, returning the log-likelihood of the most likely state sequence. More...

 
void save (Archive &ar, const uint32_t version) const
 Save the object. More...

 
void Smooth (const arma::mat &dataSeq, arma::mat &smoothSeq) const
 HMM smoothing. More...

 
double Tolerance () const
 Get the tolerance of the Baum-Welch algorithm. More...

 
double & Tolerance ()
 Modify the tolerance of the Baum-Welch algorithm. More...

 
double Train (const std::vector< arma::mat > &dataSeq)
 Train the model using the Baum-Welch algorithm, with only the given unlabeled observations. More...

 
void Train (const std::vector< arma::mat > &dataSeq, const std::vector< arma::Row< size_t > > &stateSeq)
 Train the model using the given labeled observations; the transition and emission matrices are directly estimated. More...

 
const arma::mat & Transition () const
 Return the transition matrix. More...

 
arma::mat & Transition ()
 Return a modifiable transition matrix reference. More...

 

Additional Inherited Members

- Protected Member Functions inherited from HMM< distribution::RegressionDistribution >
void Backward (const arma::mat &dataSeq, const arma::vec &logScales, arma::mat &backwardLogProb, arma::mat &logProbs) const
 The Backward algorithm (part of the Forward-Backward algorithm). More...

 
void Forward (const arma::mat &dataSeq, arma::vec &logScales, arma::mat &forwardLogProb, arma::mat &logProbs) const
 The Forward algorithm (part of the Forward-Backward algorithm). More...

 
arma::vec ForwardAtT0 (const arma::vec &emissionLogProb, double &logScales) const
 Given emission probabilities, computes forward probabilities at time t=0. More...

 
arma::vec ForwardAtTn (const arma::vec &emissionLogProb, double &logScales, const arma::vec &prevForwardLogProb) const
 Given emission probabilities, computes forward probabilities for time t>0. More...

 
- Protected Attributes inherited from HMM< distribution::RegressionDistribution >
std::vector< distribution::RegressionDistributionemission
 Set of emission probability distributions; one for each state. More...

 
arma::mat logTransition
 Transition probability matrix. No need to be mutable in mlpack 4.0. More...

 
arma::mat transitionProxy
 A proxy variable in linear space for logTransition. More...

 

Detailed Description

A class that represents a Hidden Markov Model Regression (HMMR).

HMMR is an extension of Hidden Markov Models to regression analysis. The method is described in (Fridman, 1993) https://www.ima.umn.edu/preprints/January1994/1195.pdf An HMMR is a linear regression model whose coefficients are determined by a finite-state Markov chain. The error terms are conditionally independently normally distributed with zero mean and state-dependent variance. Let Q_t be a finite-state Markov chain, X_t a vector of predictors and Y_t a response. The HMMR is $ Y_t = X_t \beta_{Q_t} + \sigma_{Q_t} \epsilon_t $

This HMMR class supports training (supervised and unsupervised), prediction of state sequences via the Viterbi algorithm, estimation of state probabilities, filtering and smoothing of responses, and calculation of the log-likelihood of a given sequence.

Usage of the HMMR class generally involves either training an HMMR or loading an already-known HMMR and using to filter a sequence. Example code for supervised training of an HMMR is given below.

// Each column is a vector of predictors for a single observation.
arma::mat predictors(5, 100, arma::fill::randn);
// Responses for each observation
arma::vec responses(100, arma::fill::randn);
// Create an untrained HMMR with 3 hidden states
RegressionDistribution rd(predictors, responses);
arma::mat transition("0.5 0.5;" "0.5 0.5;");
std::vector<RegressionDistribution> emissions(2,rd);
HMMRegression hmmr("0.9 0.1", transition, emissions);
// Train the HMM (supply a state sequence to perform supervised training)
std::vector<arma::mat> predictorsSeq(1, predictors);
std::vector< arma::vec> responsesSeq(1, responses);
hmmr.Train(predictorsSeq, responsesSeq);
hmm.Train(observations, states);

Once initialized, the HMMR can evaluate the probability of a certain sequence (with LogLikelihood()), predict the most likely sequence of hidden states (with Predict()), estimate the probabilities of each state for a sequence of observations (with Estimate()), or perform filtering or smoothing of observations.

Definition at line 69 of file hmm_regression.hpp.

Constructor & Destructor Documentation

◆ HMMRegression() [1/2]

HMMRegression ( const size_t  states,
const distribution::RegressionDistribution  emissions,
const double  tolerance = 1e-5 
)
inline

Create the Hidden Markov Model Regression with the given number of hidden states and the given default regression emission.

The dimensionality of the observations is taken from the emissions variable, so it is important that the given default emission distribution is set with the correct dimensionality. Alternately, set the dimensionality with Dimensionality(). Optionally, the tolerance for convergence of the Baum-Welch algorithm can be set.

By default, the transition matrix and initial probability vector are set to contain equal probability for each state.

Parameters
statesNumber of states.
emissionsDefault distribution for emissions.
toleranceTolerance for convergence of training algorithm (Baum-Welch).

Definition at line 89 of file hmm_regression.hpp.

◆ HMMRegression() [2/2]

HMMRegression ( const arma::vec &  initial,
const arma::mat &  transition,
const std::vector< distribution::RegressionDistribution > &  emission,
const double  tolerance = 1e-5 
)
inline

Create the Hidden Markov Model Regression with the given initial probability vector, the given transition matrix, and the given regression emission distributions.

The dimensionality of the observations of the HMMR are taken from the given emission distributions. Alternately, the dimensionality can be set with Dimensionality().

The initial state probability vector should have length equal to the number of states, and each entry represents the probability of being in the given state at time T = 0 (the beginning of a sequence).

The transition matrix should be such that T(i, j) is the probability of transition to state i from state j. The columns of the matrix should sum to 1.

Optionally, the tolerance for convergence of the Baum-Welch algorithm can be set.

Parameters
initialInitial state probabilities.
transitionTransition matrix.
emissionEmission distributions.
toleranceTolerance for convergence of training algorithm (Baum-Welch).

Definition at line 119 of file hmm_regression.hpp.

References HMMRegression::Estimate(), HMMRegression::Filter(), HMMRegression::LogLikelihood(), HMMRegression::Predict(), HMMRegression::Smooth(), and HMMRegression::Train().

Member Function Documentation

◆ Estimate() [1/2]

double Estimate ( const arma::mat &  predictors,
const arma::vec &  responses,
arma::mat &  stateProb,
arma::mat &  forwardProb,
arma::mat &  backwardProb,
arma::vec &  scales 
) const

Estimate the probabilities of each hidden state at each time step for each given data observation, using the Forward-Backward algorithm.

Each matrix which is returned has columns equal to the number of data observations, and rows equal to the number of hidden states in the model. The log-likelihood of the most probable sequence is returned.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
stateProbMatrix in which the probabilities of each state at each time interval will be stored.
forwardProbMatrix in which the forward probabilities of each state at each time interval will be stored.
backwardProbMatrix in which the backward probabilities of each state at each time interval will be stored.
scalesVector in which the scaling factors at each time interval will be stored.
Returns
Log-likelihood of most likely state sequence.

Referenced by HMMRegression::HMMRegression().

◆ Estimate() [2/2]

double Estimate ( const arma::mat &  predictors,
const arma::vec &  responses,
arma::mat &  stateProb 
) const

Estimate the probabilities of each hidden state at each time step of each given data observation, using the Forward-Backward algorithm.

The returned matrix of state probabilities has columns equal to the number of data observations, and rows equal to the number of hidden states in the model. The log-likelihood of the most probable sequence is returned.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
stateProbProbabilities of each state at each time interval.
Returns
Log-likelihood of most likely state sequence.

◆ Filter()

void Filter ( const arma::mat &  predictors,
const arma::vec &  responses,
arma::vec &  filterSeq,
size_t  ahead = 0 
) const

HMMR filtering.

Computes the k-step-ahead expected response at each time conditioned only on prior observations. That is E{ Y[t+k] | Y[0], ..., Y[t] }. The returned matrix has columns equal to the number of observations. Note that the expectation may not be meaningful for discrete emissions.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
aheadNumber of steps ahead (k) for expectations.
filterSeqVector in which the expected emission sequence will be stored.

Referenced by HMMRegression::HMMRegression().

◆ LogLikelihood()

double LogLikelihood ( const arma::mat &  predictors,
const arma::vec &  responses 
) const

Compute the log-likelihood of the given predictors and responses.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
Returns
Log-likelihood of the given sequence.

Referenced by HMMRegression::HMMRegression().

◆ Predict()

double Predict ( const arma::mat &  predictors,
const arma::vec &  responses,
arma::Row< size_t > &  stateSeq 
) const

Compute the most probable hidden state sequence for the given predictors and responses, using the Viterbi algorithm, returning the log-likelihood of the most likely state sequence.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
stateSeqVector in which the most probable state sequence will be stored.
Returns
Log-likelihood of most probable state sequence.

Referenced by HMMRegression::HMMRegression().

◆ Smooth()

void Smooth ( const arma::mat &  predictors,
const arma::vec &  responses,
arma::vec &  smoothSeq 
) const

HMM smoothing.

Computes expected emission at each time conditioned on all observations. That is E{ Y[t] | Y[0], ..., Y[T] }. The returned matrix has columns equal to the number of observations. Note that the expectation may not be meaningful for discrete emissions.

Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
smoothSeqVector in which the expected emission sequence will be stored.

Referenced by HMMRegression::HMMRegression().

◆ Train() [1/2]

void Train ( const std::vector< arma::mat > &  predictors,
const std::vector< arma::vec > &  responses 
)

Train the model using the Baum-Welch algorithm, with only the given predictors and responses.

Instead of giving a guess transition and emission here, do that in the constructor. Each matrix in the vector of predictors corresponds to an individual data sequence, and likewise for each vec in the vector of responses. The number of rows in each matrix of predictors plus one should be equal to the dimensionality of the HMM (which is set in the constructor).

It is preferable to use the other overload of Train(), with labeled data. That will produce much better results. However, if labeled data is unavailable, this will work. In addition, it is possible to use Train() with labeled data first, and then continue to train the model using this overload of Train() with unlabeled data.

The tolerance of the Baum-Welch algorithm can be set either in the constructor or with the Tolerance() method. When the change in log-likelihood of the model between iterations is less than the tolerance, the Baum-Welch algorithm terminates.

Note
Train() can be called multiple times with different sequences; each time it is called, it uses the current parameters of the HMM as a starting point for training.
Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.

Referenced by HMMRegression::HMMRegression().

◆ Train() [2/2]

void Train ( const std::vector< arma::mat > &  predictors,
const std::vector< arma::vec > &  responses,
const std::vector< arma::Row< size_t > > &  stateSeq 
)

Train the model using the given labeled observations; the transition and regression emissions are directly estimated.

Each matrix in the vector of predictors corresponds to an individual data sequence, and likewise for each vec in the vector of responses. The number of rows in each matrix of predictors plus one should be equal to the dimensionality of the HMM (which is set in the constructor).

Note
Train() can be called multiple times with different sequences; each time it is called, it uses the current parameters of the HMMR as a starting point for training.
Parameters
predictorsVector of predictor sequences.
responsesVector of response sequences.
stateSeqVector of state sequences, corresponding to each observation.

The documentation for this class was generated from the following file: