| assay-functions {RaggedExperiment} | R Documentation |
These methods transform assay() from the
default (i.e., sparseAssay()) representation to various
forms of more dense representation. compactAssay()
collapses identical ranges across samples into a single
row. disjoinAssay() creates disjoint (non-overlapping)
regions, simplifies values within each sample in a
user-specified manner, and returns a matrix of disjoint regions
x samples.
This method transforms assay() from the default
(i.e., sparseAssay()) representation to a reduced
representation summarizing each original range overlapping
ranges in query. Reduction in each cell can be tailored
to indivdual needs using the simplifyReduce functional argument.
sparseAssay(x, i = 1, withDimnames = TRUE, background = NA_integer_) compactAssay(x, i = 1, withDimnames = TRUE, background = NA_integer_) disjoinAssay(x, simplifyDisjoin, i = 1, withDimnames = TRUE, background = NA_integer_) qreduceAssay(x, query, simplifyReduce, i = 1, withDimnames = TRUE, background = NA_integer_)
x |
A |
i |
integer(1) or character(1) name of assay to be transformed. |
withDimnames |
logical(1) include dimnames on the returned
matrix. When there are no explict rownames, these are
manufactured with |
background |
A value (default NA) for the returned matrix after
|
simplifyDisjoin |
A
a
original: |-----------|
b
|----------|
a a, b b
disjoint: |----|------|---|
values <- IntegerList(a, c(a, b), b)
simplifyDisjoin(values)
|
query |
|
simplifyReduce |
A
|
sparseAssay(): A matrix() with dimensions
dim(x). Elements contain the assay value for the ith
range and jth sample.
compactAssay(): Samples with identical range are placed
in the same row. Non-disjoint ranges are NOT collapsed.
disjoinAssay(): A matrix with number of rows equal
to number of disjoint ranges across all samples. Elements of
the matrix are summarized by applying simplifyDisjoin() to
assay values of overlapping ranges
qreduceAssay(): A matrix() with dimensions
length(query) x ncol(x). Elements contain assay
values for the ith query range and jth sample, summarized
according to the function simplifyReduce.
x <- RaggedExperiment(GRangesList(
GRanges(c("A:1-3", "A:4-5", "A:10-15"), score=1:3),
GRanges(c("A:4-5", "B:1-3"), score=4:5)
))
query <- GRanges(c("A:1-2", "A:4-5", "B:1-5"))
weightedmean <- function(scores, ranges, qranges)
## weighted average score per query range
sum(scores * width(ranges)) / sum(width(ranges))
qreduceAssay(x, query, weightedmean)
## Not run:
## Extended example: non-silent mutations, summarized by genic
## region
suppressPackageStartupMessages({
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
library(GenomeInfoDb)
library(MultiAssayExperiment)
})
## TCGA Multi-assay experiment to RaggedExperiment
url <- "http://s3.amazonaws.com/multiassayexperiments/accMAEO.rds"
## download.file(url, fl <- tempfile())
## fl <- "accMAEO.rds"
mae <- readRDS(fl)[, , c("RNASeq2GeneNorm", "CNASNP", "Mutations")]
## genomic coordinates
gn <- genes(TxDb.Hsapiens.UCSC.hg19.knownGene)
gn <- keepStandardChromosomes(granges(gn), pruning.mode="coarse")
seqlevelsStyle(gn) <- "NCBI"
## reduce mutations, marking any genomic range with non-silent
## mutation as FALSE
nonsilent <- function(scores, ranges, qranges)
any(scores != "Silent")
re <- as(mae[["Mutations"]], "RaggedExperiment")
mutations <- qreduceAssay(re, gn, nonsilent, "Variant_Classification")
## reduce copy number
re <- as(mae[["CNASNP"]], "RaggedExperiment")
cn <- qreduceAssay(re, gn, weightedmean, "Segment_Mean")
## End(Not run)