quantilePositions: Determine CTSS quantile positions within clusters

Description Usage Arguments Details Value Author(s) See Also Examples

Description

Calculates the positions of “upper” and “lower” quantiles of CAGE signal along tag clusters or consensus clusters in each sample of a CAGEr object.

Usage

1
2
3
4
5
6
7
8
quantilePositions(object, clusters = c("tagClusters",
  "consensusClusters"), qLow = 0.1, qUp = 0.9, useMulticore = FALSE,
  nrCores = NULL)

## S4 method for signature 'CAGEr'
quantilePositions(object, clusters = c("tagClusters",
  "consensusClusters"), qLow = 0.1, qUp = 0.9, useMulticore = FALSE,
  nrCores = NULL)

Arguments

object

A CAGEr object.

clusters

Either tagClusters or consensusClusters.

qLow, qUp

Which “lower” or “upper” quantiles should be calculated. Numeric vector of values in range [0,1].

useMulticore

Logical, should multicore be used. useMulticore = TRUE has only effect on Unix-like platforms.

nrCores

Number of cores to use when useMulticore = TRUE. Default value NULL uses all detected cores.

Details

From the 5' end the position, the position of a quantile q is detemined as the first base in which of the cumulative expression is higher or equal to q% of the total CAGE signal of that cluster. Promoter interquantile width is defined as the distance (in base pairs) between a “lower” and an “upper” quantile position.

Value

When clusters = "tagClusters", the slots tagClustersQuantileLow and tagClustersQuantileUp of a provided CAGEset object will be occupied with the positions of specified quantiles in all tag clusters for all CAGE datasets. When clusters = "consensusClusters" the slots consensusClustersQuantileLow and consensusClustersQuantileUp will be occupied by the corresponding information for consensus clusters.

In CAGEexp objects, the positions of the quantiles are defined reliatively to the start point of their cluster, for more efficient Rle compression. The quantile data for tag clusters are stored in the TagClusters objects directly. The quantile data for consensus clusters are stored in integer matrices named “q_x”, where x represents the quantile (for instance, q_0.1), and these matrices are assays of the consensusClusters RangedSummarizedExperiment.

Author(s)

Vanja Haberle

See Also

Other CAGEr object modifiers: CTSStoGenes, CustomConsensusClusters, aggregateTagClusters, annotateCTSS, clusterCTSS, cumulativeCTSSdistribution, getCTSS, normalizeTagCount, summariseChrExpr

Other CAGEr clusters functions: CTSSclusteringMethod, CTSScumulativesTagClusters, CustomConsensusClusters, aggregateTagClusters, clusterCTSS, consensusClustersDESeq2, consensusClustersGR, cumulativeCTSSdistribution, plotInterquantileWidth, tagClusters

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
head(cbind(
  CAGEr:::tagClustersQuantileLow(exampleCAGEset, 1),
  CAGEr:::tagClustersQuantileUp (exampleCAGEset, 1)
))
quantilePositions( object = exampleCAGEset, clusters = "tagClusters"
                 , qLow = c(0.1, 0.2), qUp = c(0.8, 0.9))
head(cbind(
  CAGEr:::tagClustersQuantileLow(exampleCAGEset, 1),
  CAGEr:::tagClustersQuantileUp (exampleCAGEset,1 )
))

cumulativeCTSSdistribution(exampleCAGEset, "consensusClusters") # Defaults in object do not fit
quantilePositions( object = exampleCAGEset, clusters = "consensusClusters"
                 , qLow = c(0.1, 0.2), qUp = c(0.8, 0.9))
                 
head(cbind(
  CAGEr:::consensusClustersQuantileLow(exampleCAGEset, 1),
  CAGEr:::consensusClustersQuantileUp (exampleCAGEset , 1)
))

quantilePositions(exampleCAGEexp, "tagClusters",       qLow = c(0.1, 0.2), qUp = c(0.8, 0.9))
tagClustersGR(exampleCAGEexp)
quantilePositions(exampleCAGEexp, "consensusClusters", qLow = c(0.1, 0.2), qUp = c(0.8, 0.9))

CAGEr documentation built on Jan. 17, 2021, 2 a.m.