MosaicsHMM-class: Class "MosaicsHMM"

Description Objects from the Class Slots Methods Author(s) References See Also Examples

Description

This class represents MOSAiCS-HMM model fit.

Objects from the Class

Objects can be created by calls of the form new("MosaicsHMM", ...).

Slots

HMMfit:

Object of class "list", representing the fitted MOSAiCS-HMM model.

chrID:

Object of class "character", a vector of chromosome IDs.

coord:

Object of class "numeric", a vector of genomic coordinates.

tagCount:

Object of class "numeric", a vector of tag counts of ChIP sample.

mappability:

Object of class "numeric", a vector of mappability score.

gcContent:

Object of class "numeric", a vector of GC content score.

input:

Object of class "numeric", a vector of tag counts of matched control sample.

inputdata:

Object of class "list", representing the bin-level data.

mosaicsEst:

Object of class "MosaicsFitEst", representing estimates of MOSAiCS model fit.

init:

Object of class "character", representing the approach to initialize MOSAiCS-HMM.

initPiMat:

Object of class "numeric", representing initial transition matrix.

peakParam:

Object of class "MosaicsPeakParam", representing parameters for peak calling.

binsize:

Object of class "numeric", representing size of a bin.

nRatio:

Object of class "numeric", representing ratio of sequencing depth of ChIP vs. control.

bicMosaics:

Object of class "numeric", representing the BIC value of MOSAiCS fit.

bicMosaicsHMM:

Object of class "numeric", representing the BIC value of MOSAiCS-HMM fit.

seqDepth:

Object of class "numeric", a vector of sequencing depth of length 2, where the first and second elements correpond to sequencing depths of ChIP and control samples, respectively. If there is not control sample, the second element is set to NA.

Methods

estimates

signature(object = "MosaicsHMM"): extract estimates from MOSAiCS-HMM model fit.

=

plot

signature(x = "MosaicsHMM", y = "missing", seed=12345, parallel=FALSE, nCore=8 ): draw Goodness of Fit (GOF) plot. You can specify random seed in seed. If parallel=TRUE, parallel computing is utilized to simulate data, where nCore indicates CPUs used for parallel computing.

print

signature(x = "MosaicsHMM"): (not supported yet)

show

signature(object = "MosaicsHMM"): provide brief summary of the object.

seqDepth

signature(object = "MosaicsHMM"): provide a vector of sequencing depth.

Author(s)

Dongjun Chung, Pei Fen Kuan, Rene Welch, Sunduz Keles

References

Kuan, PF, D Chung, G Pan, JA Thomson, R Stewart, and S Keles (2011), "A Statistical Framework for the Analysis of ChIP-Seq Data", Journal of the American Statistical Association, Vol. 106, pp. 891-903.

Chung, D, Zhang Q, and Keles S (2014), "MOSAiCS-HMM: A model-based approach for detecting regions of histone modifications from ChIP-seq data", Datta S and Nettleton D (eds.), Statistical Analysis of Next Generation Sequencing Data, Springer.

See Also

mosaicsFitHMM, mosaicsPeakHMM.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
showClass("MosaicsHMM")
## Not run: 
library(mosaicsExample)

constructBins( infile=system.file( file.path("extdata","wgEncodeBroadHistoneGm12878H3k4me3StdAlnRep1_chr22_sorted.bam"), package="mosaicsExample"), 
    fileFormat="bam", outfileLoc="~/", 
    byChr=FALSE, useChrfile=FALSE, chrfile=NULL, excludeChr=NULL, 
    PET=FALSE, fragLen=200, binSize=200, capping=0 )
constructBins( infile=system.file( file.path("extdata","wgEncodeBroadHistoneGm12878ControlStdAlnRep1_chr22_sorted.bam"), package="mosaicsExample"), 
    fileFormat="bam", outfileLoc="~/", 
    byChr=FALSE, useChrfile=FALSE, chrfile=NULL, excludeChr=NULL, 
    PET=FALSE, fragLen=200, binSize=200, capping=0 )

binHM <- readBins( type=c("chip","input"),
    fileName=c( "~/wgEncodeBroadHistoneGm12878H3k4me3StdAlnRep1_chr22_sorted.bam_fragL200_bin200.txt",
    "~/wgEncodeBroadHistoneGm12878ControlStdAlnRep1_chr22_sorted.bam_fragL200_bin200.txt" ) )
binHM
plot(binHM)
plot( binHM, plotType="input" )

fitHM <- mosaicsFit( binHM, analysisType="IO", bgEst="rMOM" )
fitHM
plot(fitHM)

hmmHM <- mosaicsFitHMM( fitHM, signalModel = "2S", 
  init="mosaics", init.FDR = 0.05, parallel=TRUE, nCore=8 )
hmmHM
plot(hmmHM)

peakHM <- mosaicsPeakHMM( hmmHM, FDR = 0.05, decoding="posterior",
  thres=10, parallel=TRUE, nCore=8 )

peakHM <- extractReads( peakHM,
  chipFile=system.file( file.path("extdata","wgEncodeBroadHistoneGm12878H3k4me3StdAlnRep1_chr22_sorted.bam"), package="mosaicsExample"),
  chipFileFormat="bam", chipPET=FALSE, chipFragLen=200,
  controlFile=system.file( file.path("extdata","wgEncodeBroadHistoneGm12878ControlStdAlnRep1_chr22_sorted.bam"), package="mosaicsExample"), 
  controlFileFormat="bam", controlPET=FALSE, controlFragLen=200, parallel=TRUE, nCore=8 )
peakHM

peakHM <- findSummit( peakHM, parallel=TRUE, nCore=8 )
head(print(peakHM))
plot( peakHM, filename="~/peakplot_HM.pdf" )

peakHM <- adjustBoundary( peakHM, parallel=TRUE, nCore=8 )
peakHM
head(print(peakHM))

peakHM <- filterPeak( peakHM, parallel=TRUE, nCore=8 )
peakHM
head(print(peakHM))

export( peakHM, type = "txt", filename = "./peakHM.txt" )
export( peakHM, type = "bed", filename = "./peakHM.bed" )
export( peakHM, type = "gff", filename = "./peakHM.gff" )
export( peakHM, type = "narrowPeak", filename = "./peakHM.narrowPeak" )
export( peakHM, type = "broadPeak", filename = "./peakHM.broadPeak" )

## End(Not run)

dongjunchung/mosaics documentation built on March 1, 2020, 3:44 a.m.