Contents

1 Introduction

The mpra package provides tools for the analysis of data from massively parallel reporter assays (MPRA). Specifically, it contains the functionality described in (Myint et al. 2019). The primary analysis purpose is to enable differential analysis of activity measures, but the package can also be used to generate precision weights useful in regression analyses of activity scores on sequence features. The main workhorse of the mpra package is the mpralm() function which draws on the previously proposed voom framework for RNA-seq analysis (Law et al. 2014). In this document, we will be looking at MPRA data from a study comparing episomal and lentiviral versions of MPRA (Inoue et al. 2017). We will also look at MPRA data from a study comparing the regulatory activity of different alleles of thousands of SNPs (Tewhey et al. 2016).

1.1 How to cite

If you are using this package, please cite (Myint et al. 2019). If you are using the mpralm() function, it would be appropriate to also cite the voom framework (Law et al. 2014).

2 Dependencies

This document has the following dependencies

library(mpra)

3 Creating an MPRASet object

In this package, MPRA data are contained in MPRASet objects. Because MPRA data do not have a common prescribed format, these objects must be created manually. In this section, we demonstrate how to do this.

MPRASet objects must contain DNA and RNA count information because this is the information used to quantify activity levels of the elements being assayed. DNA and RNA count information should be specified as \(K \times S\) integer count matrices where \(K\) is the total number of barcodes over all elements if barcode-level information is being supplied or the total number of putative regulatory elements (PREs) if element-level information is being supplied. \(S\) is the number of samples (typically, the number of independent transfections).

MPRASet objects must also contain element identification information. This should be supplied as a character vector of length \(K\), the number of rows in the DNA and RNA count matrices. These are any strings used to describe/identify the unique PREs being assayed.

Optionally, the barcode sequences and PRE sequences can be specified as length \(K\) character vectors.

In the next sections we provide specific examples for how to specify this information for two common differential analysis settings: tissue and allele comparisons. Although we show simulated data, this information would typically be read from text files.

3.1 Tissue comparison

In tissue comparison studies, the same set of PREs is assayed in two or more cell types. In the following example, the experiment looks at four PREs with three barcodes each. Two tissues (liver and kidney) are studied, and each tissue has four replicates (four independent transfections each).

RNA and DNA count matrices would look as below:

E <- 4 # Number of elements
B <- 3 # Number of barcodes
s <- 4 # Samples per tissue
nt <- 2 # Number of tissues

set.seed(434)
rna <- matrix(rpois(E*B*s*nt, lambda = 30), nrow = E*B, ncol = s*nt)
dna <- matrix(rpois(E*B*s*nt, lambda = 30), nrow = E*B, ncol = s*nt)

rn <- as.character(outer(paste0("barcode_", seq_len(B), "_"), paste0("elem_", seq_len(E)), FUN = "paste0"))
cn <- c(paste0("liver_", seq_len(s)), paste0("kidney_", seq_len(s)))

rownames(rna) <- rn
rownames(dna) <- rn
colnames(rna) <- cn
colnames(dna) <- cn

rna
##                  liver_1 liver_2 liver_3 liver_4 kidney_1 kidney_2 kidney_3
## barcode_1_elem_1      31      26      36      18       27       33       28
## barcode_2_elem_1      33      28      29      32       28       33       34
## barcode_3_elem_1      22      43      25      41       28       34       30
## barcode_1_elem_2      28      27      31      20       37       33       24
## barcode_2_elem_2      32      31      32      21       38       22       27
## barcode_3_elem_2      35      36      28      38       26       23       30
## barcode_1_elem_3      30      34      26      34       40       26       26
## barcode_2_elem_3      37      27      37      29       25       29       38
## barcode_3_elem_3      34      30      20      26       30       35       33
## barcode_1_elem_4      28      27      34      26       29       30       25
## barcode_2_elem_4      31      32      28      36       31       24       29
## barcode_3_elem_4      24      28      34      31       37       29       34
##                  kidney_4
## barcode_1_elem_1       30
## barcode_2_elem_1       30
## barcode_3_elem_1       37
## barcode_1_elem_2       29
## barcode_2_elem_2       29
## barcode_3_elem_2       29
## barcode_1_elem_3       22
## barcode_2_elem_3       28
## barcode_3_elem_3       38
## barcode_1_elem_4       30
## barcode_2_elem_4       32
## barcode_3_elem_4       36
dna
##                  liver_1 liver_2 liver_3 liver_4 kidney_1 kidney_2 kidney_3
## barcode_1_elem_1      29      34      28      27       27       25       30
## barcode_2_elem_1      43      37      34      33       43       26       26
## barcode_3_elem_1      35      23      32      27       19       29       29
## barcode_1_elem_2      23      30      34      27       32       18       24
## barcode_2_elem_2      32      29      42      21       37       34       32
## barcode_3_elem_2      26      23      22      24       32       41       24
## barcode_1_elem_3      29      38      34      27       31       27       31
## barcode_2_elem_3      30      36      32      21       34       29       26
## barcode_3_elem_3      34      28      42      35       26       32       32
## barcode_1_elem_4      34      25      29      30       44       26       23
## barcode_2_elem_4      29      26      37      37       28       37       34
## barcode_3_elem_4      24      28      33      40       31       30       19
##                  kidney_4
## barcode_1_elem_1       31
## barcode_2_elem_1       33
## barcode_3_elem_1       25
## barcode_1_elem_2       25
## barcode_2_elem_2       31
## barcode_3_elem_2       23
## barcode_1_elem_3       33
## barcode_2_elem_3       39
## barcode_3_elem_3       26
## barcode_1_elem_4       33
## barcode_2_elem_4       31
## barcode_3_elem_4       26

