| selectByOverlap {methylKit} | R Documentation |
The function selects records from a methylBaseDB,
methylRawDB or methylDiffDB object
that lie inside the regions given by ranges of class GRanges
and returns an object of class
methylBase, methylRaw or methylDiff
selectByOverlap(object,ranges) ## S4 method for signature 'methylRaw' selectByOverlap(object, ranges) ## S4 method for signature 'methylRawList' selectByOverlap(object, ranges) ## S4 method for signature 'methylBase' selectByOverlap(object, ranges) ## S4 method for signature 'methylDiff' selectByOverlap(object, ranges) ## S4 method for signature 'methylRawDB' selectByOverlap(object, ranges) ## S4 method for signature 'methylRawListDB' selectByOverlap(object, ranges) ## S4 method for signature 'methylBaseDB' selectByOverlap(object, ranges) ## S4 method for signature 'methylDiffDB' selectByOverlap(object, ranges)
object |
an |
ranges |
a GRanges object specifying the regions of interest |
a methylBase,methylRaw or
methylDiff object depending on the input object.
Alexander Gosdschan
data(methylKit)
file.list=list( system.file("extdata", "test1.myCpG.txt", package = "methylKit"),
system.file("extdata", "test2.myCpG.txt", package = "methylKit"),
system.file("extdata", "control1.myCpG.txt", package = "methylKit"),
system.file("extdata", "control2.myCpG.txt", package = "methylKit") )
methylRawListDB.obj=methRead(file.list,
sample.id=list("test1","test2","ctrl1","ctrl2"),
assembly="hg18",treatment=c(1,1,0,0),
dbtype = "tabix",dbdir = "methylDB")
methylBaseDB.obj=unite(methylRawListDB.obj)
methylDiffDB.obj = calculateDiffMeth(methylBaseDB.obj)
# define the windows of interest as a GRanges object, this can be any set
# of genomic locations
library(GenomicRanges)
my.win=GRanges(seqnames="chr21",
ranges=IRanges(start=seq(from=9764513,by=10000,length.out=20),width=5000) )
# selects the records that lie inside the regions
myRaw <- selectByOverlap(methylRawListDB.obj[[1]],my.win)
# selects the records that lie inside the regions
myBase <- selectByOverlap(methylBaseDB.obj,my.win)
# selects the records that lie inside the regions
myDiff <- selectByOverlap(methylDiffDB.obj,my.win)
rm(methylRawListDB.obj)
rm(methylBaseDB.obj)
rm(methylDiffDB.obj)
unlink("methylDB",recursive=TRUE)