Description Usage Arguments Value Accessors Examples
The main object MPRAnalyze works with, contains the input data, associated annotations, model parameters and analysis results.
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 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 | MpraObject(
dnaCounts,
rnaCounts,
dnaAnnot = NULL,
rnaAnnot = NULL,
colAnnot = NULL,
controls = NA_integer_,
rowAnnot = NULL,
BPPARAM = NULL
)
## S4 method for signature 'matrix'
MpraObject(
dnaCounts,
rnaCounts,
dnaAnnot = NULL,
rnaAnnot = NULL,
colAnnot = NULL,
controls = NA_integer_,
rowAnnot = NULL,
BPPARAM = NULL
)
## S4 method for signature 'SummarizedExperiment'
MpraObject(
dnaCounts,
rnaCounts,
dnaAnnot = NULL,
rnaAnnot = NULL,
colAnnot = NULL,
controls = NA,
rowAnnot = NULL,
BPPARAM = NULL
)
dnaCounts(obj)
## S4 method for signature 'MpraObject'
dnaCounts(obj)
rnaCounts(obj)
## S4 method for signature 'MpraObject'
rnaCounts(obj)
dnaAnnot(obj)
## S4 method for signature 'MpraObject'
dnaAnnot(obj)
rnaAnnot(obj)
## S4 method for signature 'MpraObject'
rnaAnnot(obj)
rowAnnot(obj)
## S4 method for signature 'MpraObject'
rowAnnot(obj)
controls(obj)
## S4 method for signature 'MpraObject'
controls(obj)
dnaDepth(obj)
## S4 method for signature 'MpraObject'
dnaDepth(obj)
rnaDepth(obj)
## S4 method for signature 'MpraObject'
rnaDepth(obj)
model(obj)
## S4 method for signature 'MpraObject'
model(obj)
|
dnaCounts |
the DNA count matrix, or a SummarizedExperiment object containing the DNA Counts and column annotations for the DNA data. If the input is a SummarizedExperiment object, the dnaAnnot (or colAnnot) arguments will be ignored |
rnaCounts |
the RNA count matrix, or a SummarizedExperiment object containing the RNA Counts and column annotations for the RNA data. If the input is a SummarizedExperiment object, the rnaAnnot (or colAnnot) arguments will be ignored |
dnaAnnot |
data.frame with the DNA column (sample) annotations |
rnaAnnot |
data.frame with the RNA column (sample) annotations |
colAnnot |
if annotations for DNA and RNA are identical, they can be set at the same time using colAnnot instead of using both rnaAnnot and dnaAnnot |
controls |
IDs of the rows in the matrices that correspond to negative control enhancers. These are used to establish the null for quantification purposes, and to correct systemic bias in comparative analyses. Can be a character vectors (corresponding to rownames in the data matrices), logical or numeric indices. |
rowAnnot |
a data.frame with the row (candidate enhancer) annotations. The names must match the row names in the DNA and RNA count matrices. |
BPPARAM |
a parallelization backend using the BiocParallel package, see more details [here](http://bioconductor.org/packages/release/bioc/html/BiocParallel.html) |
obj |
The MpraObject to extract properties from |
an initialized MpraObject
MpraObject properties can be accessed using accessor functions
the DNA count matrix
the RNA count matrix
data.frame with the DNA column (sample) annotations
data.frame with the RNA column (sample) annotations
data.frame with the row (candidate enhancers) annotations
the distributional model used. the Gamma-Poisson convolutional
model is used by default. see setModel
The library size correction factors computed for the DNA
libraries. These are computed by the estimateDepthFactors
function and can be set manually using the setDepthFactors
function
The library size correction factors computed for the RNA
libraries These are computed by the estimateDepthFactors
function and can be set manually using the setDepthFactors
function
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 | data <- simulateMPRA(tr = rep(2,10), da=c(rep(2,5), rep(2.5,5)),
nbatch=2, nbc=20)
## use 3 of the non-active enhancers as controls
obj <- MpraObject(dnaCounts = data$obs.dna,
rnaCounts = data$obs.rna,
colAnnot = data$annot,
controls = as.integer(c(1,2,4)))
## alternatively, initialize the object with SummarizedExperiment objects:
## Not run:
se.DNA <- SummarizedExperiment(list(data$obs.dna), colData=data$annot)
se.RNA <- SummarizedExperiment(list(data$obs.rna), colData=data$annot)
obj <- MpraObject(dnaCounts = se.DNA, rnaCounts = rna.se,
controls = as.integer(c(1,2,4)))
## End(Not run)
dnaCounts <- dnaCounts(obj)
rnaCounts <- rnaCounts(obj)
dnaAnnot <- dnaAnnot(obj)
rnaAnnot <- rnaAnnot(obj)
controls <- controls(obj)
rowAnnot <- rowAnnot(obj)
model <- model(obj)
obj <- estimateDepthFactors(obj, lib.factor=c("batch", "condition"))
dnaDepth <- dnaDepth(obj)
rnaDepth <- rnaDepth(obj)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.