PRE identification strings would look as below. When counts are provided at the barcode level, the eid character vector will have repeated elements.

eid <- rep(paste0("elem_", seq_len(E)), each = B)
eid
##  [1] "elem_1" "elem_1" "elem_1" "elem_2" "elem_2" "elem_2" "elem_3" "elem_3"
##  [9] "elem_3" "elem_4" "elem_4" "elem_4"

We may also have PRE sequences as below. These sequences must be specified in a character vector of the same length as eid and the same number of rows as rna and dna.

eseq <- replicate(E, paste(sample(c("A", "T", "C", "G"), 10, replace = TRUE), collapse = ""))
eseq <- rep(eseq, each = B)
eseq
##  [1] "CAACGTTTGT" "CAACGTTTGT" "CAACGTTTGT" "TTCTCTTTGA" "TTCTCTTTGA"
##  [6] "TTCTCTTTGA" "ACCTACCTCA" "ACCTACCTCA" "ACCTACCTCA" "CCGTACTACT"
## [11] "CCGTACTACT" "CCGTACTACT"

The above pieces (rna, dna, eid, and eseq) can be supplied as arguments to the MPRASet constructor function as below. If barcode (barcode sequences) or eseq (PRE sequences) is not supplied, it must be specified as NULL.

mpraset_example <- MPRASet(DNA = dna, RNA = rna, eid = eid, eseq = eseq, barcode = NULL)
mpraset_example
## class: MPRASet 
## dim: 12 8 
## metadata(0):
## assays(2): DNA RNA
## rownames(12): barcode_1_elem_1 barcode_2_elem_1 ... barcode_2_elem_4
##   barcode_3_elem_4
## rowData names(2): eid eseq
## colnames(8): liver_1 liver_2 ... kidney_3 kidney_4
## colData names(0):
## No barcodes present

3.2 Allele comparison

In allele comparison studies, PREs that exist with two or more alleles are assayed. All allelic-versions of the PREs are assayed in the same sample. Because activity comparisons between alleles is desired, these counts must be separated into different columns. In the following example, the experiment looks at four PREs with three barcodes each. There are two alleles per PRE, and there are four replicates (four independent transfections total).

Note that because the different alleles of a single PRE are linked to different barcodes, there is not a natural way to construct the RNA and DNA count matrices as above, where a particular barcode is in a row. Further, sometimes each PRE-allele combination is paired with varying numbers of barcodes. This is yet another reason that DNA and RNA count matrices should not look as above for allelic studies. In the example below, the count matrices shown might look as follows before they are ready to be made into an MPRASet object.

E <- 2 # Number of elements
B <- 3 # Number of barcodes
s <- 4 # Total number of samples
nalleles <- 2 # Number of alleles

set.seed(434)
rna <- matrix(rpois(E*B*s*nalleles, lambda = 30), nrow = E*B*nalleles, ncol = s)
dna <- matrix(rpois(E*B*s*nalleles, lambda = 30), nrow = E*B*nalleles, ncol = s)

rn <- expand.grid(barcode = seq_len(B), allele = seq_len(nalleles), elem = seq_len(E))
rn <- paste0("barcode", rn$barcode, "_elem", rn$elem, "_allele", rn$allele)
cn <- paste0("sample", seq_len(s))

rownames(rna) <- rn
rownames(dna) <- rn
colnames(rna) <- cn
colnames(dna) <- cn

