Description Usage Arguments Objects from the Class Fields Functions and methods Author(s) Examples
Use PileupFiles() to create a reference to BAM files (and
their indicies), to be used for calculating pile-up summaries.
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 | ## Constructors
PileupFiles(files, ..., param=ApplyPileupsParam())
## S4 method for signature 'character'
PileupFiles(files, ..., param=ApplyPileupsParam())
## S4 method for signature 'list'
PileupFiles(files, ..., param=ApplyPileupsParam())
## opening / closing
## open(con, ...)
## close(con, ...)
## accessors; also path()
## S4 method for signature 'PileupFiles'
isOpen(con, rw="")
plpFiles(object)
plpParam(object)
## actions
## S4 method for signature 'PileupFiles,missing'
applyPileups(files, FUN, ..., param)
## S4 method for signature 'PileupFiles,ApplyPileupsParam'
applyPileups(files, FUN, ..., param)
## display
## S4 method for signature 'PileupFiles'
show(object)
|
files |
For For |
... |
Additional arguments, currently ignored. |
con, object |
An instance of |
FUN |
A function of one argument; see |
param |
An instance of |
rw |
character() indicating mode of file; not used for
|
Objects are created by calls of the form PileupFiles().
The PileupFiles class is implemented as an S4 reference
class. It has the following fields:
A list of BamFile instances.
An instance of ApplyPileupsParam.
Opening / closing:
Opens the (local or remote) path and
index of each file in the PileupFiles
instance. Returns a PileupFiles instance.
Closes each file in the PileupFiles
instance; returning (invisibly) the updated
PileupFiles. The instance may be re-opened with
open.PileupFiles.
Accessors:
Returns the list of the files in the
PileupFiles instance.
Returns the ApplyPileupsParam content of the
PileupFiles instance.
Methods:
Calculate the pileup across all files in
files according to criteria in param (or
plpParam(files) if param is missing), invoking
FUN on each range or collection of positions. See
applyPileups.
Compactly display the object.
Martin Morgan
1 |
Loading required package: GenomeInfoDb
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, sd, var, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, basename, cbind, colMeans, colSums, colnames,
dirname, 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,
rowMeans, rowSums, rownames, sapply, setdiff, sort, table, tapply,
union, unique, unsplit, which, which.max, which.min
Loading required package: S4Vectors
Loading required package: stats4
Attaching package: 'S4Vectors'
The following object is masked from 'package:base':
expand.grid
Loading required package: IRanges
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector
Attaching package: 'Biostrings'
The following object is masked from 'package:base':
strsplit
Attaching package: 'Rsamtools'
The following object is masked from 'package:BiocGenerics':
path
applyP> fl <- system.file("extdata", "ex1.bam", package="Rsamtools",
applyP+ mustWork=TRUE)
applyP> fls <- PileupFiles(c(fl, fl))
applyP> calcInfo <-
applyP+ function(x)
applyP+ {
applyP+ ## information at each pile-up position
applyP+ info <- apply(x[["seq"]], 2, function(y) {
applyP+ y <- y[c("A", "C", "G", "T"),,drop=FALSE]
applyP+ y <- y + 1L # continuity
applyP+ cvg <- colSums(y)
applyP+ p <- y / cvg[col(y)]
applyP+ h <- -colSums(p * log(p))
applyP+ ifelse(cvg == 4L, NA, h)
applyP+ })
applyP+ list(seqnames=x[["seqnames"]], pos=x[["pos"]], info=info)
applyP+ }
applyP> which <- GRanges(c("seq1", "seq2"), IRanges(c(1000, 1000), 2000))
applyP> param <- ApplyPileupsParam(which=which, what="seq")
applyP> res <- applyPileups(fls, calcInfo, param=param)
applyP> str(res)
List of 2
$ :List of 3
..$ seqnames: Named int 570
.. ..- attr(*, "names")= chr "seq1"
..$ pos : int [1:570] 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 ...
..$ info : num [1:570, 1:2] 0.356 0.356 0.336 0.336 0.33 ...
.. ..- attr(*, "dimnames")=List of 2
.. .. ..$ : NULL
.. .. ..$ : chr [1:2] "ex1.bam" "ex1.bam"
$ :List of 3
..$ seqnames: Named int 568
.. ..- attr(*, "names")= chr "seq2"
..$ pos : int [1:568] 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 ...
..$ info : num [1:568, 1:2] 0.238 0.235 0.295 0.241 0.235 ...
.. ..- attr(*, "dimnames")=List of 2
.. .. ..$ : NULL
.. .. ..$ : chr [1:2] "ex1.bam" "ex1.bam"
applyP> head(res[[1]][["pos"]]) # positions matching param
[1] 1000 1001 1002 1003 1004 1005
applyP> head(res[[1]][["info"]]) # inforamtion in each file
ex1.bam ex1.bam
[1,] 0.3556980 0.3556980
[2,] 0.3556980 0.3556980
[3,] 0.3357909 0.3357909
[4,] 0.3357909 0.3357909
[5,] 0.3296843 0.3296843
[6,] 0.3296843 0.3296843
applyP> ## 'param' as part of 'files'
applyP> fls1 <- PileupFiles(c(fl, fl), param=param)
applyP> res1 <- applyPileups(fls1, calcInfo)
applyP> identical(res, res1)
[1] TRUE
applyP> ## yield by position, across ranges
applyP> param <- ApplyPileupsParam(which=which, yieldSize=500L,
applyP+ yieldBy="position", what="seq")
applyP> res <- applyPileups(fls, calcInfo, param=param)
applyP> sapply(res, "[[", "seqnames")
[[1]]
seq1
500
[[2]]
seq1 seq2
70 430
[[3]]
seq2
138
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.