Model

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

HMM

HMMBase.HMMType
HMM([a, ]A, B) -> HMM

Build an HMM with transition matrix A and observation distributions B. If the initial state distribution a is not specified, a uniform distribution is assumed.

Observations distributions can be of different types (for example Normal and Exponential), but they must be of the same dimension.

Arguments

  • a::AbstractVector{T}: initial probabilities vector.
  • A::AbstractMatrix{T}: transition matrix.
  • B::AbstractVector{<:Distribution{F}}: observations distributions.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
Base.randFunction
rand([rng, ]hmm, T; init, seq) -> Array | (Vector, Array)

Sample a trajectory of T timesteps from hmm.

Keyword Arguments

  • init::Integer = rand(Categorical(hmm.a)): initial state.
  • seq::Bool = false: whether to return the hidden state sequence or not.

Output

  • Vector{Int} (if seq == true): hidden state sequence.
  • Vector{Float64} (for Univariate HMMs): observations (T).
  • Matrix{Float64} (for Multivariate HMMs): observations (T x dim(obs)).

Examples

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000) # or
z, y = rand(hmm, 1000, seq = true)
size(y) # (1000,)
using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [MvNormal(ones(2)), MvNormal(ones(2))])
y = rand(hmm, 1000) # or
z, y = rand(hmm, 1000, seq = true)
size(y) # (1000, 2)
rand([rng, ]hmm, z) -> Array

Sample observations from hmm according to trajectory z.

Output

  • Vector{Float64} (for Univariate HMMs): observations (T).
  • Matrix{Float64} (for Multivariate HMMs): observations (T x dim(obs)).

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, [1, 1, 2, 2, 1])
Base.sizeFunction
size(hmm, [dim]) -> Int | Tuple

Return the number of states in hmm and the dimension of the observations.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
size(hmm)
# output
(2, 1)
HMMBase.nparamsFunction
nparams(hmm) -> Int

Return the number of free parameters in hmm, without counting the observation distributions parameters.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
nparams(hmm)
# output
3
HMMBase.permuteFunction
permute(hmm, perm) -> HMM

Permute the states of hmm according to perm.

Arguments

  • perm::Vector{<:Integer}: permutation of the states.

Example

using Distributions, HMMBase
hmm = HMM([0.8 0.2; 0.1 0.9], [Normal(0,1), Normal(10,1)])
hmm = permute(hmm, [2, 1])
hmm.A # [0.9 0.1; 0.2 0.8]
hmm.B # [Normal(10,1), Normal(0,1)]
HMMBase.istransmatFunction
istransmat(A) -> Bool

Return true if A is square and its rows sums to 1.

Base.copyFunction
copy(hmm) -> HMM

Return a copy of hmm.

Observations Likelihoods

HMMBase.likelihoodsFunction
likelihoods(hmm, observations; logl) -> Matrix

Return the likelihood per-state and per-observation.

Output

  • Matrix{Float64}: likelihoods matrix (T x K).

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
y = rand(hmm, 1000)
L = likelihoods(hmm, y)
LL = likelihoods(hmm, y, logl = true)
StatsBase.loglikelihoodFunction
loglikelihood(hmm, observations; logl, robust) -> Float64

Compute the log-likelihood of the observations under the model. This is defined as the sum of the log of the normalization coefficients in the forward filter.

Output

  • Float64: log-likelihood of the observations sequence under the model.

Example

using Distributions, HMMBase
hmm = HMM([0.9 0.1; 0.1 0.9], [Normal(0,1), Normal(10,1)])
loglikelihood(hmm, [0.15, 0.10, 1.35])
# output
-4.588183811489616

Stationnary Distribution

HMMBase.statdistsFunction
statdists(hmm) -> Vector{Vector}

Return the stationnary distribution(s) of hmm. That is, the eigenvectors of transpose(hmm.A) with eigenvalues 1.

Abstract Type

HMMBase.AbstractHMMType
AbstractHMM{F<:VariateForm}

A custom HMM type must at-least implement the following interface:

struct CustomHMM{F,T} <: AbstractHMM{F}
    a::AbstractVector{T}               # Initial state distribution
    A::AbstractMatrix{T}               # Transition matrix
    B::AbstractVector{Distribution{F}} # Observations distributions
    # Optional, custom, fields ....
end