Skip to contents

Counts

Usage

counts(object, ...)
counts(object, ...) <- value

# S4 method for bcbioRNASeq
counts(object, normalized = FALSE)

Arguments

object

Object.

normalized

character(1) or logical(1). Normalization method to apply:

  • FALSE: Raw counts. When using a tximport-compatible caller, these are length scaled by default (see countsFromAbundance argument). When using a featureCounts-compatible caller, these are integer.

tximport caller-specific normalizations:

  • "tpm": Transcripts per million.

Additional gene-level-specific normalizations:

  • TRUE / "sf": Size factor (i.e. library size) normalized counts.
    See DESeq2::sizeFactors for details.

  • "fpkm": Fragments per kilobase per million mapped fragments.
    Requires fast = FALSE in bcbioRNASeq() call and gene annotations in rowRanges() with defined width().
    See DESeq2::fpkm() for details.

  • "vst": Variance-stabilizing transformation (log2).
    Requires fast = FALSE to be set during bcbioRNASeq() call.
    See DESeq2::varianceStabilizingTransformation() for more information.

  • "tmm": Trimmed mean of M-values.
    Calculated on the fly.
    See edgeR::calcNormFactors() for details.

  • "rle": Relative log expression transformation.
    Calculated on the fly.
    See relativeLogExpression() for details.

  • "rlog": Deprecated. Regularized log transformation (log2).
    No longer calculated automatically during bcbioRNASeq() call, but may be defined in legacy objects.
    See DESeq2::rlog() for details.
    Note that VST is more performant and now recommended by default instead.

Note that logical(1) support only applies to counts(). Other functions in the package require character(1) and use match.arg() internally.

...

Additional arguments.

value

Value to assign.

Value

matrix.

Details

By default, return the raw counts. Normalized counts in a variety of formats can be accessed using the normalized argument.

Note

Updated 2022-03-07.

References

  • TMM: Robinson and Oshlack (2010).

  • RLE: Anders and Huber (2010).

Author

Michael Steinbaugh, Lorena Pantano

Examples

## bcbioRNASeq ====
data(bcb)
summary(counts(bcb))
#>   control_rep1      control_rep2       control_rep3      fa_day7_rep1    
#>  Min.   :   0.00   Min.   :    0.00   Min.   :   0.00   Min.   :   0.00  
#>  1st Qu.:   0.00   1st Qu.:    0.00   1st Qu.:   0.00   1st Qu.:   0.00  
#>  Median :  12.12   Median :   14.17   Median :  10.72   Median :  33.48  
#>  Mean   : 176.59   Mean   :  602.45   Mean   : 161.05   Mean   : 294.82  
#>  3rd Qu.: 120.48   3rd Qu.:  419.69   3rd Qu.: 102.85   3rd Qu.: 302.74  
#>  Max.   :2865.68   Max.   :13680.79   Max.   :2154.71   Max.   :4227.90  
#>   fa_day7_rep2       fa_day7_rep3     
#>  Min.   :   0.000   Min.   :   0.000  
#>  1st Qu.:   0.761   1st Qu.:   0.005  
#>  Median :  30.355   Median :  27.788  
#>  Mean   : 262.722   Mean   : 298.832  
#>  3rd Qu.: 238.995   3rd Qu.: 215.032  
#>  Max.   :5004.905   Max.   :4540.914