create_rse | R Documentation |
Once you have identified a project you want to work with, you can use this
function to construct a recount3
RangedSummarizedExperiment-class
(RSE) object at the gene or exon expression feature level. This function will
retrieve the data, cache it, then assemble the RSE object.
create_rse(
project_info,
type = c("gene", "exon", "jxn"),
annotation = annotation_options(project_info$organism),
bfc = recount3_cache(),
jxn_format = c("ALL", "UNIQUE"),
recount3_url = getOption("recount3_url", "http://duffel.rail.bio/recount3"),
verbose = getOption("recount3_verbose", TRUE)
)
project_info |
A |
type |
A |
annotation |
A |
bfc |
A BiocFileCache-class
object where the files will be cached to, typically created by
|
jxn_format |
A |
recount3_url |
A |
verbose |
A |
A RangedSummarizedExperiment-class object.
## Find all available human projects
human_projects <- available_projects()
## Find the project you are interested in
proj_info <- subset(
human_projects,
project == "SRP009615" & project_type == "data_sources"
)
## Create a RSE object at the gene level
rse_gene_SRP009615 <- create_rse(proj_info)
## Explore the resulting RSE gene object
rse_gene_SRP009615
## Information about how this RSE object was made
metadata(rse_gene_SRP009615)
## Number of genes by number of samples
dim(rse_gene_SRP009615)
## Information about the genes
rowRanges(rse_gene_SRP009615)
## Sample metadata
colnames(colData(rse_gene_SRP009615))
## Check how much memory this RSE object uses
pryr::object_size(rse_gene_SRP009615)
## Create an RSE object using gencode_v29 instead of gencode_v26
rse_gene_SRP009615_gencode_v29 <- create_rse(
proj_info,
annotation = "gencode_v29",
verbose = FALSE
)
rowRanges(rse_gene_SRP009615_gencode_v29)
## Create an RSE object using FANTOM6_CAT instead of gencode_v26
rse_gene_SRP009615_fantom6_cat <- create_rse(
proj_info,
annotation = "fantom6_cat"
)
rowRanges(rse_gene_SRP009615_fantom6_cat)
## Create an RSE object using RefSeq instead of gencode_v26
rse_gene_SRP009615_refseq <- create_rse(
proj_info,
annotation = "refseq"
)
rowRanges(rse_gene_SRP009615_refseq)
## Create an RSE object using ERCC instead of gencode_v26
rse_gene_SRP009615_ercc <- create_rse(
proj_info,
annotation = "ercc"
)
rowRanges(rse_gene_SRP009615_ercc)
## Create an RSE object using SIRV instead of gencode_v26
rse_gene_SRP009615_sirv <- create_rse(
proj_info,
annotation = "sirv"
)
rowRanges(rse_gene_SRP009615_sirv)
## Obtain a list of RSE objects for all gene annotations
rses_gene <- lapply(annotation_options(), function(x) {
create_rse(proj_info, type = "gene", annotation = x)
})
names(rses_gene) <- annotation_options()
rses_gene
## Create a RSE object at the exon level
rse_exon_SRP009615 <- create_rse(
proj_info,
type = "exon"
)
## Explore the resulting RSE exon object
rse_exon_SRP009615
dim(rse_exon_SRP009615)
rowRanges(rse_exon_SRP009615)
pryr::object_size(rse_exon_SRP009615)
## Create a RSE object at the exon-exon junction level
rse_jxn_SRP009615 <- create_rse(
proj_info,
type = "jxn"
)
## Explore the resulting RSE exon-exon junctions object
rse_jxn_SRP009615
dim(rse_jxn_SRP009615)
rowRanges(rse_jxn_SRP009615)
pryr::object_size(rse_jxn_SRP009615)
## Obtain a list of RSE objects for all exon annotations
## Not run:
rses_exon <- lapply(annotation_options(), function(x) {
create_rse(proj_info, type = "exon", annotation = x, verbose = FALSE)
})
names(rses_exon) <- annotation_options()
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.