Title: | Representation of Sparse Experiments and Assays Across Samples |
---|---|
Description: | This package provides a flexible representation of copy number, mutation, and other data that fit into the ragged array schema for genomic location data. The basic representation of such data provides a rectangular flat table interface to the user with range information in the rows and samples/specimen in the columns. The RaggedExperiment class derives from a GRangesList representation and provides a semblance of a rectangular dataset. |
Authors: | Martin Morgan [aut], Marcel Ramos [aut, cre] , Lydia King [ctb] |
Maintainer: | Marcel Ramos <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.31.0 |
Built: | 2024-10-31 03:41:59 UTC |
Source: | https://github.com/bioc/RaggedExperiment |
RaggedExperiment allows the user to represent, copy number, mutation, and other types of range-based data formats where optional information about samples can be provided. At the backbone of this package is the GRangesList class. The RaggedExperiment class uses this representation and presents the data in a couple of different ways:
rowRanges
colData
The rowRanges method will return the internal
GRangesList
representation of the dataset. A distinction
between the SummarizedExperiment and the
RaggedExperiment
classes is that the RaggedExperiment
class allows for ragged ranges, meaning that there may be a
different number of ranges or rows per sample.
Maintainer: Marcel Ramos [email protected] (ORCID)
Authors:
Martin Morgan [email protected]
Other contributors:
Lydia King [email protected] [contributor]
Useful links:
Report bugs at https://github.com/Bioconductor/RaggedExperiment/issues
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_, sparse = FALSE ) compactAssay( x, i = 1, withDimnames = TRUE, background = NA_integer_, sparse = FALSE ) disjoinAssay( x, simplifyDisjoin, i = 1, withDimnames = TRUE, background = NA_integer_ ) qreduceAssay( x, query, simplifyReduce, i = 1, withDimnames = TRUE, background = NA_integer_ )
sparseAssay( x, i = 1, withDimnames = TRUE, background = NA_integer_, sparse = FALSE ) compactAssay( x, i = 1, withDimnames = TRUE, background = NA_integer_, sparse = FALSE ) 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
|
sparse |
logical(1) whether to return a
|
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. Use 'sparse=TRUE' to obtain
a sparseMatrix
assay representation.
compactAssay()
: Samples with identical range are placed
in the same row. Non-disjoint ranges are NOT collapsed. Use
'sparse=TRUE' to obtain a sparseMatrix
assay
representation.
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
.
re4 <- RaggedExperiment(GRangesList( GRanges(c(A = "chr1:1-10:-", B = "chr1:8-14:-", C = "chr2:15-18:+"), score = 3:5), GRanges(c(D = "chr1:1-10:-", E = "chr2:11-18:+"), score = 1:2) ), colData = DataFrame(id = 1:2)) query <- GRanges(c("chr1:1-14:-", "chr2:11-18:+")) weightedmean <- function(scores, ranges, qranges) { ## weighted average score per query range ## the weight corresponds to the size of the overlap of each ## overlapping subject range with the corresponding query range isects <- pintersect(ranges, qranges) sum(scores * width(isects)) / sum(width(isects)) } qreduceAssay(re4, 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) library(curatedTCGAData) library(TCGAutils) }) ## TCGA MultiAssayExperiment with RaggedExperiment data mae <- curatedTCGAData("ACC", c("RNASeq2GeneNorm", "CNASNP", "Mutation"), version = "1.1.38", dry.run = FALSE) ## genomic coordinates gn <- genes(TxDb.Hsapiens.UCSC.hg19.knownGene) gn <- keepStandardChromosomes(granges(gn), pruning.mode="coarse") seqlevelsStyle(gn) <- "NCBI" genome(gn) gn <- unstrand(gn) ## reduce mutations, marking any genomic range with non-silent ## mutation as FALSE nonsilent <- function(scores, ranges, qranges) any(scores != "Silent") mre <- mae[["ACC_Mutation-20160128"]] seqlevelsStyle(rowRanges(mre)) <- "NCBI" ## hack to make genomes match genome(mre) <- paste0(correctBuild(unique(genome(mre)), "NCBI"), ".p13") mutations <- qreduceAssay(mre, gn, nonsilent, "Variant_Classification") genome(mre) <- correctBuild(unique(genome(mre)), "NCBI") ## reduce copy number re <- mae[["ACC_CNASNP-20160128"]] class(re) ## [1] "RaggedExperiment" seqlevelsStyle(re) <- "NCBI" genome(re) <- "GRCh37.p13" cn <- qreduceAssay(re, gn, weightedmean, "Segment_Mean") genome(re) <- "GRCh37" ## ALTERNATIVE ## ## TCGAutils helper function to convert RaggedExperiment objects to ## RangedSummarizedExperiment based on annotated gene ranges mae2 <- mae mae2[[1L]] <- re mae2[[2L]] <- mre qreduceTCGA(mae2) ## End(Not run)
re4 <- RaggedExperiment(GRangesList( GRanges(c(A = "chr1:1-10:-", B = "chr1:8-14:-", C = "chr2:15-18:+"), score = 3:5), GRanges(c(D = "chr1:1-10:-", E = "chr2:11-18:+"), score = 1:2) ), colData = DataFrame(id = 1:2)) query <- GRanges(c("chr1:1-14:-", "chr2:11-18:+")) weightedmean <- function(scores, ranges, qranges) { ## weighted average score per query range ## the weight corresponds to the size of the overlap of each ## overlapping subject range with the corresponding query range isects <- pintersect(ranges, qranges) sum(scores * width(isects)) / sum(width(isects)) } qreduceAssay(re4, 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) library(curatedTCGAData) library(TCGAutils) }) ## TCGA MultiAssayExperiment with RaggedExperiment data mae <- curatedTCGAData("ACC", c("RNASeq2GeneNorm", "CNASNP", "Mutation"), version = "1.1.38", dry.run = FALSE) ## genomic coordinates gn <- genes(TxDb.Hsapiens.UCSC.hg19.knownGene) gn <- keepStandardChromosomes(granges(gn), pruning.mode="coarse") seqlevelsStyle(gn) <- "NCBI" genome(gn) gn <- unstrand(gn) ## reduce mutations, marking any genomic range with non-silent ## mutation as FALSE nonsilent <- function(scores, ranges, qranges) any(scores != "Silent") mre <- mae[["ACC_Mutation-20160128"]] seqlevelsStyle(rowRanges(mre)) <- "NCBI" ## hack to make genomes match genome(mre) <- paste0(correctBuild(unique(genome(mre)), "NCBI"), ".p13") mutations <- qreduceAssay(mre, gn, nonsilent, "Variant_Classification") genome(mre) <- correctBuild(unique(genome(mre)), "NCBI") ## reduce copy number re <- mae[["ACC_CNASNP-20160128"]] class(re) ## [1] "RaggedExperiment" seqlevelsStyle(re) <- "NCBI" genome(re) <- "GRCh37.p13" cn <- qreduceAssay(re, gn, weightedmean, "Segment_Mean") genome(re) <- "GRCh37" ## ALTERNATIVE ## ## TCGAutils helper function to convert RaggedExperiment objects to ## RangedSummarizedExperiment based on annotated gene ranges mae2 <- mae mae2[[1L]] <- re mae2[[2L]] <- mre qreduceTCGA(mae2) ## End(Not run)
The RaggedExperiment
class is a container for
storing range-based data, including but not limited to copy
number data, and mutation data. It can store a collection of
GRanges
objects, as it is derived from the
GenomicRangesList
.
RaggedExperiment(..., colData = DataFrame(), metadata = list()) ## S4 method for signature 'RaggedExperiment' seqinfo(x) ## S4 replacement method for signature 'RaggedExperiment' seqinfo(x, new2old = NULL, pruning.mode = c("error", "coarse", "fine", "tidy")) <- value ## S4 method for signature 'RaggedExperiment' rowRanges(x, ...) ## S4 replacement method for signature 'RaggedExperiment,GRanges' rowRanges(x, ...) <- value ## S4 method for signature 'RaggedExperiment' mcols(x, use.names = FALSE, ...) ## S4 replacement method for signature 'RaggedExperiment' mcols(x, ...) <- value ## S4 method for signature 'RaggedExperiment' rowData(x, use.names = TRUE, ...) ## S4 replacement method for signature 'RaggedExperiment' rowData(x, ...) <- value ## S4 method for signature 'RaggedExperiment' dim(x) ## S4 method for signature 'RaggedExperiment' dimnames(x) ## S4 replacement method for signature 'RaggedExperiment,list' dimnames(x) <- value ## S4 replacement method for signature 'RaggedExperiment,ANY' dimnames(x) <- value ## S4 method for signature 'RaggedExperiment' length(x) ## S4 method for signature 'RaggedExperiment' colData(x, ...) ## S4 replacement method for signature 'RaggedExperiment,DataFrame' colData(x) <- value ## S4 method for signature 'RaggedExperiment,missing' assay(x, i, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment,ANY' assay(x, i, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment' assays(x, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment' assayNames(x, ...) ## S4 method for signature 'RaggedExperiment' show(object) ## S4 method for signature 'RaggedExperiment' as.list(x, ...) ## S4 method for signature 'RaggedExperiment' as.data.frame(x, row.names = NULL, optional = FALSE, ...) ## S4 method for signature 'RaggedExperiment' x$name ## S4 method for signature 'RaggedExperiment,ANY,ANY,ANY' x[i, j, ..., drop = TRUE] ## S4 method for signature 'RaggedExperiment,Vector' overlapsAny( query, subject, maxgap = 0L, minoverlap = 1L, type = c("any", "start", "end", "within", "equal"), ... ) ## S4 method for signature 'RaggedExperiment,Vector' subsetByOverlaps( x, ranges, maxgap = -1L, minoverlap = 0L, type = c("any", "start", "end", "within", "equal"), invert = FALSE, ... ) ## S4 method for signature 'RaggedExperiment' subset(x, subset, select, ...)
RaggedExperiment(..., colData = DataFrame(), metadata = list()) ## S4 method for signature 'RaggedExperiment' seqinfo(x) ## S4 replacement method for signature 'RaggedExperiment' seqinfo(x, new2old = NULL, pruning.mode = c("error", "coarse", "fine", "tidy")) <- value ## S4 method for signature 'RaggedExperiment' rowRanges(x, ...) ## S4 replacement method for signature 'RaggedExperiment,GRanges' rowRanges(x, ...) <- value ## S4 method for signature 'RaggedExperiment' mcols(x, use.names = FALSE, ...) ## S4 replacement method for signature 'RaggedExperiment' mcols(x, ...) <- value ## S4 method for signature 'RaggedExperiment' rowData(x, use.names = TRUE, ...) ## S4 replacement method for signature 'RaggedExperiment' rowData(x, ...) <- value ## S4 method for signature 'RaggedExperiment' dim(x) ## S4 method for signature 'RaggedExperiment' dimnames(x) ## S4 replacement method for signature 'RaggedExperiment,list' dimnames(x) <- value ## S4 replacement method for signature 'RaggedExperiment,ANY' dimnames(x) <- value ## S4 method for signature 'RaggedExperiment' length(x) ## S4 method for signature 'RaggedExperiment' colData(x, ...) ## S4 replacement method for signature 'RaggedExperiment,DataFrame' colData(x) <- value ## S4 method for signature 'RaggedExperiment,missing' assay(x, i, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment,ANY' assay(x, i, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment' assays(x, withDimnames = TRUE, ...) ## S4 method for signature 'RaggedExperiment' assayNames(x, ...) ## S4 method for signature 'RaggedExperiment' show(object) ## S4 method for signature 'RaggedExperiment' as.list(x, ...) ## S4 method for signature 'RaggedExperiment' as.data.frame(x, row.names = NULL, optional = FALSE, ...) ## S4 method for signature 'RaggedExperiment' x$name ## S4 method for signature 'RaggedExperiment,ANY,ANY,ANY' x[i, j, ..., drop = TRUE] ## S4 method for signature 'RaggedExperiment,Vector' overlapsAny( query, subject, maxgap = 0L, minoverlap = 1L, type = c("any", "start", "end", "within", "equal"), ... ) ## S4 method for signature 'RaggedExperiment,Vector' subsetByOverlaps( x, ranges, maxgap = -1L, minoverlap = 0L, type = c("any", "start", "end", "within", "equal"), invert = FALSE, ... ) ## S4 method for signature 'RaggedExperiment' subset(x, subset, select, ...)
... |
Constructor: GRanges, list of GRanges, or GRangesList OR assay: Additional arguments for assay. See details for more information. |
colData |
A |
metadata |
A |
x |
A RaggedExperiment object. |
new2old |
The
If Note that most of the times it's easier to proceed in 2 steps:
This 2-step approach will typically look like this: seqlevels(x) <- seqlevels(value) # align seqlevels seqinfo(x) <- seqinfo(value) # guaranteed to work Or, if seqlevels(x, pruning.mode="coarse") <- seqlevels(value) seqinfo(x) <- seqinfo(value) # guaranteed to work The |
pruning.mode |
When some of the seqlevels to drop from
See the "B. DROP SEQLEVELS FROM A LIST-LIKE OBJECT" section in the examples below for an extensive illustration of these pruning modes. |
value |
|
use.names |
(logical default FALSE) whether to propagate rownames from
the object to rownames of metadata |
i |
logical(1), integer(1), or character(1) indicating the
assay to be reported. For |
withDimnames |
logical (default TRUE) whether to use dimension names in the resulting object |
object |
A RaggedExperiment object. |
row.names |
|
optional |
logical. If |
name |
a literal character string or a name (possibly backtick
quoted). For extraction, this is normally (see under
‘Environments’) partially matched to the |
j |
integer(), character(), or logical() index selecting columns from RaggedExperiment |
drop |
logical (default TRUE) whether to drop empty samples |
query |
A RaggedExperiment instance. |
subject , ranges
|
Each of them can be an IntegerRanges (e.g.
IRanges, Views) or
IntegerRangesList (e.g. IRangesList,
ViewsList) derivative.
In addition, if If If both arguments are list-like objects with names, each list element from the 2nd argument is paired with the list element from the 1st argument with the matching name, if any. Otherwise, list elements are paired by position. The overlap is then computed between the pairs as described below. If |
maxgap |
A single integer >= -1. If If |
minoverlap |
A single non-negative integer. Only ranges with a minimum of When |
type |
By default, any overlap is accepted. By specifying the The |
invert |
If |
subset |
logical expression indicating elements or rows to keep: missing values are taken as false. |
select |
If If |
constructor returns a RaggedExperiment
object
'rowRanges' returns a GRanges
object
summarizing ranges corresponding to assay()
rows.
'rowRanges<-' returns a RaggedExperiment
object
with replaced ranges
'mcols' returns a DataFrame
object
of the metadata columns
'assays' returns a SimpleList
'overlapsAny' returns a logical vector of length equal
to the number of rows in the query
; TRUE
when the
copy number region overlaps the subject
.
'subsetByOverlaps' returns a RaggedExperiment containing
only copy number regions overlapping subject
.
seqinfo(RaggedExperiment)
: seqinfo accessor
seqinfo(RaggedExperiment) <- value
: Replace seqinfo metadata of the ranges
rowRanges(RaggedExperiment)
: rowRanges accessor
rowRanges(x = RaggedExperiment) <- value
: rowRanges replacement
mcols(RaggedExperiment)
: get the metadata columns of the ranges,
rectangular representation of the 'assays'
mcols(RaggedExperiment) <- value
: set the metadata columns of the ranges
corresponding to the assays
rowData(RaggedExperiment)
: get the rowData or metadata for the ranges
rowData(RaggedExperiment) <- value
: set the rowData or metadata for the ranges
dim(RaggedExperiment)
: get dimensions (number of sample-specific row
ranges by number of samples)
dimnames(RaggedExperiment)
: get row (sample-specific) range names
and sample names
dimnames(x = RaggedExperiment) <- value
: set row (sample-specific) range names
and sample names
dimnames(x = RaggedExperiment) <- value
: set row range names and sample names to NULL
length(RaggedExperiment)
: get the length of row vectors in the object,
similar to SummarizedExperiment
colData(RaggedExperiment)
: get column data
colData(x = RaggedExperiment) <- value
: change the colData
assay(x = RaggedExperiment, i = missing)
: assay missing method uses first metadata column
assay(x = RaggedExperiment, i = ANY)
: assay numeric method.
assays(RaggedExperiment)
: assays
assayNames(RaggedExperiment)
: names in each assay
show(RaggedExperiment)
: show method
as.list(RaggedExperiment)
: Allow extraction of metadata columns as a plain
list
as.data.frame(RaggedExperiment)
: Allow conversion to plain data.frame
$
: Easily access the colData
columns with
the dollar sign operator
x[i
: Subset a RaggedExperiment object
overlapsAny(query = RaggedExperiment, subject = Vector)
: Determine whether copy number ranges
defined by query
overlap ranges of subject
.
subsetByOverlaps(x = RaggedExperiment, ranges = Vector)
: Subset the RaggedExperiment to contain only
copy number ranges overlapping ranges of subject
.
subset(RaggedExperiment)
: subset helper function for dividing by rowData
and / or colData values
RaggedExperiment(..., colData=DataFrame())
: Creates a
RaggedExperiment object using multiple GRanges
objects or a list
of GRanges
objects. Additional column data may be provided
as a DataFrame
object.
In the following, 'x' represents a RaggedExperiment
object:
rowRanges(x)
:
Get the ranged data. Value is a GenomicRanges
object.
assays(x)
:
Get the assays. Value is a SimpleList
.
assay(x, i)
:
An alternative to assays(x)[[i]]
to get the ith
(default first) assay element.
mcols(x), mcols(x) <- value
:
Get or set the metadata columns. For RaggedExperiment
, the
columns correspond to the assay ith elements.
rowData(x), rowData(x) <- value
:
Get or set the row data. Value is a DataFrame
object. Also corresponds to the mcols
data.
Note for advanced users and developers. Both
mcols
and rowData
setters may reduce the size of the
internal RaggedExperiment
data representation. Particularly after
subsetting, the internal row index is modified and such setter
operations will use the index to subset the data and reduce the
"rows" of the internal data representation.
x[i, j]
:
Get ranges or elements (i
and j
, respectively) with
optional metadata columns where i
or j
can be missing,
an NA-free logical, numeric, or character vector.
In the following, 'object' represents a RaggedExperiment
object:
as(object, "GRangesList")
:
Creates a GRangesList object from a RaggedExperiment
.
as(from, "RaggedExperiment")
:
Creates a RaggedExperiment
object from a GRangesList,
or GRanges object.
## Create an empty RaggedExperiment instance re0 <- RaggedExperiment() re0 ## Create a couple of GRanges objects with row ranges names sample1 <- GRanges( c(a = "chr1:1-10:-", b = "chr1:11-18:+"), score = 1:2) sample2 <- GRanges( c(c = "chr2:1-10:-", d = "chr2:11-18:+"), score = 3:4) ## Include column data colDat <- DataFrame(id = 1:2) ## Create a RaggedExperiment object from a couple of GRanges re1 <- RaggedExperiment(sample1=sample1, sample2=sample2, colData = colDat) re1 ## With list of GRanges lgr <- list(sample1 = sample1, sample2 = sample2) ## Create a RaggedExperiment from a list of GRanges re2 <- RaggedExperiment(lgr, colData = colDat) grl <- GRangesList(sample1 = sample1, sample2 = sample2) ## Create a RaggedExperiment from a GRangesList re3 <- RaggedExperiment(grl, colData = colDat) ## Subset a RaggedExperiment assay(re3[c(1, 3),]) subsetByOverlaps(re3, GRanges("chr1:1-5")) # by ranges
## Create an empty RaggedExperiment instance re0 <- RaggedExperiment() re0 ## Create a couple of GRanges objects with row ranges names sample1 <- GRanges( c(a = "chr1:1-10:-", b = "chr1:11-18:+"), score = 1:2) sample2 <- GRanges( c(c = "chr2:1-10:-", d = "chr2:11-18:+"), score = 3:4) ## Include column data colDat <- DataFrame(id = 1:2) ## Create a RaggedExperiment object from a couple of GRanges re1 <- RaggedExperiment(sample1=sample1, sample2=sample2, colData = colDat) re1 ## With list of GRanges lgr <- list(sample1 = sample1, sample2 = sample2) ## Create a RaggedExperiment from a list of GRanges re2 <- RaggedExperiment(lgr, colData = colDat) grl <- GRangesList(sample1 = sample1, sample2 = sample2) ## Create a RaggedExperiment from a GRangesList re3 <- RaggedExperiment(grl, colData = colDat) ## Subset a RaggedExperiment assay(re3[c(1, 3),]) subsetByOverlaps(re3, GRanges("chr1:1-5")) # by ranges
These methods transform RaggedExperiment
objects to similar SummarizedExperiment
objects. They do
so by transforming assay data to more rectangular
representations, following the rules outlined for similarly
names transformations sparseAssay()
,
compactAssay()
, disjoinAssay()
, and
qreduceAssay()
. Because of the complexity of the
transformation, ti usually only makes sense transform
RaggedExperiment
objects with a single assay; this is
currently enforced at time of coercion.
sparseSummarizedExperiment(x, i = 1, withDimnames = TRUE, sparse = FALSE) compactSummarizedExperiment(x, i = 1L, withDimnames = TRUE, sparse = FALSE) disjoinSummarizedExperiment(x, simplifyDisjoin, i = 1L, withDimnames = TRUE) qreduceSummarizedExperiment( x, query, simplifyReduce, i = 1L, withDimnames = TRUE )
sparseSummarizedExperiment(x, i = 1, withDimnames = TRUE, sparse = FALSE) compactSummarizedExperiment(x, i = 1L, withDimnames = TRUE, sparse = FALSE) disjoinSummarizedExperiment(x, simplifyDisjoin, i = 1L, withDimnames = TRUE) qreduceSummarizedExperiment( x, query, simplifyReduce, i = 1L, withDimnames = TRUE )
x |
|
i |
|
withDimnames |
|
sparse |
logical(1) whether to return a
|
simplifyDisjoin |
|
query |
|
simplifyReduce |
|
All functions return RangedSummarizedExperiment
.
sparseSummarizedExperiment
has rowRanges()
identical to the row ranges of x
, and assay()
data as sparseAssay()
. This is very space-inefficient
representation of ragged data. Use 'sparse=TRUE' to obtain
a sparseMatrix
assay representation.
compactSummarizedExperiment
has rowRanges()
identical to the row ranges of x
, and assay()
data as compactAssay()
. This is space-inefficient
representation of ragged data when samples are primarily
composed of different ranges. Use 'sparse=TRUE' to obtain
a sparseMatrix
assay representation.
disjoinSummarizedExperiment
has rowRanges()
identical to the disjoint row ranges of x
,
disjoint(rowRanges(x))
, and assay()
data as
disjoinAssay()
.
qreduceSummarizedExperiment
has rowRanges()
identical to query
, and assay()
data as
qreduceAssay()
.
Convert a dgCMatrix
to a RaggedExperiment
given that the rownames
are coercible to GRanges
.
In the following example, x
is a dgCMatrix
from the Matrix
package.
`as(x, "RaggedExperiment")`
x <- RaggedExperiment(GRangesList( GRanges(c("A:1-5", "A:4-6", "A:10-15"), score=1:3), GRanges(c("A:1-5", "B:1-3"), score=4:5) )) ## sparseSummarizedExperiment sse <- sparseSummarizedExperiment(x) assay(sse) rowRanges(sse) ## compactSummarizedExperiment cse <- compactSummarizedExperiment(x) assay(cse) rowRanges(cse) ## disjoinSummarizedExperiment disjoinAssay(x, lengths) dse <- disjoinSummarizedExperiment(x, lengths) assay(dse) rowRanges(dse) ## qreduceSummarizedExperiment 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 ## the weight corresponds to the size of the overlap of each ## overlapping subject range with the corresponding query range isects <- pintersect(ranges, qranges) sum(scores * width(isects)) / sum(width(isects)) } qreduceAssay(x, query, weightedmean) qse <- qreduceSummarizedExperiment(x, query, weightedmean) assay(qse) rowRanges(qse) sm <- Matrix::sparseMatrix( i = c(2, 3, 4, 3, 4, 3, 4), j = c(1, 1, 1, 3, 3, 4, 4), x = c(2L, 4L, 2L, 2L, 2L, 4L, 2L), dims = c(4, 4), dimnames = list( c("chr2:1-10", "chr2:2-10", "chr2:3-10", "chr2:4-10"), LETTERS[1:4] ) ) as(sm, "RaggedExperiment")
x <- RaggedExperiment(GRangesList( GRanges(c("A:1-5", "A:4-6", "A:10-15"), score=1:3), GRanges(c("A:1-5", "B:1-3"), score=4:5) )) ## sparseSummarizedExperiment sse <- sparseSummarizedExperiment(x) assay(sse) rowRanges(sse) ## compactSummarizedExperiment cse <- compactSummarizedExperiment(x) assay(cse) rowRanges(cse) ## disjoinSummarizedExperiment disjoinAssay(x, lengths) dse <- disjoinSummarizedExperiment(x, lengths) assay(dse) rowRanges(dse) ## qreduceSummarizedExperiment 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 ## the weight corresponds to the size of the overlap of each ## overlapping subject range with the corresponding query range isects <- pintersect(ranges, qranges) sum(scores * width(isects)) / sum(width(isects)) } qreduceAssay(x, query, weightedmean) qse <- qreduceSummarizedExperiment(x, query, weightedmean) assay(qse) rowRanges(qse) sm <- Matrix::sparseMatrix( i = c(2, 3, 4, 3, 4, 3, 4), j = c(1, 1, 1, 3, 3, 4, 4), x = c(2L, 4L, 2L, 2L, 2L, 4L, 2L), dims = c(4, 4), dimnames = list( c("chr2:1-10", "chr2:2-10", "chr2:3-10", "chr2:4-10"), LETTERS[1:4] ) ) as(sm, "RaggedExperiment")