rna
##                        sample1 sample2 sample3 sample4
## barcode1_elem1_allele1      31      26      36      18
## barcode2_elem1_allele1      33      28      29      32
## barcode3_elem1_allele1      22      43      25      41
## barcode1_elem1_allele2      28      27      31      20
## barcode2_elem1_allele2      32      31      32      21
## barcode3_elem1_allele2      35      36      28      38
## barcode1_elem2_allele1      30      34      26      34
## barcode2_elem2_allele1      37      27      37      29
## barcode3_elem2_allele1      34      30      20      26
## barcode1_elem2_allele2      28      27      34      26
## barcode2_elem2_allele2      31      32      28      36
## barcode3_elem2_allele2      24      28      34      31
dna
##                        sample1 sample2 sample3 sample4
## barcode1_elem1_allele1      27      33      28      30
## barcode2_elem1_allele1      28      33      34      30
## barcode3_elem1_allele1      28      34      30      37
## barcode1_elem1_allele2      37      33      24      29
## barcode2_elem1_allele2      38      22      27      29
## barcode3_elem1_allele2      26      23      30      29
## barcode1_elem2_allele1      40      26      26      22
## barcode2_elem2_allele1      25      29      38      28
## barcode3_elem2_allele1      30      35      33      38
## barcode1_elem2_allele2      29      30      25      30
## barcode2_elem2_allele2      31      24      29      32
## barcode3_elem2_allele2      37      29      34      36

Most often with allelic studies, we will want to aggregate counts over barcodes to have summarized counts for each PRE-allele combination in each sample as below:

agg_output <- lapply(seq_len(E), function(elem_id) {
    pattern1 <- paste0(paste0("elem", elem_id), "_allele1")
    bool_rna_allele1 <- grepl(pattern1, rownames(rna))
    pattern2 <- paste0(paste0("elem", elem_id), "_allele2")
    bool_rna_allele2 <- grepl(pattern2, rownames(rna))
    agg_rna <- c(
        colSums(rna[bool_rna_allele1,]),
        colSums(rna[bool_rna_allele2,])
    )
    names(agg_rna) <- paste0(rep(c("allele1", "allele2"), each = s), "_", names(agg_rna))
    bool_dna_allele1 <- grepl(pattern1, rownames(dna))
    bool_dna_allele2 <- grepl(pattern2, rownames(dna))
    agg_dna <- c(
        colSums(dna[bool_dna_allele1,]),
        colSums(dna[bool_dna_allele2,])
    )
    names(agg_dna) <- paste0(rep(c("allele1", "allele2"), each = s), "_", names(agg_dna))
    list(rna = agg_rna, dna = agg_dna)
})
agg_rna <- do.call(rbind, lapply(agg_output, "[[", "rna"))
agg_dna <- do.call(rbind, lapply(agg_output, "[[", "dna"))
eid <- paste0("elem", seq_len(E))
rownames(agg_rna) <- eid
rownames(agg_dna) <- eid
eseq <- replicate(E, paste(sample(c("A", "T", "C", "G"), 10, replace = TRUE), collapse = ""))

With the relevant information defined, we can use the MPRASet constructor function as in the first example:

mpraset_example2 <- MPRASet(DNA = agg_dna, RNA = agg_rna, eid = eid, eseq = eseq, barcode = NULL)
mpraset_example2
## class: MPRASet 
## dim: 2 8 
## metadata(0):
## assays(2): DNA RNA
## rownames(2): elem1 elem2
## rowData names(2): eid eseq
## colnames(8): allele1_sample1 allele1_sample2 ... allele2_sample3
##   allele2_sample4
## colData names(0):
## No barcodes present

4 Analysis

4.1 Tissue comparison

While the above section demonstrated how to create MPRASet objects, we will use preconstructed objects containing data from a comparison of episomal and lentiviral versions of MPRA (Inoue et al. 2017).

data(mpraSetExample)

We create the design matrix with an indicator for the episomal (mutant integrase) samples and fit the precision-weighted linear model with mpralm. In MPRA experiments, activity measures are quantified as the log2 ratio of RNA counts over DNA counts. When there is barcode level information (as in this experiment), there are various ways to summarize information over barcodes to compute the final element- and sample-specific log ratios that are used for subsequent statistical modeling.

We have specified aggregate = "mean" to indicate that the element- and sample-specific log ratios will be computed by first computing the log ratio of RNA counts over DNA counts for each barcode, then taking the mean over barcodes in a particular element and sample. This is termed the “average estimator” in (Myint et al. 2019). In contrast, specifying aggregate = "sum" would indicate the use of the “aggregate estimator”, in which counts are first summed over barcodes to create total RNA and DNA counts, and the log ratio activity measure is the ratio of these total counts.

