1 Introduction

The RaggedExperiment package provides a flexible data representation for copy number, mutation and other ragged array schema for genomic location data. It aims to provide a framework for a set of samples that have differing numbers of genomic ranges.

The RaggedExperiment class derives from a GRangesList representation and provides a semblance of a rectangular dataset. The row and column dimensions of the RaggedExperiment correspond to the number of ranges in the entire dataset and the number of samples represented in the data, respectively.

2 Installation

source("https://bioconductor.org/biocLite.R")
BiocInstaller::biocLite("RaggedExperiment")

Loading the package:

library(RaggedExperiment)

3 Constructing a RaggedExperiment object

We start with a couple of GRanges objects, each representing an individual sample:

sample1 <- GRanges(
    c(A = "chr1:1-10:-", B = "chr1:8-14:+", C = "chr2:15-18:+"),
    score = 3:5)
sample2 <- GRanges(
    c(D = "chr1:1-10:-", E = "chr2:11-18:+"),
    score = 1:2)

Include column data colData to describe the samples:

colDat <- DataFrame(id = 1:2)

3.1 Using GRanges objects

ragexp <- RaggedExperiment(sample1 = sample1,
                           sample2 = sample2,
                           colData = colDat)
ragexp
## class: RaggedExperiment 
## dim: 5 2 
## assays(1): score
## rownames(5): A B C D E
## colnames(2): sample1 sample2
## colData names(1): id

3.2 Using a GRangesList instance

grl <- GRangesList(sample1 = sample1, sample2 = sample2)
RaggedExperiment(grl, colData = colDat)
## class: RaggedExperiment 
## dim: 5 2 
## assays(1): score
## rownames(5): A B C D E
## colnames(2): sample1 sample2
## colData names(1): id

3.3 Using a list of GRanges

rangeList <- list(sample1 = sample1, sample2 = sample2)
RaggedExperiment(rangeList, colData = colDat)
## class: RaggedExperiment 
## dim: 5 2 
## assays(1): score
## rownames(5): A B C D E
## colnames(2): sample1 sample2
## colData names(1): id

3.4 Using a List of GRanges with metadata

Note: In cases where a SimpleGenomicRangesList is provided along with accompanying metadata (accessed by mcols), the metadata is used as the colData for the RaggedExperiment.

grList <- List(sample1 = sample1, sample2 = sample2)
mcols(grList) <- colDat
RaggedExperiment(grList)
## class: RaggedExperiment 
## dim: 5 2 
## assays(1): score
## rownames(5): A B C D E
## colnames(2): sample1 sample2
## colData names(1): id

4 Accessors

4.1 Range data

rowRanges(ragexp)
## GRanges object with 5 ranges and 0 metadata columns:
##     seqnames    ranges strand
##        <Rle> <IRanges>  <Rle>
##   A     chr1  [ 1, 10]      -
##   B     chr1  [ 8, 14]      +
##   C     chr2  [15, 18]      +
##   D     chr1  [ 1, 10]      -
##   E     chr2  [11, 18]      +
##   -------
##   seqinfo: 2 sequences from an unspecified genome; no seqlengths

4.2 Dimension names

dimnames(ragexp)
## [[1]]
## [1] "A" "B" "C" "D" "E"
## 
## [[2]]
## [1] "sample1" "sample2"

4.3 colData

colData(ragexp)
## DataFrame with 2 rows and 1 column
##                id
##         <integer>
## sample1         1
## sample2         2

5 Subsetting

5.1 by dimension

Subsetting a RaggedExperiment is akin to subsetting a matrix object. Rows correspond to genomic ranges and columns to samples or specimen. It is possible to subset using integer, character, and logical indices.

5.2 by genomic ranges

The overlapsAny and subsetByOverlaps functionalities are available for use for RaggedExperiment. Please see the corresponding documentation in RaggedExperiment and GenomicRanges.

6 *Assay functions

RaggedExperiment package provides several different functions for representing ranged data in a rectangular matrix via the *Assay methods.

6.1 sparseAssay

The most straightforward matrix representation of a RaggedExperiment will return a matrix of dimensions equal to the product of the number of ranges and samples.

dim(ragexp)
## [1] 5 2
Reduce(`*`, dim(ragexp))
## [1] 10
sparseAssay(ragexp)
##   sample1 sample2
## A       3      NA
## B       4      NA
## C       5      NA
## D      NA       1
## E      NA       2
length(sparseAssay(ragexp))
## [1] 10

6.2 compactAssay

Samples with identical ranges are placed in the same row. Non-disjoint ranges are not collapsed.

compactAssay(ragexp)
##              sample1 sample2
## chr1:8-14:+        4      NA
## chr1:1-10:-        3       1
## chr2:11-18:+      NA       2
## chr2:15-18:+       5      NA

6.3 disjoinAssay

