API reference
HiddenMarkovModels — ModuleHiddenMarkovModelsA Julia package for HMM modeling, simulation, inference and learning.
Exports
Sequence formatting
Most algorithms below ingest the data with two positional arguments obs_seq (mandatory) and control_seq (optional), and a keyword argument seq_ends (optional).
- If the data consists of a single sequence,
obs_seqandcontrol_seqare the corresponding vectors of observations and controls, and you don't need to provideseq_ends. - If the data consists of multiple sequences,
obs_seqandcontrol_seqare concatenations of several vectors, whose end indices are given byseq_ends. Starting from separate sequencesobs_seqsandcontrol_seqs, you can run the following snippet:
obs_seq = reduce(vcat, obs_seqs)
control_seq = reduce(vcat, control_seqs)
seq_ends = cumsum(length.(obs_seqs))Types
HiddenMarkovModels.AbstractHMM — TypeAbstractHMMAbstract supertype for an HMM amenable to simulation, inference and learning.
Interface
To create your own subtype of AbstractHMM, you need to implement the following methods:
initializationtransition_matrixobs_distributionsfit!(for learning)
Applicable functions
Any AbstractHMM which satisfies the interface can be given to the following functions:
randlogdensityofforwardviterbiforward_backwardbaum_welch(if[fit!](@ref)is implemented)
HiddenMarkovModels.HMM — Typestruct HMM{V<:(AbstractVector), M<:(AbstractMatrix), VD<:(AbstractVector), Vl<:(AbstractVector), Ml<:(AbstractMatrix)} <: AbstractHMMBasic implementation of an HMM.
Fields
init::AbstractVector: initial state probabilitiestrans::AbstractMatrix: state transition probabilitiesdists::AbstractVector: observation distributionsloginit::AbstractVector: logarithms of initial state probabilitieslogtrans::AbstractMatrix: logarithms of state transition probabilities
Interface
HiddenMarkovModels.initialization — Functioninitialization(hmm)Return the vector of initial state probabilities for hmm.
HiddenMarkovModels.transition_matrix — Functiontransition_matrix(hmm)
transition_matrix(hmm, control)Return the matrix of state transition probabilities for hmm (possibly when control is applied).
When processing sequences, the control at time t influences the transition from time t to t+1 (and not from time t-1 to t).
HiddenMarkovModels.obs_distributions — Functionobs_distributions(hmm)
obs_distributions(hmm, control)Return a vector of observation distributions, one for each state of hmm (possibly when control is applied).
These distribution objects should implement
Random.rand(rng, dist)for samplingDensityInterface.logdensityof(dist, obs)for inferenceStatsAPI.fit!(dist, obs_seq, weight_seq)for learning
Utils
Base.length — Functionlength(hmm)Return the number of states of hmm.
Base.rand — Functionrand([rng,] hmm, T)
rand([rng,] hmm, control_seq)Simulate hmm for T time steps, or when the sequence control_seq is applied.
Return a named tuple (; state_seq, obs_seq).
Base.eltype — Functioneltype(hmm, obs, control)Return a type that can accommodate forward-backward computations for hmm on observations similar to obs.
It is typically a promotion between the element type of the initialization, the element type of the transition matrix, and the type of an observation logdensity evaluated at obs.
HiddenMarkovModels.seq_limits — Functionseq_limits(seq_ends, k)
Return a tuple (t1, t2) giving the begin and end indices of subsequence k within a set of sequences ending at seq_ends.
Inference
DensityInterface.logdensityof — Functionlogdensityof(hmm)Return the prior loglikelihood associated with the parameters of hmm.
logdensityof(hmm, obs_seq; ...)
logdensityof(hmm, obs_seq, control_seq; seq_ends)
Run the forward algorithm to compute the loglikelihood of obs_seq for hmm, integrating over all possible state sequences.
HiddenMarkovModels.joint_logdensityof — Functionjoint_logdensityof(hmm, obs_seq, state_seq; ...)
joint_logdensityof(
hmm,
obs_seq,
state_seq,
control_seq;
seq_ends
)
Run the forward algorithm to compute the the joint loglikelihood of obs_seq and state_seq for hmm.
HiddenMarkovModels.forward — Functionforward(hmm, obs_seq; ...)
forward(
hmm,
obs_seq,
control_seq;
seq_ends,
error_if_not_finite
)
Apply the forward algorithm to infer the current state after sequence obs_seq for hmm.
Return a tuple (storage.α, storage.logL) where storage is of type ForwardStorage.
HiddenMarkovModels.viterbi — Functionviterbi(hmm, obs_seq; ...)
viterbi(hmm, obs_seq, control_seq; seq_ends)
Apply the Viterbi algorithm to infer the most likely state sequence corresponding to obs_seq for hmm.
Return a tuple (storage.q, storage.logL) where storage is of type ViterbiStorage.
HiddenMarkovModels.forward_backward — Functionforward_backward(hmm, obs_seq; ...)
forward_backward(hmm, obs_seq, control_seq; seq_ends)
Apply the forward-backward algorithm to infer the posterior state and transition marginals during sequence obs_seq for hmm.
Return a tuple (storage.γ, storage.logL) where storage is of type ForwardBackwardStorage.
Learning
HiddenMarkovModels.baum_welch — Functionbaum_welch(hmm_guess, obs_seq; ...)
baum_welch(
hmm_guess,
obs_seq,
control_seq;
seq_ends,
atol,
max_iterations,
loglikelihood_increasing
)
Apply the Baum-Welch algorithm to estimate the parameters of an HMM on obs_seq, starting from hmm_guess.
Return a tuple (hmm_est, loglikelihood_evolution) where hmm_est is the estimated HMM and loglikelihood_evolution is a vector of loglikelihood values, one per iteration of the algorithm.
Keyword arguments
atol: minimum loglikelihood increase at an iteration of the algorithm (otherwise the algorithm is deemed to have converged)max_iterations: maximum number of iterations of the algorithmloglikelihood_increasing: whether to throw an error if the loglikelihood decreases
StatsAPI.fit! — FunctionStatsAPI.fit!(
hmm, fb_storage::ForwardBackwardStorage,
obs_seq, [control_seq]; seq_ends,
)Update hmm in-place based on information generated during forward-backward.
This function is allowed to reuse fb_storage as a scratch space, so its contents should not be trusted afterwards.
In-place versions
Forward
HiddenMarkovModels.ForwardStorage — Typestruct ForwardStorage{R}Fields
Only the fields with a description are part of the public API.
α::Matrix: posterior last state marginalsα[i] = ℙ(X[T]=i | Y[1:T])logL::Vector: one loglikelihood per observation sequenceB::Matrixc::Vector
HiddenMarkovModels.initialize_forward — Functioninitialize_forward(hmm, obs_seq, control_seq; seq_ends)
HiddenMarkovModels.forward! — Functionforward!(
storage,
hmm,
obs_seq,
control_seq;
seq_ends,
error_if_not_finite
)
Viterbi
HiddenMarkovModels.ViterbiStorage — Typestruct ViterbiStorage{R}Fields
Only the fields with a description are part of the public API.
q::Vector{Int64}: most likely state sequenceq[t] = argmaxᵢ ℙ(X[t]=i | Y[1:T])logL::Vector: one joint loglikelihood per pair of observation sequence and most likely state sequencelogB::Matrixϕ::Matrixψ::Matrix{Int64}
HiddenMarkovModels.initialize_viterbi — Functioninitialize_viterbi(hmm, obs_seq, control_seq; seq_ends)
HiddenMarkovModels.viterbi! — Functionviterbi!(storage, hmm, obs_seq, control_seq; seq_ends)
Forward-backward
HiddenMarkovModels.ForwardBackwardStorage — Typestruct ForwardBackwardStorage{R, M<:AbstractArray{R, 2}}Fields
Only the fields with a description are part of the public API.
γ::Matrix: posterior state marginalsγ[i,t] = ℙ(X[t]=i | Y[1:T])ξ::Vector{M} where {R, M<:AbstractMatrix{R}}: posterior transition marginalsξ[t][i,j] = ℙ(X[t]=i, X[t+1]=j | Y[1:T])logL::Vector: one loglikelihood per observation sequenceB::Matrixα::Matrixc::Vectorβ::MatrixBβ::Matrix
HiddenMarkovModels.initialize_forward_backward — Functioninitialize_forward_backward(
hmm,
obs_seq,
control_seq;
seq_ends,
transition_marginals
)
HiddenMarkovModels.forward_backward! — Functionforward_backward!(
storage,
hmm,
obs_seq,
control_seq;
seq_ends,
transition_marginals
)
Baum-Welch
HiddenMarkovModels.baum_welch! — Functionbaum_welch!(
fb_storage,
logL_evolution,
hmm,
obs_seq,
control_seq;
seq_ends,
atol,
max_iterations,
loglikelihood_increasing
)
Miscellaneous
HiddenMarkovModels.valid_hmm — Functionvalid_hmm(hmm)Perform some checks to rule out obvious inconsistencies with an AbstractHMM object.
HiddenMarkovModels.rand_prob_vec — Functionrand_prob_vec([rng, ::Type{R},] N)Generate a random probability distribution of size N with normalized uniform entries.
HiddenMarkovModels.rand_trans_mat — Functionrand_trans_mat([rng, ::Type{R},] N)Generate a random transition matrix of size (N, N) with normalized uniform entries.
HiddenMarkovModels.fit_in_sequence! — Functionfit_in_sequence!(dists, i, x, w)
Modify the i-th element of dists by fitting it to an observation sequence x with associated weight sequence w.
Default behavior:
fit!(dists[i], x, w)Override for Distributions.jl (in the package extension)
dists[i] = fit(eltype(dists), turn_into_vector(x), w)Internals
HiddenMarkovModels.LightDiagNormal — Typestruct LightDiagNormal{T1, T2, T3, V1<:AbstractArray{T1, 1}, V2<:AbstractArray{T2, 1}, V3<:AbstractArray{T3, 1}}An HMMs-compatible implementation of a multivariate normal distribution with diagonal covariance, enabling allocation-free in-place estimation.
This is not part of the public API and is expected to change.
Fields
μ::AbstractVector: meansσ::AbstractVector: standard deviationslogσ::AbstractVector: log standard deviations
HiddenMarkovModels.LightCategorical — Typestruct LightCategorical{T1, T2, V1<:AbstractArray{T1, 1}, V2<:AbstractArray{T2, 1}}An HMMs-compatible implementation of a discrete categorical distribution, enabling allocation-free in-place estimation.
This is not part of the public API and is expected to change.
Fields
p::AbstractVector: class probabilitieslogp::AbstractVector: log class probabilities
HiddenMarkovModels.log_initialization — Functionlog_initialization(hmm)Return the vector of initial state log-probabilities for hmm.
Falls back on initialization.
HiddenMarkovModels.log_transition_matrix — Functionlog_transition_matrix(hmm)
log_transition_matrix(hmm, control)Return the matrix of state transition log-probabilities for hmm (possibly when control is applied).
Falls back on transition_matrix.
When processing sequences, the control at time t influences the transition from time t-1 to t (since version 0.7 of the package).
HiddenMarkovModels.mul_rows_cols! — Functionmul_rows_cols!(B, l, A, r)Perform the in-place operation B .= l .* A .* transpose(r).
HiddenMarkovModels.argmaxplus_transmul! — Functionargmaxplus_transmul!(y, ind, A, x)Perform the in-place multiplication transpose(A) * x in the sense of max-plus algebra, store the result in y, and store the index of the maximum for each component of y in ind.
Index
HiddenMarkovModelsHiddenMarkovModels.AbstractHMMHiddenMarkovModels.ForwardBackwardStorageHiddenMarkovModels.ForwardStorageHiddenMarkovModels.HMMHiddenMarkovModels.LightCategoricalHiddenMarkovModels.LightDiagNormalHiddenMarkovModels.ViterbiStorageBase.eltypeBase.lengthBase.randDensityInterface.logdensityofHiddenMarkovModels.argmaxplus_transmul!HiddenMarkovModels.baum_welchHiddenMarkovModels.baum_welch!HiddenMarkovModels.fit_in_sequence!HiddenMarkovModels.forwardHiddenMarkovModels.forward!HiddenMarkovModels.forward_backwardHiddenMarkovModels.forward_backward!HiddenMarkovModels.initializationHiddenMarkovModels.initialize_forwardHiddenMarkovModels.initialize_forward_backwardHiddenMarkovModels.initialize_viterbiHiddenMarkovModels.joint_logdensityofHiddenMarkovModels.log_initializationHiddenMarkovModels.log_transition_matrixHiddenMarkovModels.mul_rows_cols!HiddenMarkovModels.obs_distributionsHiddenMarkovModels.rand_prob_vecHiddenMarkovModels.rand_trans_matHiddenMarkovModels.seq_limitsHiddenMarkovModels.transition_matrixHiddenMarkovModels.valid_hmmHiddenMarkovModels.viterbiHiddenMarkovModels.viterbi!StatsAPI.fit!