Description Usage Arguments Details Value See Also Examples
View source: R/transform_counts.R
In preparation for a differential expression analysis, you will have to
choose how to scale the raw counts provided by the recount3 project. These
raw counts are similar to those provided by the recount2 project, except
that they were generated with a different aligner and a modified counting
approach. The raw coverage counts for recount2 are described with
illustrative figures at https://doi.org/10.12688/f1000research.12223.1.
Note that the raw counts are the sum of the base level coverage so you have
to take into account the total base-pair coverage for the given sample
(default option) by using the area under the coverage (AUC), or alternatively
use the mapped read lengths. You might want to do some further scaling to
take into account the gene or exon lengths. If you prefer to calculate read
counts without scaling check the function compute_read_counts()
.
1 2 3 4 5 6 7 8 |
rse |
A
RangedSummarizedExperiment-class
created by |
by |
Either |
targetSize |
A |
L |
A |
round |
A |
... |
Further arguments passed to |
This function is similar to
recount::scale_counts()
but more general and with a different name to
avoid NAMESPACE conflicts.
To compute RPKMs, use recount::getRPKM(length_var = "score")
. Similarly,
for TPMs, use recount::getTPM(length_var = "score")
.
A matrix()
with the transformed (scaled) counts.
Other count transformation functions:
compute_read_counts()
,
compute_scale_factors()
,
is_paired_end()
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 | ## Create a RSE object at the gene level
rse_gene_SRP009615 <- create_rse_manual("SRP009615")
## Scale the counts using the AUC
assays(rse_gene_SRP009615)$counts <- transform_counts(rse_gene_SRP009615)
## See that now we have two assayNames()
rse_gene_SRP009615
assayNames(rse_gene_SRP009615)
## You can compare the scaled counts against those from
## recount::scale_counts() from the recount2 project
## which used a different RNA-seq aligner
## If needed, install recount, the R/Bioconductor package for recount2:
# BiocManager::install("recount")
recount2_sizes <- colSums(assay(recount::scale_counts(
recount::rse_gene_SRP009615,
by = "auc"
), "counts")) / 1e6
recount3_sizes <- colSums(assay(rse_gene_SRP009615, "counts")) / 1e6
recount_sizes <- data.frame(
recount2 = recount2_sizes[order(names(recount2_sizes))],
recount3 = recount3_sizes[order(names(recount3_sizes))]
)
plot(recount2 ~ recount3, data = recount_sizes)
abline(a = 0, b = 1, col = "purple", lwd = 2, lty = 2)
## Compute RPKMs
assays(rse_gene_SRP009615)$RPKM <- recount::getRPKM(rse_gene_SRP009615, length_var = "score")
colSums(assay(rse_gene_SRP009615, "RPKM"))
## Compute TPMs
assays(rse_gene_SRP009615)$TPM <- recount::getTPM(rse_gene_SRP009615, length_var = "score")
colSums(assay(rse_gene_SRP009615, "TPM")) / 1e6 ## Should all be equal to 1
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.