Skip to contents

Extract genes by row and samples by column.

Usage

# S4 method for bcbioRNASeq,ANY,ANY,ANY
[(x, i, j, recalculate = TRUE, drop = FALSE)

Arguments

x

Object.

i

Indices specifying elements to extract or replace. Indices are numeric or character vectors, empty (missing), or NULL.

For more information:

help(topic = "Extract", package = "base")

j

Indices specifying elements to extract or replace. Indices are numeric or character vectors, empty (missing), or NULL.

For more information:

help(topic = "Extract", package = "base")

recalculate

logical(1). Recalculate DESeq2 normalized counts and variance-stabilizing transformations defined in assays.
Recommended by default, but can take a long time for large datasets.
\ If FALSE, these assays will be removed automatically: normalized, rlog, vst.

drop

For matrices and arrays. If TRUE the result is coerced to the lowest possible dimension (see the examples). This only works for extracting elements, not for the replacement. See drop for further details.

Value

bcbioRNASeq.

Details

Internal count transformations are rescaled automatically, if defined. DESeq2 transformations will only be updated when recalculate = TRUE.

Note

Updated 2024-03-27.

References

Becker, R. A., Chambers, J. M. and Wilks, A. R. (1988) The New S Language. Wadsworth & Brooks/Cole.

Author

Michael Steinbaugh, Lorena Pantano

Examples

## bcbioRNASeq ====
data(bcb)
object <- bcb

## Minimum of 100 genes, 2 samples.
genes <- head(rownames(object), 100L)
head(genes)
#> [1] "ENSMUSG00000000001" "ENSMUSG00000000003" "ENSMUSG00000000028"
#> [4] "ENSMUSG00000000049" "ENSMUSG00000000058" "ENSMUSG00000000078"
samples <- head(colnames(object), 2L)
head(samples)
#> [1] "control_rep1" "control_rep2"

## Extract by sample name.
object[, samples]
#> → Recalculating DESeq2 normalizations.
#> → Generating <DESeqDataSet> with DESeq2 1.42.1.
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> → `varianceStabilizingTransformation()`
#> bcbioRNASeq 0.6.0 with 100 rows and 2 columns
#> uploadDir: /Users/mike/git/monorepo/r-packages/bcbioRNASeq/inst/extdata/bcbio
#> dates(2): 2018-03-18 2023-10-05
#> level: genes
#> caller: salmon
#> organism: Mus musculus
#> ensemblRelease: 90
#> interestingGroups(2): treatment day
#> class: RangedSummarizedExperiment 
#> dim: 100 2 
#> metadata(28): allSamples bcbioCommandsLog ... yaml subset
#> assays(7): counts aligned ... vst fpkm
#> rownames(100): ENSMUSG00000000001 ENSMUSG00000000003 ...
#>   ENSMUSG00000062661 ENSMUSG00000074340
#> rowData names(9): broadClass description ... ncbiGeneId seqCoordSystem
#> colnames(2): control_rep1 control_rep2
#> colData names(26): averageInsertSize averageReadLength ... treatment
#>   x5x3Bias

## Extract by gene list.
object[genes, ]
#> bcbioRNASeq 0.6.0 with 100 rows and 6 columns
#> uploadDir: /Users/mike/git/monorepo/r-packages/bcbioRNASeq/inst/extdata/bcbio
#> dates(2): 2018-03-18 2023-10-05
#> level: genes
#> caller: salmon
#> organism: Mus musculus
#> ensemblRelease: 90
#> interestingGroups(2): treatment day
#> class: RangedSummarizedExperiment 
#> dim: 100 6 
#> metadata(27): allSamples bcbioCommandsLog ... wd yaml
#> assays(7): counts aligned ... vst fpkm
#> rownames(100): ENSMUSG00000000001 ENSMUSG00000000003 ...
#>   ENSMUSG00000062661 ENSMUSG00000074340
#> rowData names(9): broadClass description ... ncbiGeneId seqCoordSystem
#> colnames(6): control_rep1 control_rep2 ... fa_day7_rep2 fa_day7_rep3
#> colData names(26): averageInsertSize averageReadLength ... treatment
#>   x5x3Bias

## Extract by both genes and samples.
x <- object[genes, samples]
#> → Recalculating DESeq2 normalizations.
#> → Generating <DESeqDataSet> with DESeq2 1.42.1.
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> → `varianceStabilizingTransformation()`
print(x)
#> bcbioRNASeq 0.6.0 with 100 rows and 2 columns
#> uploadDir: /Users/mike/git/monorepo/r-packages/bcbioRNASeq/inst/extdata/bcbio
#> dates(2): 2018-03-18 2023-10-05
#> level: genes
#> caller: salmon
#> organism: Mus musculus
#> ensemblRelease: 90
#> interestingGroups(2): treatment day
#> class: RangedSummarizedExperiment 
#> dim: 100 2 
#> metadata(28): allSamples bcbioCommandsLog ... yaml subset
#> assays(7): counts aligned ... vst fpkm
#> rownames(100): ENSMUSG00000000001 ENSMUSG00000000003 ...
#>   ENSMUSG00000062661 ENSMUSG00000074340
#> rowData names(9): broadClass description ... ncbiGeneId seqCoordSystem
#> colnames(2): control_rep1 control_rep2
#> colData names(26): averageInsertSize averageReadLength ... treatment
#>   x5x3Bias

## Fast subsetting, by skipping DESeq2 recalculations.
## Note that `normalized`, `rlog`, and `vst` assays will be removed.
x <- object[, samples, recalculate = FALSE]
#> ! Skipping `DESeq2()` normalizations.
print(x)
#> bcbioRNASeq 0.6.0 with 100 rows and 2 columns
#> uploadDir: /Users/mike/git/monorepo/r-packages/bcbioRNASeq/inst/extdata/bcbio
#> dates(2): 2018-03-18 2023-10-05
#> level: genes
#> caller: salmon
#> organism: Mus musculus
#> ensemblRelease: 90
#> interestingGroups(2): treatment day
#> class: RangedSummarizedExperiment 
#> dim: 100 2 
#> metadata(28): allSamples bcbioCommandsLog ... yaml subset
#> assays(4): counts aligned avgTxLength tpm
#> rownames(100): ENSMUSG00000000001 ENSMUSG00000000003 ...
#>   ENSMUSG00000062661 ENSMUSG00000074340
#> rowData names(9): broadClass description ... ncbiGeneId seqCoordSystem
#> colnames(2): control_rep1 control_rep2
#> colData names(26): averageInsertSize averageReadLength ... treatment
#>   x5x3Bias