ApplyPileupsParam-class: Parameters for creating pileups from BAM files

Description Usage Arguments Objects from the Class Slots Functions and methods Author(s) See Also Examples

Description

Use ApplyPileupsParam() to create a parameter object influencing what fields and which records are used to calculate pile-ups, and to influence the values returned.

Usage

 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
# Constructor
ApplyPileupsParam(flag = scanBamFlag(),
    minBaseQuality = 13L, minMapQuality = 0L,
    minDepth = 0L, maxDepth = 250L,
    yieldSize = 1L, yieldBy = c("range", "position"), yieldAll = FALSE,
    which = GRanges(), what = c("seq", "qual"))

# Accessors
plpFlag(object)
plpFlag(object) <- value
plpMaxDepth(object)
plpMaxDepth(object) <- value
plpMinBaseQuality(object)
plpMinBaseQuality(object) <- value
plpMinDepth(object)
plpMinDepth(object) <- value
plpMinMapQuality(object)
plpMinMapQuality(object) <- value
plpWhat(object)
plpWhat(object) <- value
plpWhich(object)
plpWhich(object) <- value
plpYieldAll(object)
plpYieldAll(object) <- value
plpYieldBy(object)
plpYieldBy(object) <- value
plpYieldSize(object)
plpYieldSize(object) <- value

## S4 method for signature 'ApplyPileupsParam'
show(object)

Arguments

flag

An instance of the object returned by scanBamFlag, restricting various aspects of reads to be included or excluded.

minBaseQuality

The minimum read base quality below which the base is ignored when summarizing pileup information.

minMapQuality

The minimum mapping quality below which the entire read is ignored.

minDepth

The minimum depth of the pile-up below which the position is ignored.

maxDepth

The maximum depth of reads considered at any position; this can be used to limit memory consumption.

yieldSize

The number of records to include in each call to FUN.

yieldBy

How records are to be counted. By range (in which case yieldSize must equal 1) means that FUN is invoked once for each range in which. By position means that FUN is invoked whenever pile-ups have been accumulated for yieldSize positions, regardless of ranges in which.

yieldAll

Whether to report all positions (yieldAll=TRUE), or just those passing the filtering criteria of flag, minBaseQuality, etc. When yieldAll=TRUE, positions not passing filter criteria have ‘0’ entries in seq or qual.

which

A GRanges or IntegerRangesList instance restricting pileup calculations to the corresponding genomic locations.

what

A character() instance indicating what values are to be returned. One or more of c("seq", "qual").

object

An instace of class ApplyPileupsParam.

value

An instance to be assigned to the corresponding slot of the ApplyPileupsParam instance.

Objects from the Class

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

Slots

Slot interpretation is as described in the ‘Arguments’ section.

flag

Object of class integer encoding flags to be kept when they have their '0' (keep0) or '1' (keep1) bit set.

minBaseQuality

An integer(1).

minMapQuality

An integer(1).

minDepth

An integer(1).

maxDepth

An integer(1).

yieldSize

An integer(1).

yieldBy

An character(1).

yieldAll

A logical(1).

which

A GRanges or IntegerRangesList object.

what

A character().

Functions and methods

See 'Usage' for details on invocation.

Constructor:

ApplyPileupsParam:

Returns a ApplyPileupsParam object.

Accessors: get or set corresponding slot values; for setters, value is coerced to the type of the corresponding slot.

plpFlag, plpFlag<-

Returns or sets the named integer vector of flags; see scanBamFlag.

plpMinBaseQuality, plpMinBaseQuality<-

Returns or sets an integer(1) vector of miminum base qualities.

plpMinMapQuality, plpMinMapQuality<-

Returns or sets an integer(1) vector of miminum map qualities.

plpMinDepth, plpMinDepth<-

Returns or sets an integer(1) vector of miminum pileup depth.

plpMaxDepth, plpMaxDepth<-

Returns or sets an integer(1) vector of the maximum depth to which pileups are calculated.

plpYieldSize, plpYieldSize<-

Returns or sets an integer(1) vector of yield size.

plpYieldBy, plpYieldBy<-

Returns or sets an character(1) vector determining how pileups will be returned.

plpYieldAll, plpYieldAll<-

Returns or sets an logical(1) vector indicating whether all positions, or just those satisfying pileup positions, are to be returned.

plpWhich, plpWhich<-

Returns or sets the object influencing which locations pileups are calculated over.

plpWhat, plpWhat<-

Returns or sets the character vector describing what summaries are returned by pileup.

Methods:

show

Compactly display the object.

Author(s)

Martin Morgan

See Also

applyPileups.

Examples

1

Example output

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 

Rsamtools documentation built on Nov. 8, 2020, 8:11 p.m.