Add gene lengths to an SingleCellExperiment object
addGeneLengths(sce, method = c("generate", "sample"), loc = 7.9, scale = 0.7, lengths = NULL)
sce | SingleCellExperiment to add gene lengths to. |
---|---|
method | Method to use for creating lengths. |
loc | Location parameter for the generate method. |
scale | Scale parameter for the generate method. |
lengths | Vector of lengths for the sample method. |
SingleCellExperiment with added gene lengths
This function adds simulated gene lengths to the
rowData
slot of a
SingleCellExperiment
object that can be
used for calculating length normalised expression values such as TPM or FPKM.
The generate
method simulates lengths using a (rounded) log-normal
distribution, with the default loc
and scale
parameters based
on human protein-coding genes. Alternatively the sample
method can be
used which randomly samples lengths (with replacement) from a supplied
vector.
#>#>#>#> DataFrame with 6 rows and 3 columns #> Gene GeneMean Length #> <factor> <numeric> <numeric> #> Gene1 Gene1 0.0138665247174492 6630 #> Gene2 Gene2 1.25529432999633 4894 #> Gene3 Gene3 0.0375780248572044 1723 #> Gene4 Gene4 1.69560565875577 1886 #> Gene5 Gene5 1.12938171645951 438 #> Gene6 Gene6 0.14592367445476 1508# Sample method (human coding genes)# NOT RUN { library(TxDb.Hsapiens.UCSC.hg19.knownGene) library(GenomicFeatures) txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene tx.lens <- transcriptLengths(txdb, with.cds_len = TRUE) tx.lens <- tx.lens[tx.lens$cds_len > 0, ] gene.lens <- max(splitAsList(tx.lens$tx_len, tx.lens$gene_id)) sce <- addGeneLengths(sce, method = "sample", lengths = gene.lens) # }