We have specifed normalize = TRUE to perform total count normalization on the RNA and DNA libraries. This scales all libraries to have a common size of 10 million reads.

Because this experiment looks at a set of PREs in two different cellular conditions, the different samples (columns of the MPRASet object) are independent. Thus we specify model_type = "indep_groups" to perform an unpaired analysis. In contrast, if we were performing an allele comparison, we would specify model_type = "corr_groups" to performed a paired analysis (indicating that different columns of the MPRASet object are linked).

Finally, we specify plot = TRUE to plot the relationship between log ratio variability versus element copy number.

design <- data.frame(intcpt = 1, episomal = grepl("MT", colnames(mpraSetExample)))
mpralm_fit <- mpralm(object = mpraSetExample, design = design, aggregate = "mean", normalize = TRUE, model_type = "indep_groups", plot = TRUE)

The resulting fit object can be used with topTable from the limma package.

toptab <- topTable(mpralm_fit, coef = 2, number = Inf)
toptab6 <- head(toptab)

Because the element codes are rather long for this dataset, we do some tricks to print the top differential elements:

rownames(toptab6)
## [1] "C:SLEA_hg18:chr2:210861483-210861650|5:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;26:V_HNF1_C:AGTTAATGATTAACCAA;45:V_HNF1_C:AGTTAATGATTAACCAA;64:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;85:V_HNF1_C:AGTTAATGATTAACCAA;104:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;125:V_HNF1_C:AGTTAATGATTAACCAA;144:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC"      
## [2] "R:EP300-NoMod_chr9:12814543-12814714_[chr9:12814543-12814714]"                                                                                                                                                                                                                                                       
## [3] "C:SLEA_hg18:chr2:210861483-210861650|4:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;25:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;46:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;67:V_HNF1_C:AGTTAATGATTAACCAA;86:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;107:V_HNF1_C:AGTTAATGATTAACCAA;126:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;147:V_HNF1_C:AGTTAATGATTAACCAA"
## [4] "C:SLEA_hg18:chr2:210861483-210861650|6:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;27:V_HNF1_C:AGTTAATGATTAACCAA;46:V_HNF1_C:AGTTAATGATTAACCAA;65:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC;86:V_HNF1_C:AGTTAATGATTAACCAA;105:V_HNF1_C:AGTTAATGATTAACCAA;124:V_HNF1_C:AGTTAATGATTAACCAA;143:V_AHRARNT_02:GGGGATCGCGTGCCAGCCC"            
## [5] "R:FOXA1_FOXA2-ChMod_chr1:38000211-38000353_[chr1:38000196-38000367]"                                                                                                                                                                                                                                                 
## [6] "R:EP300-NoMod_chr16:3070958-3071129_[chr16:3070958-3071129]"
rownames(toptab6) <- NULL
toptab6
##        logFC  AveExpr         t      P.Value    adj.P.Val        B
## 1 -1.2076627 1.284363 -26.99214 1.450176e-36 3.538429e-33 73.05825
## 2 -1.7316415 1.668562 -23.50738 4.357972e-33 5.316726e-30 64.87553
## 3 -1.1603766 1.497371 -23.11093 1.150733e-32 9.359292e-30 64.12459
## 4 -0.9560848 1.016586 -18.57491 2.136470e-27 1.303247e-24 52.09394
## 5 -1.1713718 1.730330 -18.25135 5.494116e-27 2.681128e-24 51.09648
## 6 -0.7639104 0.733524 -17.71599 2.689472e-26 1.093719e-23 49.54888

4.2 Allelic comparison

We will also demonstrate an allelic comparison analysis using the data in (Tewhey et al. 2016). In this study, the investigators compare reference and alternate allele versions of sequences containing SNPs. These sequences were believed to be eQTLs based on previous work.

We create a design matrix with an indicator for whether the counts come from the “B” allele (as opposed to the “A” allele). We also create an integer block_vector to indicate the actual sample that each column in the DNA and RNA count matrices comes from. In this case, columns 1 and 6 of the DNA and RNA count matrices come from sample 1 (the A and B alleles measured in transfection 1), columns 2 and 7 from sample 2, and so on. This information is needed to accurately model the within-sample correlation for a paired analysis.

data(mpraSetAllelicExample)

