Description Usage Arguments Objects from the Class Slots Functions and methods Author(s) Examples
Use BamViews()
to reference a set of disk-based BAM files to be
processed (e.g., queried using scanBam
) as a single
‘experiment’.
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 | ## Constructor
BamViews(bamPaths=character(0),
bamIndicies=bamPaths,
bamSamples=DataFrame(row.names=make.unique(basename(bamPaths))),
bamRanges, bamExperiment = list(), ...)
## S4 method for signature 'missing'
BamViews(bamPaths=character(0),
bamIndicies=bamPaths,
bamSamples=DataFrame(row.names=make.unique(basename(bamPaths))),
bamRanges, bamExperiment = list(), ..., auto.range=FALSE)
## Accessors
bamPaths(x)
bamSamples(x)
bamSamples(x) <- value
bamRanges(x)
bamRanges(x) <- value
bamExperiment(x)
## S4 method for signature 'BamViews'
names(x)
## S4 replacement method for signature 'BamViews'
names(x) <- value
## S4 method for signature 'BamViews'
dimnames(x)
## S4 replacement method for signature 'BamViews,ANY'
dimnames(x) <- value
bamDirname(x, ...) <- value
## Subset
## S4 method for signature 'BamViews,ANY,ANY'
x[i, j, ..., drop=TRUE]
## S4 method for signature 'BamViews,ANY,missing'
x[i, j, ..., drop=TRUE]
## S4 method for signature 'BamViews,missing,ANY'
x[i, j, ..., drop=TRUE]
## Input
## S4 method for signature 'BamViews'
scanBam(file, index = file, ..., param = ScanBamParam(what=scanBamWhat()))
## S4 method for signature 'BamViews'
countBam(file, index = file, ..., param = ScanBamParam())
## Show
## S4 method for signature 'BamViews'
show(object)
|
bamPaths |
A character() vector of BAM path names. |
bamIndicies |
A character() vector of BAM index file path names, without the ‘.bai’ extension. |
bamSamples |
A |
bamRanges |
Missing or a |
bamExperiment |
A list() containing additional information about the experiment. |
auto.range |
If |
... |
Additional arguments. |
x |
An instance of |
object |
An instance of |
value |
An object of appropriate type to replace content. |
i |
During subsetting, a logical or numeric index into
|
j |
During subsetting, a logical or numeric index into
|
drop |
A logical(1), ignored by all |
file |
An instance of |
index |
A character vector of indices, corresponding to the
|
param |
An optional |
Objects are created by calls of the form BamViews()
.
A character() vector of BAM path names.
A character() vector of BAM index path names.
A DataFrame
instance with as
many rows as length(bamPaths)
, containing sample information
associated with each path.
A GRanges
instance with
ranges defined on the spaces of the BAM files. Ranges are not
validated against the BAM files.
A list() containing additional information about the experiment.
See 'Usage' for details on invocation.
Constructor:
Returns a BamViews
object.
Accessors:
Returns a character() vector of BAM path names.
Returns a character() vector of BAM index path names.
Returns a DataFrame
instance
with as many rows as length(bamPaths)
, containing sample
information associated with each path.
Assign a DataFrame
instance
with as many rows as length(bamPaths)
, containing sample
information associated with each path.
Returns a GRanges
instance
with ranges defined on the spaces of the BAM files. Ranges are
not validated against the BAM files.
Assign a GRanges
instance
with ranges defined on the spaces of the BAM files. Ranges are
not validated against the BAM files.
Returns a list() containing additional information about the experiment.
Return the column names of the BamViews
instance; same as names(bamSamples(x))
.
Assign the column names of the BamViews
instance.
Return the row and column names of the
BamViews
instance.
Assign the row and column names of the
BamViews
instance.
Methods:
Subset the object by bamRanges
or bamSamples
.
Visit each path in bamPaths(file)
, returning
the result of scanBam
applied to the specified
path. bamRanges(file)
takes precedence over
bamWhich(param)
.
Visit each path in bamPaths(file)
, returning
the result of countBam
applied to the specified
path. bamRanges(file)
takes precedence over
bamWhich(param)
.
Compactly display the object.
Martin Morgan
1 2 3 4 5 6 7 8 9 10 11 12 | fls <- system.file("extdata", "ex1.bam", package="Rsamtools",
mustWork=TRUE)
rngs <- GRanges(seqnames = Rle(c("chr1", "chr2"), c(9, 9)),
ranges = c(IRanges(seq(10000, 90000, 10000), width=500),
IRanges(seq(100000, 900000, 100000), width=5000)),
Count = seq_len(18L))
v <- BamViews(fls, bamRanges=rngs)
v
v[1:5,]
bamRanges(v[c(1:5, 11:15),])
bamDirname(v) <- getwd()
v
|
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’:
anyDuplicated, append, as.data.frame, basename, cbind, colnames,
dirname, do.call, duplicated, eval, evalq, Filter, Find, get, grep,
grepl, intersect, is.unsorted, lapply, Map, mapply, match, mget,
order, paste, pmax, pmax.int, pmin, pmin.int, Position, rank,
rbind, Reduce, rownames, sapply, setdiff, sort, table, tapply,
union, unique, unsplit, 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
BamViews dim: 18 ranges x 1 samples
names: ex1.bam
detail: use bamPaths(), bamSamples(), bamRanges(), ...
BamViews dim: 5 ranges x 1 samples
names: ex1.bam
detail: use bamPaths(), bamSamples(), bamRanges(), ...
GRanges object with 10 ranges and 1 metadata column:
seqnames ranges strand | Count
<Rle> <IRanges> <Rle> | <integer>
[1] chr1 10000-10499 * | 1
[2] chr1 20000-20499 * | 2
[3] chr1 30000-30499 * | 3
[4] chr1 40000-40499 * | 4
[5] chr1 50000-50499 * | 5
[6] chr2 200000-204999 * | 11
[7] chr2 300000-304999 * | 12
[8] chr2 400000-404999 * | 13
[9] chr2 500000-504999 * | 14
[10] chr2 600000-604999 * | 15
-------
seqinfo: 2 sequences from an unspecified genome; no seqlengths
BamViews dim: 18 ranges x 1 samples
names: ex1.bam
detail: use bamPaths(), bamSamples(), bamRanges(), ...
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.