Extract genes by row and samples by column.

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

Arguments

x

object from which to extract element(s) or in which to replace element(s).

i

indices specifying elements to extract or replace. Indices are numeric or character vectors or empty (missing) or NULL. Numeric values are coerced to integer as by as.integer (and hence truncated towards zero). Character vectors will be matched to the names of the object (or for matrices/arrays, the dimnames): see ‘Character indices’ below for further details.

For [-indexing only: i, j, ... can be logical vectors, indicating elements/slices to select. Such vectors are recycled if necessary to match the corresponding extent. i, j, ... can also be negative integers, indicating elements/slices to leave out of the selection.

When indexing arrays by [ a single argument i can be a matrix with as many columns as there are dimensions of x; the result is then a vector with elements corresponding to the sets of indices in each row of i.

An index value of NULL is treated as if it were integer(0).

j

indices specifying elements to extract or replace. Indices are numeric or character vectors or empty (missing) or NULL. Numeric values are coerced to integer as by as.integer (and hence truncated towards zero). Character vectors will be matched to the names of the object (or for matrices/arrays, the dimnames): see ‘Character indices’ below for further details.

For [-indexing only: i, j, ... can be logical vectors, indicating elements/slices to select. Such vectors are recycled if necessary to match the corresponding extent. i, j, ... can also be negative integers, indicating elements/slices to leave out of the selection.

When indexing arrays by [ a single argument i can be a matrix with as many columns as there are dimensions of x; the result is then a vector with elements corresponding to the sets of indices in each row of i.

An index value of NULL is treated as if it were integer(0).

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.

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.

Value

bcbioRNASeq.

Details

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

Note

Updated 2019-08-20.

References

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

Examples

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.25.17.
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> Applying variance-stabilizing transformation.
#> bcbioRNASeq 0.3.27 #> uploadDir: /data00/draco/acidbase/packages/bcbioRNASeq/inst/extdata/bcbio #> dates(2): [bcbio] 2018-03-18; [R] 2019-09-16 #> level: genes #> caller: salmon #> organism: Mus musculus #> 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(8): broadClass description ... geneName seqCoordSystem #> colnames(2): control_rep1 control_rep2 #> colData names(26): averageInsertSize averageReadLength ... treatment #> x5x3Bias
## Extract by gene list. object[genes, ]
#> bcbioRNASeq 0.3.27 #> uploadDir: /data00/draco/acidbase/packages/bcbioRNASeq/inst/extdata/bcbio #> dates(2): [bcbio] 2018-03-18; [R] 2019-09-16 #> level: genes #> caller: salmon #> organism: Mus musculus #> 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(8): broadClass description ... geneName 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.25.17.
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> Applying variance-stabilizing transformation.
#> bcbioRNASeq 0.3.27 #> uploadDir: /data00/draco/acidbase/packages/bcbioRNASeq/inst/extdata/bcbio #> dates(2): [bcbio] 2018-03-18; [R] 2019-09-16 #> level: genes #> caller: salmon #> organism: Mus musculus #> 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(8): broadClass description ... geneName seqCoordSystem #> colnames(2): control_rep1 control_rep2 #> colData names(26): averageInsertSize averageReadLength ... treatment #> x5x3Bias
assayNames(x)
#> [1] "counts" "aligned" "avgTxLength" "tpm" "normalized" #> [6] "vst" "fpkm"
## Fast subsetting, by skipping DESeq2 recalculations. ## Note that `normalized`, `rlog`, and `vst` assays will be removed. x <- object[, samples, recalculate = FALSE]
#> Skipping DESeq2 normalizations.
#> bcbioRNASeq 0.3.27 #> uploadDir: /data00/draco/acidbase/packages/bcbioRNASeq/inst/extdata/bcbio #> dates(2): [bcbio] 2018-03-18; [R] 2019-09-16 #> level: genes #> caller: salmon #> organism: Mus musculus #> 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(8): broadClass description ... geneName seqCoordSystem #> colnames(2): control_rep1 control_rep2 #> colData names(26): averageInsertSize averageReadLength ... treatment #> x5x3Bias
names(assays(x))
#> [1] "counts" "aligned" "avgTxLength" "tpm"