[][src]Module bio::stats::hmm

An implementation of Hidden Markov Models in Rust.

Examples

Discrete Emission Distribution

We construct the example from Borodovsky & Ekisheva (2006), pp. 80 (also see these slides.

#[macro_use] extern crate approx;
#[macro_use] extern crate ndarray;

use bio::stats::Prob;
use bio::stats::hmm::discrete_emission::Model as DiscreteEmissionHMM;
use bio::stats::hmm::viterbi;

let transition = array![[0.5, 0.5], [0.4, 0.6]];
let observation = array![[0.2, 0.3, 0.3, 0.2], [0.3, 0.2, 0.2, 0.3]];
let initial = array![0.5, 0.5];

let hmm = DiscreteEmissionHMM::with_float(&transition, &observation, &initial)
    .expect("Dimensions should be consistent");
let (path, log_prob) = viterbi(&hmm, &vec![2, 2, 1, 0, 1, 3, 2, 0, 0]);
let prob = Prob::from(log_prob);
assert_relative_eq!(4.25e-8_f64, *prob, epsilon = 1e-9_f64);

Continuous (Gaussian) Emission Distribution

#[macro_use] extern crate approx;
#[macro_use] extern crate ndarray;

use bio::stats::Prob;
use bio::stats::hmm::univariate_continuous_emission::GaussianModel as GaussianHMM;
use statrs::distribution::Normal;
use bio::stats::hmm::viterbi;

let transition = array![[0.5, 0.5], [0.4, 0.6]];
let observation = vec![
    Normal::new(0.0, 1.0).unwrap(),
    Normal::new(2.0, 1.0).unwrap(),
];
let initial = array![0.5, 0.5];

let hmm = GaussianHMM::with_float(&transition, observation, &initial)
    .expect("Dimensions should be consistent");
let (path, log_prob) = viterbi(
    &hmm,
    &vec![-0.1, 0.1, -0.2, 0.5, 0.8, 1.1, 1.2, 1.5, 0.5, 0.2],
);
let prob = Prob::from(log_prob);
assert_relative_eq!(2.64e-8_f64, *prob, epsilon = 1e-9_f64);

Numeric Stability

The implementation uses log-scale probabilities for numeric stability.

Limitations

Currently, only discrete and single-variate Gaussian continuous HMMs are implemented. Also, only dense transition matrices are supported.

References

  • Rabiner, Lawrence R. "A tutorial on hidden Markov models and selected applications in speech recognition." Proceedings of the IEEE 77, no. 2 (1989): 257-286.

Modules

discrete_emission

Implementation of Hidden Markov Model with emission values from discrete distributions.

univariate_continuous_emission

Implementation of Hidden Markov Models with emission values from univariate continuous distributions.

Structs

State

A newtype for HMM states.

StateIter

Iterate over the states of a Model.

StateTransition

Transition between two states in a Model.

StateTransitionIter

Iterate over all state transitions of a Model.

Enums

HMMError

Traits

Model

A trait for Hidden Markov Models (HMM) with generic Observation type.

Functions

backward

Execute the backward algorithm and return the backward probabilities as LogProb values and the resulting backward probability.

forward

Execute the forward algorithm and return the forward probabilites as LogProb values and the resulting forward probability.

viterbi

Execute Viterbi algorithm on the given slice of Observation values to get the maximum a posteriori (MAP) probability.