srnadiffExp | R Documentation |
srnadiffExp
is an S4 class providing the infrastructure (slots)
to store the input data, methods parameters, intermediate calculations
and results of a sRNA-diff approach.
srnadiffExp( bamFiles = NULL, sampleInfo = NULL, annotReg = NULL, diffMethod = NULL, normFactors = NULL ) ## S4 method for signature 'srnadiffExp' show(object)
bamFiles |
A vector with the full paths to the BAM files. |
sampleInfo |
A |
annotReg |
Optional annotation information. Annotated regions as a
|
diffMethod |
A character storing the name of the method used
to compute the p-values. Can be |
normFactors |
A numeric vector, one size factor for each sample in the data. |
object |
An |
srnadiffExp
load and summarize sample BAM
files into base-resolution coverage and estimate the size factors (the
effective library size) from the coverage data.
To facilitate programming pipelines, NULL
values are input
for regions
, parameters
and countMatrix
slots,
in which case the default value is used as if the argument had been
missing. These slots will be updated after differential expression
(srnadiff
) approach.
srnadiffExp
constructor returns an srnadiffExp
object of class S4.
The show
method informatively display object contents.
bamFiles
A BamFileList
object
with the full paths to the BAM files.
sampleInfo
A data.frame
with sample and experimental
design information. Each row describes one sample.
annotReg
A GRanges
with annotation information.
diffMethod
A character storing the name of the method used
to compute the p-values. Can be DESeq2
,
edgeR
, or baySeq
.
chromosomeSizes
A named vector with the sizes of the chromosomes.
coverages
The sample coverages, a named
RleList
object.
normFactors
A vector of normalization factors.
regions
A GenomicRanges
of the candidate
differentially expressed regions.
countMatrix
A matrix of non-negative integer count values, one row per region and one column per sample.
parameters
An named list
. The parameters for the
segmentation methods. See parameters
.
basedir <- system.file("extdata", package="srnadiff", mustWork = TRUE) sampleInfo <- read.csv(file.path(basedir, "dataInfo.csv")) gtfFile <- file.path(basedir, "Homo_sapiens.GRCh38.76.gtf.gz") annotReg <- readAnnotation(gtfFile, feature="gene", source="miRNA") bamFiles <- file.path(basedir, sampleInfo$FileName) srnaExp <- srnadiffExp(bamFiles, sampleInfo, annotReg) srnaExp
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.