Last data update: 2014.03.03

R: Generalized negative log likelihood and Viterbi algorithms
negLogLikeR Documentation

Generalized negative log likelihood and Viterbi algorithms

Description

negLogLike: Returns the negative log likelihood calculated with the forward equations.

viterbiPath: Calculates the most likely sequence of hidden states for the Markov model given the current parameters.

Usage

  negLogLike(par,fx.par,data,nstates,stFn,trFn,emFn)
  viterbiPath(par,fx.par,data,nstates,stFn,trFn,emFn)

Arguments

par

A list of parameters, over which the likelihood will be optimized.

fx.par

A list of fixed parameters.

data

A list of data objects, which must contain a vector O, which represents the observed sequence of the HMM.

nstates

The number of states of the HMM.

stFn

A function which takes arguments par, fx.par, data, and nstates, and returns a vector of length nstates of starting probabilities.

trFn

A function which takes arguments par, fx.par, data, and nstates, and returns a matrix of dimension (nstates,nstates) of the transition probabilities.

emFn

A function which takes arguments par, fx.par, data, and nstates, and returns a matrix of dimension (nstates,length(O)) of the emission probabilities.

Value

negLogLike: The negative log likelihood of the HMM. The likelihood is slightly modified to account for ranges with read counts which have zero probability of originating from any of the states. In this case the likelihood is lowered and the range is skipped.

viterbiPath: The Viterbi path through the states given the parameters.

References

On the forward equations and the Viterbi algorithm:

Rabiner, L. R. (1989): "A tutorial on hidden Markov models and selected applications in speech recognition," Proceedings of the IEEE, 77, 257, 286, http://dx.doi.org/10.1109/5.18626.

Examples

## functions for starting, transition, and emission probabilities
stFn <- function(par,fx.par,data,nstates) rep(1/nstates,nstates)
trFn <- function(par,fx.par,data,nstates) {
  A <- matrix(1/(nstates*10),ncol=nstates,nrow=nstates)
  diag(A) <- 1 - rowSums(A)
  A
}
emFn <- function(par,fx.par,data,nstates) {
  t(sapply(1:nstates,function(j) dnorm(data$O,par$means[j],fx.par$sdev)))
}

## simulate some observations from two states
Q <- c(rep(1,100),rep(2,100),rep(1,100),rep(2,100))
T <- length(Q)
means <- c(-0.5,0.5)
sdev <- 1
O <- rnorm(T,means[Q],sdev)

## use viterbiPath() to recover the state chain using parameters
viterbi.path <- viterbiPath(par=list(means=means),
fx.par=list(sdev=sdev), data=list(O=O), nstates=2,stFn,trFn,emFn) 
plot(O,pch=Q,col=c("darkgreen","orange")[viterbi.path])

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(exomeCopy)
Loading required package: IRanges
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: GenomicRanges
Loading required package: GenomeInfoDb
Loading required package: Rsamtools
Loading required package: Biostrings
Loading required package: XVector
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/exomeCopy/negLogLike.Rd_%03d_medium.png", width=480, height=480)
> ### Name: negLogLike
> ### Title: Generalized negative log likelihood and Viterbi algorithms
> ### Aliases: negLogLike viterbiPath
> 
> ### ** Examples
> 
> ## functions for starting, transition, and emission probabilities
> stFn <- function(par,fx.par,data,nstates) rep(1/nstates,nstates)
> trFn <- function(par,fx.par,data,nstates) {
+   A <- matrix(1/(nstates*10),ncol=nstates,nrow=nstates)
+   diag(A) <- 1 - rowSums(A)
+   A
+ }
> emFn <- function(par,fx.par,data,nstates) {
+   t(sapply(1:nstates,function(j) dnorm(data$O,par$means[j],fx.par$sdev)))
+ }
> 
> ## simulate some observations from two states
> Q <- c(rep(1,100),rep(2,100),rep(1,100),rep(2,100))
> T <- length(Q)
> means <- c(-0.5,0.5)
> sdev <- 1
> O <- rnorm(T,means[Q],sdev)
> 
> ## use viterbiPath() to recover the state chain using parameters
> viterbi.path <- viterbiPath(par=list(means=means),
+ fx.par=list(sdev=sdev), data=list(O=O), nstates=2,stFn,trFn,emFn) 
> plot(O,pch=Q,col=c("darkgreen","orange")[viterbi.path])
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>