View source: R/transform_counts.R
compute_scale_factors | R Documentation |
This function computes the count scaling factors used by
transform_counts()
. This function is similar to
recount::scale_counts(factor_only = TRUE)
, but it is more general.
compute_scale_factors(
x,
by = c("auc", "mapped_reads"),
targetSize = 4e+07,
L = 100,
auc = "recount_qc.bc_auc.all_reads_all_bases",
avg_mapped_read_length = "recount_qc.star.average_mapped_length",
mapped_reads = "recount_qc.star.all_mapped_reads",
paired_end = is_paired_end(x, avg_mapped_read_length)
)
x |
Either a
RangedSummarizedExperiment-class
created by |
by |
Either |
targetSize |
A |
L |
A |
auc |
A |
avg_mapped_read_length |
A |
mapped_reads |
A |
paired_end |
A |
A numeric()
with the sample scale factors that are used by
transform_counts()
.
Other count transformation functions:
compute_read_counts()
,
is_paired_end()
,
transform_counts()
## Download the metadata for SRP009615, a single-end study
SRP009615_meta <- read_metadata(
metadata_files = file_retrieve(
locate_url(
"SRP009615",
"data_sources/sra",
)
)
)
## Compute the scaling factors
compute_scale_factors(SRP009615_meta, by = "auc")
compute_scale_factors(SRP009615_meta, by = "mapped_reads")
## Download the metadata for DRP000499, a paired-end study
DRP000499_meta <- read_metadata(
metadata_files = file_retrieve(
locate_url(
"DRP000499",
"data_sources/sra",
)
)
)
## Compute the scaling factors
compute_scale_factors(DRP000499_meta, by = "auc")
compute_scale_factors(DRP000499_meta, by = "mapped_reads")
## You can compare the factors 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_factors <- recount::scale_counts(
recount::rse_gene_SRP009615,
by = "auc", factor_only = TRUE
)
recount3_factors <- compute_scale_factors(SRP009615_meta, by = "auc")
recount_factors <- data.frame(
recount2 = recount2_factors[order(names(recount2_factors))],
recount3 = recount3_factors[order(names(recount3_factors))]
)
plot(recount2 ~ recount3, data = recount_factors)
abline(a = 0, b = 1, col = "purple", lwd = 2, lty = 2)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.