Last data update: 2014.03.03

R: Views into a set of BAM files
BamViewsR Documentation

Views into a set of BAM files

Description

Use BamViews() to reference a set of disk-based BAM files to be processed (e.g., queried using scanBam) as a single ‘experiment’.

Usage


## 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)

Arguments

bamPaths

A character() vector of BAM path names.

bamIndicies

A character() vector of BAM index file path names, without the ‘.bai’ extension.

bamSamples

A DataFrame instance with as many rows as length(bamPaths), containing sample information associated with each path.

bamRanges

Missing or a GRanges instance with ranges defined on the reference chromosomes of the BAM files. Ranges are not validated against the BAM files.

bamExperiment

A list() containing additional information about the experiment.

auto.range

If TRUE and all bamPaths exist, populate the ranges with the union of ranges returned in the target element of scanBamHeader.

...

Additional arguments.

x

An instance of BamViews.

object

An instance of BamViews.

value

An object of appropriate type to replace content.

i

During subsetting, a logical or numeric index into bamRanges.

j

During subsetting, a logical or numeric index into bamSamples and bamPaths.

drop

A logical(1), ignored by all BamViews subsetting methods.

file

An instance of BamViews.

index

A character vector of indices, corresponding to the bamPaths(file).

param

An optional ScanBamParam instance to further influence scanning or counting.

Objects from the Class

Objects are created by calls of the form BamViews().

Slots

bamPaths

A character() vector of BAM path names.

bamIndicies

A character() vector of BAM index path names.

bamSamples

A DataFrame instance with as many rows as length(bamPaths), containing sample information associated with each path.

bamRanges

A GRanges instance with ranges defined on the spaces of the BAM files. Ranges are not validated against the BAM files.

bamExperiment

A list() containing additional information about the experiment.

Functions and methods

See 'Usage' for details on invocation.

Constructor:

BamViews:

Returns a BamViews object.

Accessors:

bamPaths

Returns a character() vector of BAM path names.

bamIndicies

Returns a character() vector of BAM index path names.

bamSamples

Returns a DataFrame instance with as many rows as length(bamPaths), containing sample information associated with each path.

bamSamples<-

Assign a DataFrame instance with as many rows as length(bamPaths), containing sample information associated with each path.

bamRanges

Returns a GRanges instance with ranges defined on the spaces of the BAM files. Ranges are not validated against the BAM files.

bamRanges<-

Assign a GRanges instance with ranges defined on the spaces of the BAM files. Ranges are not validated against the BAM files.

bamExperiment

Returns a list() containing additional information about the experiment.

names

Return the column names of the BamViews instance; same as names(bamSamples(x)).

names<-

Assign the column names of the BamViews instance.

dimnames

Return the row and column names of the BamViews instance.

dimnames<-

Assign the row and column names of the BamViews instance.

Methods:

"["

Subset the object by bamRanges or bamSamples.

scanBam

Visit each path in bamPaths(file), returning the result of scanBam applied to the specified path. bamRanges(file) takes precedence over bamWhich(param).

countBam

Visit each path in bamPaths(file), returning the result of countBam applied to the specified path. bamRanges(file) takes precedence over bamWhich(param).

show

Compactly display the object.

Author(s)

Martin Morgan

Examples

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

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(Rsamtools)
Loading required package: GenomeInfoDb
Loading required package: stats4
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, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, 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, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: S4Vectors

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: GenomicRanges
Loading required package: Biostrings
Loading required package: XVector
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/Rsamtools/BamViews-class.Rd_%03d_medium.png", width=480, height=480)
> ### Name: BamViews
> ### Title: Views into a set of BAM files
> ### Aliases: BamViews-class BamViews BamViews,missing-method
> ###   BamViews,GRanges-method bamPaths bamIndicies bamSamples bamRanges
> ###   bamExperiment bamDirname<- bamSamples<- bamRanges<-
> ###   names,BamViews-method names<-,BamViews-method dim,BamViews-method
> ###   dimnames,BamViews-method dimnames<-,BamViews-method
> ###   dimnames<-,BamViews,ANY-method [,BamViews,ANY,missing-method
> ###   [,BamViews,missing,ANY-method [,BamViews,ANY,ANY-method
> ###   [,BamViews,ANY,missing-method scanBam,BamViews-method
> ###   countBam,BamViews-method show,BamViews-method
> ### Keywords: classes
> 
> ### ** Examples
> 
> 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
BamViews dim: 18 ranges x 1 samples 
names: ex1.bam 
detail: use bamPaths(), bamSamples(), bamRanges(), ... 
> v[1:5,]
BamViews dim: 5 ranges x 1 samples 
names: ex1.bam 
detail: use bamPaths(), bamSamples(), bamRanges(), ... 
> bamRanges(v[c(1:5, 11:15),])
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
> bamDirname(v) <- getwd()
> v
BamViews dim: 18 ranges x 1 samples 
names: ex1.bam 
detail: use bamPaths(), bamSamples(), bamRanges(), ... 
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>