| simpleRNASeq,BamFileList,RnaSeqParam-method {easyRNASeq} | R Documentation |
This function is a wrapper around the more low level functionalities of the
package. It is the simplest way to get a RangedSummarizedExperiment
object from a set of bam files. RangedSummarizedExperiment are
containers meant to hold any Next-Generation Sequencing experiment results and
metadata. The simpleRNASeq method replaces the
easyRNASeq function to
simplify the usability. It does the following:
use GenomicAlignments
for reading/pre-processing the BAM files.
get the annotations
depending on the selected parameters
calculate the coverage from the provided file(s)
summarizes the
read counts according to the selected summarization
returns a RangedSummarizedExperiment object.
## S4 method for signature 'BamFileList,RnaSeqParam' simpleRNASeq(bamFiles = BamFileList(), param = RnaSeqParam(), nnodes = 1, verbose = TRUE, override = FALSE)
bamFiles |
a |
param |
RnaSeqParam a |
nnodes |
The number of CPU cores to use in parallel |
verbose |
a logical to be report progress or not. |
override |
Should the provided parameters override the detected ones |
returns a RangedSummarizedExperiment object.
Nicolas Delhomme
For the input:
For the output:
RangedSummarizedExperiment
For related functions:
# the data
library(curl)
# get the example data files - we retrieve a set of example bam files
# from GitHub using curl, as well as their index.
invisible(sapply(c("ACACTG","ACTAGC"),function(bam){
curl_download(paste0("https://github.com/UPSCb/UPSCb/raw/",
"master/tutorial/easyRNASeq/",bam,".bam"),paste0(bam,".bam"))
curl_download(paste0("https://github.com/UPSCb/UPSCb/raw/",
"master/tutorial/easyRNASeq/",bam,".bam.bai"),paste0(bam,".bam.bai"))
}))
# and some annotation
invisible(curl_download(paste0("https://github.com/UPSCb/UPSCb/raw/",
"master/tutorial/easyRNASeq/Drosophila_melanogaster.BDGP5.77.with-chr.gtf.gz"),
"Drosophila_melanogaster.BDGP5.77.with-chr.gtf.gz"))
# get the BamFileList
bamFiles <- getBamFileList(dir(".",pattern="^[A,T].*\\.bam$",full.names=TRUE))
# create the AnnotParam
annotParam <- AnnotParam("Drosophila_melanogaster.BDGP5.77.with-chr.gtf.gz",
type="gtf")
# create the RnaSeqParam
rnaSeqParam <- RnaSeqParam(annotParam=annotParam)
# get a RangedSummarizedExperiment containing the counts table
sexp <- simpleRNASeq(
bamFiles=bamFiles,
param=rnaSeqParam,
verbose=TRUE
)
# get the counts
assays(sexp)$exons