Algorithms

See common options for the documentation of a, A, L, ...

Baum–Welch

Distributions.fit_mleFunction
fit_mle(hmm, observations; ...) -> AbstractHMM

Estimate the HMM parameters using the EM (Baum-Welch) algorithm, with hmm as the initial state.

Keyword Arguments

  • display::Symbol = :none: when to display convergence logs, can be set to :iter or :final.
  • init::Symbol = :none: if set to :kmeans the HMM parameters will be initialized using a K-means clustering.
  • maxiter::Integer = 100: maximum number of iterations to perform.
  • tol::Integer = 1e-3: stop the algorithm when the improvement in the log-likelihood is less than tol.

Output

  • <:AbstractHMM: a copy of the original HMM with the updated parameters.

Forward-Backward

HMMBase.forwardFunction
forward(a, A, L; logl) -> (Vector, Float)

Compute forward probabilities using samples likelihoods. See Forward-backward algorithm.

Output

  • Vector{Float64}: forward probabilities.
  • Float64: log-likelihood of the observed sequence.
forward(hmm, observations; logl, robust) -> (Vector, Float)

Compute forward probabilities of the observations given the hmm model.

Output

  • Vector{Float64}: forward probabilities.
  • Float64: log-likelihood of the observed sequence.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000)
probs, tot = forward(hmm, y)
HMMBase.backwardFunction
backward(a, A, L; logl) -> (Vector, Float)

Compute backward probabilities using samples likelihoods. See Forward-backward algorithm.

Output

  • Vector{Float64}: backward probabilities.
  • Float64: log-likelihood of the observed sequence.
backward(hmm, observations; logl, robust) -> (Vector, Float)

Compute backward probabilities of the observations given the hmm model.

Output

  • Vector{Float64}: backward probabilities.
  • Float64: log-likelihood of the observed sequence.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000)
probs, tot = backward(hmm, y)
HMMBase.posteriorsFunction
posteriors(α, β) -> Vector

Compute posterior probabilities from α and β.

Arguments

  • α::AbstractVector: forward probabilities.
  • β::AbstractVector: backward probabilities.
posteriors(a, A, L; logl) -> Vector

Compute posterior probabilities using samples likelihoods.

posteriors(hmm, observations; logl, robust) -> Vector

Compute posterior probabilities using samples likelihoods.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000)
γ = posteriors(hmm, y)

Viterbi

HMMBase.viterbiFunction
viterbi(a, A, L; logl) -> Vector

Find the most likely hidden state sequence, see Viterbi algorithm.

viterbi(hmm, observations; logl, robust) -> Vector

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000)
zv = viterbi(hmm, y)