design <- data.frame(intcpt = 1, alleleB = grepl("allele_B", colnames(mpraSetAllelicExample)))
block_vector <- rep(1:5, 2)
mpralm_allele_fit <- mpralm(object = mpraSetAllelicExample, design = design, aggregate = "none", normalize = TRUE, block = block_vector, model_type = "corr_groups", plot = TRUE)
## Warning in regularize.values(x, y, ties, missing(ties), na.rm = na.rm):
## collapsing to unique 'x' values
## Warning: Zero sample variances detected, have been offset away from zero

toptab_allele <- topTable(mpralm_allele_fit, coef = 2, number = Inf)
head(toptab_allele)
##                 logFC  AveExpr         t      P.Value    adj.P.Val        B
## rs11080327   3.113780 1.694106  79.29544 4.261359e-16 1.596305e-12 26.21112
## rs9661285   -2.109272 1.299262 -49.94253 5.760237e-14 1.078892e-10 22.28452
## rs71535706  -1.298922 1.201426 -37.78349 1.104021e-12 1.378555e-09 19.54272
## rs7257930    1.705642 1.851310  35.15861 2.361939e-12 2.034619e-09 18.85088
## rs112372623 -1.829474 2.458348 -34.69678 2.715722e-12 2.034619e-09 18.72500
## rs2191501    1.891782 1.436918  32.00768 6.359293e-12 3.716120e-09 18.01258

5 Session Info

R version 4.2.1 (2022-06-23) Platform: x86_64-pc-linux-gnu (64-bit) Running under: Ubuntu 20.04.5 LTS

Matrix products: default BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so

locale: [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_GB 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] stats4 stats graphics grDevices utils datasets methods
[8] base

other attached packages: [1] mpra_1.20.0 limma_3.54.0
[3] SummarizedExperiment_1.28.0 Biobase_2.58.0
[5] GenomicRanges_1.50.0 GenomeInfoDb_1.34.0
[7] IRanges_2.32.0 S4Vectors_0.36.0
[9] MatrixGenerics_1.10.0 matrixStats_0.62.0
[11] BiocGenerics_0.44.0 BiocStyle_2.26.0

loaded via a namespace (and not attached): [1] Rcpp_1.0.9 highr_0.9 bslib_0.4.0
[4] compiler_4.2.1 BiocManager_1.30.19 jquerylib_0.1.4
[7] XVector_0.38.0 bitops_1.0-7 tools_4.2.1
[10] zlibbioc_1.44.0 statmod_1.4.37 digest_0.6.30
[13] lifecycle_1.0.3 lattice_0.20-45 jsonlite_1.8.3
[16] evaluate_0.17 rlang_1.0.6 Matrix_1.5-1
[19] DelayedArray_0.24.0 cli_3.4.1 magick_2.7.3
[22] yaml_2.3.6 xfun_0.34 fastmap_1.1.0
[25] GenomeInfoDbData_1.2.9 stringr_1.4.1 knitr_1.40
[28] sass_0.4.2 grid_4.2.1 R6_2.5.1
[31] rmarkdown_2.17 bookdown_0.29 farver_2.1.1
[34] magrittr_2.0.3 scales_1.2.1 htmltools_0.5.3
[37] colorspace_2.0-3 stringi_1.7.8 munsell_0.5.0
[40] RCurl_1.98-1.9 cachem_1.0.6

References

Inoue, Fumitaka, Martin Kircher, Beth Martin, Gregory M Cooper, Daniela M Witten, Michael T McManus, Nadav Ahituv, and Jay Shendure. 2017. “A Systematic Comparison Reveals Substantial Differences in Chromosomal Versus Episomal Encoding of Enhancer Activity.” Genome Research 27: 38–52. https://doi.org/10.1101/gr.212092.116.

Law, Charity W, Yunshun Chen, Wei Shi, and Gordon K Smyth. 2014. “Voom: Precision Weights Unlock Linear Model Analysis Tools for RNA-seq Read Counts.” Genome Biology 15: R29. https://doi.org/10.1186/gb-2014-15-2-r29.

Myint, Leslie, Dimitrios G Avramopoulos, Loyal A Goff, and Kasper D Hansen. 2019. “Linear Models Enable Powerful Differential Activity Analysis in Massively Parallel Reporter Assays.” BMC Genomics 20: 209. https://doi.org/10.1186/s12864-019-5556-x.

Tewhey, Ryan, Dylan Kotliar, Daniel S Park, Brandon Liu, Sarah Winnicki, Steven K Reilly, Kristian G Andersen, et al. 2016. “Direct Identification of Hundreds of Expression-Modulating Variants Using a Multiplexed Reporter Assay.” Cell 165: 1519–29. https://doi.org/10.1016/j.cell.2016.04.027.