This function returns a matrix of disjoint ranges across all samples. Elements of the matrix are summarized by applying the simplifyDisjoin functional argument to assay values of overlapping ranges.

disjoinAssay(ragexp, simplifyDisjoin = mean)
##              sample1 sample2
## chr1:8-14:+        4      NA
## chr1:1-10:-        3       1
## chr2:11-14:+      NA       2
## chr2:15-18:+       5       2

6.4 qreduceAssay

The qreduceAssay function works with a query parameter that highlights a window of ranges for the resulting matrix. The returned matrix will have dimensions length(query) by ncol(x). Elements contain assay values for the i th query range and the j th sample, summarized according to the simplifyReduce functional argument.

For demonstration purposes, we can have a look at the original GRangesList and the associated scores from which the current ragexp object is derived:

unlist(grl, use.names = FALSE)
## GRanges object with 5 ranges and 1 metadata column:
##     seqnames    ranges strand |     score
##        <Rle> <IRanges>  <Rle> | <integer>
##   A     chr1  [ 1, 10]      - |         3
##   B     chr1  [ 8, 14]      + |         4
##   C     chr2  [15, 18]      + |         5
##   D     chr1  [ 1, 10]      - |         1
##   E     chr2  [11, 18]      + |         2
##   -------
##   seqinfo: 2 sequences from an unspecified genome; no seqlengths

This data is represented as rowRanges and assays in RaggedExperiment:

rowRanges(ragexp)
## GRanges object with 5 ranges and 0 metadata columns:
##     seqnames    ranges strand
##        <Rle> <IRanges>  <Rle>
##   A     chr1  [ 1, 10]      -
##   B     chr1  [ 8, 14]      +
##   C     chr2  [15, 18]      +
##   D     chr1  [ 1, 10]      -
##   E     chr2  [11, 18]      +
##   -------
##   seqinfo: 2 sequences from an unspecified genome; no seqlengths
assay(ragexp, "score")
##   sample1 sample2
## A       3      NA
## B       4      NA
## C       5      NA
## D      NA       1
## E      NA       2

Here we provide the “query” region of interest:

(query <- GRanges(c("chr1:1-14:-", "chr2:11-18:+")))
## GRanges object with 2 ranges and 0 metadata columns:
##       seqnames    ranges strand
##          <Rle> <IRanges>  <Rle>
##   [1]     chr1  [ 1, 14]      -
##   [2]     chr2  [11, 18]      +
##   -------
##   seqinfo: 2 sequences from an unspecified genome; no seqlengths

The simplifyReduce argument in qreduceAssay allows the user to summarize overlapping regions with a custom method for the given “query” region of interest. We provide one for calculating a weighted average score per query range, where the weight is proportional to the overlap widths between overlapping ranges and a query range.

Note that there are three arguments to this function. See the documentation for additional details.

weightedmean <- function(scores, ranges, qranges)
{
    isects <- pintersect(ranges, qranges)
    sum(scores * width(isects)) / sum(width(isects))
}

A call to qreduceAssay involves the RaggedExperiment, the GRanges query and the simplifyReduce functional argument.

qreduceAssay(ragexp, query, simplifyReduce = weightedmean)
##              sample1 sample2
## chr1:1-14:-        3       1
## chr2:11-18:+       5       2

7 Coercion

The RaggedExperiment provides a family of parallel functions for coercing to the SummarizedExperiment class. By selecting a particular assay index (i), a parallel assay coercion method can be achieved.

Here is the list of function names:

  • sparseSummarizedExperiment
  • compactSummarizedExperiment
  • disjoinSummarizedExperiment
  • qreduceSummarizedExperiment

See the documentation for details.

8 Session Information

sessionInfo()
## R version 3.4.3 (2017-11-30)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.3 LTS
## 
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.6-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.6-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
## [1] RaggedExperiment_1.2.5 GenomicRanges_1.30.0   GenomeInfoDb_1.14.0   
## [4] IRanges_2.12.0         S4Vectors_0.16.0       BiocGenerics_0.24.0   
## [7] BiocStyle_2.6.1       
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_0.12.14               knitr_1.17                
##  [3] XVector_0.18.0             magrittr_1.5              
##  [5] zlibbioc_1.24.0            lattice_0.20-35           
##  [7] stringr_1.2.0              tools_3.4.3               
##  [9] grid_3.4.3                 SummarizedExperiment_1.8.0
## [11] Biobase_2.38.0             matrixStats_0.52.2        
## [13] htmltools_0.3.6            yaml_2.1.16               
## [15] rprojroot_1.2              digest_0.6.12             
## [17] bookdown_0.5               Matrix_1.2-12             
## [19] GenomeInfoDbData_0.99.1    bitops_1.0-6              
## [21] RCurl_1.95-4.8             evaluate_0.10.1           
## [23] rmarkdown_1.8              DelayedArray_0.4.1        
## [25] stringi_1.1.6              compiler_3.4.3            
## [27] backports_1.1.2