Force an object to belong to a class.

Note

Updated 2019-08-07.

bcbioRNASeq to DESeqDataSet

  1. Coerce to RangedSummarizedExperiment.

  2. Round raw counts to integer matrix.

  3. Subset colData() to include only clean factor columns. See sampleData() for details.

  4. Simplify metadata() to include only relevant information and updates sessionInfo.

Note that gene-level counts are required. Alternatively, summarizeToGene() can be called to convert transcript-level counts to gene-level. By default, we're using length-scaled TPM, so a corresponding average transcript length matrix isn't necessary. The average transcript length matrix is only necessary when raw counts matrix isn't scaled during tximport call (see countsFromAbundance in tximport() documentation).

bcbioRNASeq to DESeqTransform

  1. Coerce to DESeqDataSet.

  2. Call DESeq2::DESeq().

  3. Call DESeq2::varianceStabilizingTransformation().

bcbioRNASeq to DGEList

  1. Obtain per-observation scaling factors for length, adjusted to avoid changing the magnitude of the counts.

  2. Computing effective library sizes from scaled counts, to account for composition biases between samples.

  3. Combine effective library sizes with the length factors, and calculate offsets for a log-link GLM.

  4. Call edgeR::DGEList()

  5. Apply offset matrix using edgeR::scaleOffset().

See also

Examples

data(bcb) ## bcbioRNASeq to DESeqDataSet ==== x <- as(bcb, "DESeqDataSet")
#> Generating DESeqDataSet with DESeq2 1.25.17.
names(S4Vectors::mcols(x))
#> [1] "broadClass" "description" "entrezID" "geneBiotype" #> [5] "geneID" "geneIDVersion" "geneName" "seqCoordSystem"
#> [1] "DESeqDataSet" #> attr(,"package") #> [1] "DESeq2"
show(x)
#> class: DESeqDataSet #> dim: 100 6 #> metadata(11): caller countsFromAbundance ... date sessionInfo #> assays(1): counts #> 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
## bcbioRNASeq to RangedSummarizedExperiment ==== x <- as(bcb, "RangedSummarizedExperiment") slotNames(x)
#> [1] "rowRanges" "colData" "assays" "NAMES" #> [5] "elementMetadata" "metadata"
show(x)
#> 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
## bcbioRNASeq to SummarizedExperiment ==== ## Coerce to RangedSummarizedExperiment first. x <- as(bcb, "RangedSummarizedExperiment") x <- as(x, "SummarizedExperiment") class(x)
#> [1] "SummarizedExperiment" #> attr(,"package") #> [1] "SummarizedExperiment"
slotNames(x)
#> [1] "colData" "assays" "NAMES" "elementMetadata" #> [5] "metadata"
show(x)
#> class: SummarizedExperiment #> 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