coerce-functions: Create SummarizedExperiment representations by transforming...

sparseSummarizedExperimentR Documentation

Create SummarizedExperiment representations by transforming ragged assays to rectangular form.

Description

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.

Usage

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
)

Arguments

x

RaggedExperiment

i

integer(1), character(1), or logical() selecting the assay to be transformed.

withDimnames

logical(1) default TRUE. propagate dimnames to SummarizedExperiment.

sparse

logical(1) whether to return a sparseMatrix representation

simplifyDisjoin

function of 1 argument, used to transform assays. See assay-functions.

query

GRanges provding regions over which reduction is to occur.

simplifyReduce

function of 3 arguments used to transform assays. See assay-functions.

Value

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().

sparseMatrix

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")`

Examples

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")


Bioconductor/RaggedExperiment documentation built on Feb. 8, 2024, 11:11 p.m.