selectByOverlap-methods: selects records of methylDB objects lying inside a GRanges...

selectByOverlapR Documentation

selects records of methylDB objects lying inside a GRanges range

Description

The function selects records from any methylKit object that lie inside the regions given by ranges of class GRanges and returns an in-memory equivalent of this object

Usage

selectByOverlap(object,ranges)

## S4 method for signature 'methylRaw,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylRawList,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylBase,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylDiff,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylRawDB,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylRawListDB,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylBaseDB,GRanges'
selectByOverlap(object, ranges)

## S4 method for signature 'methylDiffDB,GRanges'
selectByOverlap(object, ranges)

Arguments

object

an methylRaw,methylRawDB, methylRawList, methylRawListDB, methylBase, methylBaseDB, methylDiff or methylDiffDB object

ranges

a GRanges object specifying the regions of interest

Value

a methylBase,methylRaw, methylRawList or methylDiff object depending on the input object.

Author(s)

Alexander Gosdschan

Examples

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)

# selects the records that lie inside the regions
myRaw2 <- selectByOverlap(methylRawList.obj[[1]],my.win)

# selects the records that lie inside the regions
myRawList2 <- selectByOverlap(methylRawList.obj,my.win)

# selects the records that lie inside the regions
myBase2 <- selectByOverlap(methylBase.obj,my.win)

# selects the records that lie inside the regions
myDiff2 <- selectByOverlap(methylDiff.obj,my.win)


rm(methylRawListDB.obj)
rm(methylBaseDB.obj)
rm(methylDiffDB.obj)
unlink("methylDB",recursive=TRUE)


al2na/methylKit documentation built on Feb. 1, 2024, 4:42 p.m.