Example Workflow For Processing a Single Pooled Screen

Russell Bainer

2023-10-24

Example Workflow For Processing a Single Screen

This is an example workflow for processing a pooled screening eperiment using the provided sample data. See the various manpages for additional visualization options and algorithmic details.

Note that what follows describes a very basic analysis. If you are considering integrating the results of many different screen contrasts or even different experiments and/or technologies, refer to the Advanced Screen Analysis: Contrast Comparisons vignette.

Load dependencies and data

suppressPackageStartupMessages(library(Biobase))
suppressPackageStartupMessages(library(limma))
suppressPackageStartupMessages(library(gCrisprTools))

data("es", package = "gCrisprTools")
data("ann", package = "gCrisprTools")
data("aln", package = "gCrisprTools")
knitr::opts_chunk$set(message = FALSE, fig.width = 8, fig.height = 8, warning = FALSE)

Make a sample key, structured as a factor with control samples in the first level

sk <- relevel(as.factor(pData(es)$TREATMENT_NAME), "ControlReference")
names(sk) <- row.names(pData(es))

Generate a contrast of interest using voom/limma; pairing replicates is a good idea if that information is available.

design <- model.matrix(~ 0 + REPLICATE_POOL + TREATMENT_NAME, pData(es))
colnames(design) <- gsub('TREATMENT_NAME', '', colnames(design))
contrasts <-makeContrasts(DeathExpansion - ControlExpansion, levels = design)

Optionally, trim of trace reads from the unnormalized object (see man page for details)

es <- ct.filterReads(es, trim = 1000, sampleKey = sk)

plot of chunk unnamed-chunk-3

Normalize, convert to a voom object, and generate a contrast

es <- ct.normalizeGuides(es, method = "scale", plot.it = TRUE) #See man page for other options
vm <- voom(exprs(es), design)

fit <- lmFit(vm, design)
fit <- contrasts.fit(fit, contrasts)
fit <- eBayes(fit)

Edit the annotation file if you used ct.filterReads above

ann <- ct.prepareAnnotation(ann, fit, controls = "NoTarget")

Summarize gRNA signals to identify target genes of interest

resultsDF <-
  ct.generateResults(
    fit,
    annotation = ann,
    RRAalphaCutoff = 0.1,
    permutations = 1000,
    scoring = "combined", 
    permutation.seed = 2
  )

Alternative Annotations

In some cases, reagents might target multiple known elements (e.g., gRNAs in a CRISPRi library that target multiple promoters of the same gene). In such cases, you can specify this via the alt.annotation argument to ct.generateResults(). Alternative annotations are supplied as a list of character vectors named for the reagents.

# Create random alternative target associations 

altann <- sapply(ann$ID, 
                 function(x){
                   out <- as.character(ann$geneSymbol)[ann$ID %in% x]
                   if(runif(1) < 0.01){out <- c(out, sample(as.character(ann$geneSymbol), size = 1))}
                   return(out)
                 }, simplify = FALSE)

resultsDF <-
  ct.generateResults(
    fit,
    annotation = ann,
    RRAalphaCutoff = 0.1,
    permutations = 1000,
    scoring = "combined", 
    alt.annotation = altann,
    permutation.seed = 2
  )

Optionally, just load an example results object for testing purposes (trimming out reads as necessary)

data("fit", package = "gCrisprTools")
data("resultsDF", package = "gCrisprTools")

fit <- fit[(row.names(fit) %in% row.names(ann)),]
resultsDF <- resultsDF[(row.names(resultsDF) %in% row.names(ann)),]

targetResultDF <- ct.simpleResult(resultsDF) #For a simpler target-level result object

Quality Control

gCrisprTools contains a variety of pooled screen-specific quality control and visualization tools (see man pages for details):

ct.alignmentChart(aln, sk)

plot of chunk unnamed-chunk-9

ct.rawCountDensities(es, sk)

plot of chunk unnamed-chunk-9

Visualize gRNA abundance distributions

ct.gRNARankByReplicate(es, sk) 

plot of chunk unnamed-chunk-10

ct.gRNARankByReplicate(es, sk, annotation = ann, geneSymb = "NoTarget")  #Show locations of NTC gRNAs

plot of chunk unnamed-chunk-10

Visualize control guide behavior across conditions

ct.viewControls(es, ann, sk, normalize = FALSE)

plot of chunk unnamed-chunk-11

ct.viewControls(es, ann, sk, normalize = TRUE)

plot of chunk unnamed-chunk-11

Visualize GC bias across samples, or within an experimental contrast

ct.GCbias(es, ann, sk)

plot of chunk unnamed-chunk-12

ct.GCbias(fit, ann, sk)

plot of chunk unnamed-chunk-12

View most variable gRNAs/Genes (as % of sequencing library)

ct.stackGuides(es,
               sk,
               plotType = "gRNA",
               annotation = ann,
               nguides = 40)

plot of chunk unnamed-chunk-13

ct.stackGuides(es, 
               sk, 
               plotType = "Target", 
               annotation = ann)

plot of chunk unnamed-chunk-14

ct.stackGuides(es,
               sk,
               plotType = "Target",
               annotation = ann,
               subset = names(sk)[grep('Expansion', sk)])

plot of chunk unnamed-chunk-15

View a CDF of genes/guides

ct.guideCDF(es, sk, plotType = "gRNA")

plot of chunk unnamed-chunk-16

ct.guideCDF(es, sk, plotType = "Target", annotation = ann)

plot of chunk unnamed-chunk-16

Target-Level Visualization and Analysis

View the overall enrichment and depletion trends identified in the screen:

ct.contrastBarchart(resultsDF)

plot of chunk unnamed-chunk-17

View top enriched/depleted candidates

ct.topTargets(fit,
              resultsDF,
              ann,
              targets = 10,
              enrich = TRUE)

plot of chunk unnamed-chunk-18

ct.topTargets(fit,
              resultsDF,
              ann,
              targets = 10,
              enrich = FALSE)

plot of chunk unnamed-chunk-18

View the behavior of reagents targeting a particular gene of interest

ct.viewGuides("Target1633", fit, ann)

plot of chunk unnamed-chunk-19

ct.gRNARankByReplicate(es, sk, annotation = ann, geneSymb = "Target1633")

plot of chunk unnamed-chunk-19

Observe the effects detected for sets of targets within a screen contrast

ct.signalSummary(resultsDF,
                 targets = list('TargetSetA' = c(sample(unique(resultsDF$geneSymbol), 3)),
                                'TargetSetB' = c(sample(unique(resultsDF$geneSymbol), 2))))

plot of chunk unnamed-chunk-20

You could test a known gene set for enrichment within target candidates:

data("essential.genes", package = "gCrisprTools")
ct.targetSetEnrichment(resultsDF, essential.genes)
## $targets
##  [1] "8877"   "26999"  "10061"  "4050"   "11054"  "79586"  "51171"  "5260"  
##  [9] "127544" "8645"   "3702"   "4621"   "2000"   "7248"   "8518"   "844"   
## [17] "6950"   "845"    "10237"  "11044"  "462"    "54739"  "1856"   "7077"  
## [25] "10645"  "3759"   "3455"   "140465" "80273"  "87178"  "84913"  "5313"  
## [33] "9592"   "6696"   "6426"   "26119"  "9510"   "8993"   "4843"   "6383"  
## [41] "5133"   "1289"   "5395"   "10478"  "4129"   "4897"   "10538"  "1737"  
## [49] "715"    "5739"   "3112"   "5880"   "567"    "3651"   "9971"   "10856" 
## [57] "7096"   "6258"   "1642"   "16"     "7127"   "3764"   "4316"   "6518"  
## [65] "3489"   "1984"   "290"    "1345"   "7475"   "5650"   "7298"   "1823"  
## [73] "633"    "657"    "9456"   "39"     "8476"   "5641"   "3659"   "6881"  
## [81] "10840"  "79109"  "6820"   "597"    "8630"   "10134"  "4778"   "967"   
## [89] "6370"   "4504"   "4482"   "6150"   "26275"  "55109"  "9935"   "4209"  
## [97] "51148"  "134429" "1213"  
## 
## $P.values
##      Cutoff Sig Hypergeometric_P
## [1,]  0e+00   0      1.000000000
## [2,]  1e-05   0      1.000000000
## [3,]  1e-04   0      1.000000000
## [4,]  1e-03   0      1.000000000
## [5,]  1e-02   0      1.000000000
## [6,]  1e-01   5      0.810437224
## [7,]  5e-01  13      0.434952850
## [8,]  1e+00  99      0.009988581
## 
## $Q.values
##      Cutoff Sig Hypergeometric_P
## [1,]  0e+00   0      1.000000000
## [2,]  1e-05   0      1.000000000
## [3,]  1e-04   0      1.000000000
## [4,]  1e-03   0      1.000000000
## [5,]  1e-02   0      1.000000000
## [6,]  1e-01   0      1.000000000
## [7,]  5e-01   0      1.000000000
## [8,]  1e+00  99      0.009988581

Or optionally add a visualization:

ROC <- ct.ROC(resultsDF, essential.genes, direction = "deplete")

plot of chunk rocprc

PRC <- ct.PRC(resultsDF, essential.genes, direction = "deplete")

plot of chunk rocprc

show(ROC) # show(PRC) is equivalent for the PRC analysis
## $sensitivity
##   [1] 0.02020202 0.02020202 0.02020202 0.02020202 0.03030303 0.03030303
##   [7] 0.03030303 0.03030303 0.03030303 0.03030303 0.04040404 0.04040404
##  [13] 0.04040404 0.04040404 0.04040404 0.04040404 0.04040404 0.04040404
##  [19] 0.04040404 0.05050505 0.05050505 0.05050505 0.05050505 0.05050505
##  [25] 0.06060606 0.06060606 0.07070707 0.07070707 0.07070707 0.07070707
##  [31] 0.07070707 0.07070707 0.07070707 0.07070707 0.07070707 0.08080808
##  [37] 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808
##  [43] 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808
##  [49] 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808
##  [55] 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808 0.08080808
##  [61] 0.08080808 0.08080808 0.08080808 0.09090909 0.09090909 0.09090909
##  [67] 0.09090909 0.10101010 0.10101010 0.10101010 0.10101010 0.10101010
##  [73] 0.10101010 0.10101010 0.10101010 0.10101010 0.10101010 0.10101010
##  [79] 0.10101010 0.10101010 0.10101010 0.11111111 0.11111111 0.11111111
##  [85] 0.11111111 0.11111111 0.11111111 0.11111111 0.11111111 0.11111111
##  [91] 0.11111111 0.11111111 0.11111111 0.11111111 0.11111111 0.11111111
##  [97] 0.11111111 0.11111111 0.11111111 0.12121212 0.12121212 0.12121212
## [103] 0.12121212 0.12121212 0.12121212 0.12121212 0.12121212 0.12121212
## [109] 0.12121212 0.12121212 0.13131313 0.13131313 0.13131313 0.13131313
## [115] 0.13131313 0.13131313 0.13131313 0.13131313 0.13131313 0.14141414
## [121] 0.14141414 0.14141414 0.14141414 0.14141414 0.14141414 0.14141414
## [127] 0.15151515 0.15151515 0.15151515 0.15151515 0.15151515 0.15151515
## [133] 0.15151515 0.15151515 0.16161616 0.16161616 0.16161616 0.16161616
## [139] 0.16161616 0.16161616 0.16161616 0.17171717 0.17171717 0.17171717
## [145] 0.17171717 0.17171717 0.17171717 0.18181818 0.18181818 0.18181818
## [151] 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818
## [157] 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818
## [163] 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818
## [169] 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818
## [175] 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818 0.18181818
## [181] 0.19191919 0.19191919 0.19191919 0.19191919 0.19191919 0.19191919
## [187] 0.19191919 0.20202020 0.20202020 0.20202020 0.20202020 0.20202020
## [193] 0.20202020 0.20202020 0.20202020 0.20202020 0.20202020 0.20202020
## [199] 0.20202020 0.20202020 0.20202020 0.20202020 0.20202020 0.20202020
## [205] 0.20202020 0.21212121 0.21212121 0.21212121 0.21212121 0.21212121
## [211] 0.21212121 0.21212121 0.21212121 0.21212121 0.21212121 0.21212121
## [217] 0.21212121 0.21212121 0.21212121 0.21212121 0.21212121 0.21212121
## [223] 0.21212121 0.21212121 0.21212121 0.22222222 0.22222222 0.23232323
## [229] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [235] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [241] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [247] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [253] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [259] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323
## [265] 0.23232323 0.23232323 0.23232323 0.23232323 0.23232323 0.24242424
## [271] 0.24242424 0.24242424 0.24242424 0.24242424 0.24242424 0.24242424
## [277] 0.25252525 0.25252525 0.25252525 0.25252525 0.25252525 0.25252525
## [283] 0.25252525 0.26262626 0.26262626 0.26262626 0.26262626 0.26262626
## [289] 0.26262626 0.26262626 0.26262626 0.26262626 0.26262626 0.26262626
## [295] 0.26262626 0.26262626 0.26262626 0.26262626 0.26262626 0.27272727
## [301] 0.27272727 0.27272727 0.27272727 0.27272727 0.27272727 0.27272727
## [307] 0.27272727 0.27272727 0.27272727 0.27272727 0.27272727 0.27272727
## [313] 0.27272727 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828
## [319] 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828
## [325] 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828
## [331] 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828 0.28282828
## [337] 0.28282828 0.28282828 0.28282828 0.28282828 0.29292929 0.29292929
## [343] 0.29292929 0.29292929 0.29292929 0.29292929 0.29292929 0.29292929
## [349] 0.29292929 0.29292929 0.29292929 0.30303030 0.30303030 0.31313131
## [355] 0.31313131 0.31313131 0.31313131 0.31313131 0.31313131 0.31313131
## [361] 0.31313131 0.31313131 0.32323232 0.32323232 0.32323232 0.32323232
## [367] 0.32323232 0.32323232 0.32323232 0.32323232 0.32323232 0.32323232
## [373] 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333
## [379] 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333
## [385] 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333 0.33333333
## [391] 0.33333333 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [397] 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [403] 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [409] 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [415] 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [421] 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434 0.34343434
## [427] 0.34343434 0.35353535 0.35353535 0.35353535 0.35353535 0.35353535
## [433] 0.35353535 0.36363636 0.36363636 0.36363636 0.36363636 0.36363636
## [439] 0.38383838 0.38383838 0.38383838 0.39393939 0.39393939 0.39393939
## [445] 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939
## [451] 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939
## [457] 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939 0.39393939
## [463] 0.39393939 0.39393939 0.40404040 0.40404040 0.40404040 0.41414141
## [469] 0.42424242 0.42424242 0.42424242 0.42424242 0.42424242 0.42424242
## [475] 0.42424242 0.42424242 0.42424242 0.42424242 0.42424242 0.42424242
## [481] 0.42424242 0.42424242 0.42424242 0.42424242 0.42424242 0.43434343
## [487] 0.44444444 0.45454545 0.45454545 0.45454545 0.45454545 0.45454545
## [493] 0.46464646 0.46464646 0.46464646 0.46464646 0.46464646 0.46464646
## [499] 0.47474747 0.48484848 0.48484848 0.48484848 0.48484848 0.48484848
## [505] 0.48484848 0.48484848 0.48484848 0.48484848 0.48484848 0.48484848
## [511] 0.48484848 0.48484848 0.48484848 0.48484848 0.48484848 0.48484848
## [517] 0.48484848 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [523] 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [529] 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [535] 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [541] 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [547] 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949 0.49494949
## [553] 0.49494949 0.49494949 0.50505051 0.50505051 0.50505051 0.50505051
## [559] 0.50505051 0.51515152 0.51515152 0.51515152 0.51515152 0.52525253
## [565] 0.52525253 0.52525253 0.52525253 0.52525253 0.52525253 0.52525253
## [571] 0.52525253 0.52525253 0.52525253 0.53535354 0.53535354 0.53535354
## [577] 0.53535354 0.54545455 0.54545455 0.54545455 0.54545455 0.54545455
## [583] 0.54545455 0.54545455 0.54545455 0.54545455 0.54545455 0.54545455
## [589] 0.54545455 0.54545455 0.54545455 0.55555556 0.55555556 0.55555556
## [595] 0.55555556 0.55555556 0.55555556 0.55555556 0.55555556 0.56565657
## [601] 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657
## [607] 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657
## [613] 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657
## [619] 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657 0.56565657
## [625] 0.56565657 0.57575758 0.57575758 0.57575758 0.57575758 0.57575758
## [631] 0.57575758 0.57575758 0.57575758 0.57575758 0.57575758 0.57575758
## [637] 0.57575758 0.57575758 0.57575758 0.57575758 0.57575758 0.57575758
## [643] 0.57575758 0.57575758 0.58585859 0.58585859 0.58585859 0.58585859
## [649] 0.58585859 0.58585859 0.58585859 0.58585859 0.58585859 0.58585859
## [655] 0.58585859 0.58585859 0.58585859 0.59595960 0.59595960 0.59595960
## [661] 0.59595960 0.60606061 0.60606061 0.60606061 0.60606061 0.63636364
## [667] 0.63636364 0.63636364 0.63636364 0.63636364 0.63636364 0.63636364
## [673] 0.63636364 0.63636364 0.63636364 0.63636364 0.63636364 0.63636364
## [679] 0.64646465 0.64646465 0.64646465 0.64646465 0.64646465 0.64646465
## [685] 0.64646465 0.64646465 0.64646465 0.64646465 0.64646465 0.64646465
## [691] 0.66666667 0.66666667 0.66666667 0.66666667 0.67676768 0.67676768
## [697] 0.67676768 0.67676768 0.67676768 0.67676768 0.67676768 0.68686869
## [703] 0.68686869 0.68686869 0.68686869 0.68686869 0.68686869 0.68686869
## [709] 0.68686869 0.68686869 0.68686869 0.68686869 0.68686869 0.68686869
## [715] 0.68686869 0.68686869 0.68686869 0.69696970 0.69696970 0.69696970
## [721] 0.70707071 0.70707071 0.70707071 0.71717172 0.71717172 0.71717172
## [727] 0.71717172 0.71717172 0.71717172 0.71717172 0.71717172 0.71717172
## [733] 0.71717172 0.71717172 0.71717172 0.71717172 0.71717172 0.71717172
## [739] 0.71717172 0.71717172 0.71717172 0.71717172 0.72727273 0.72727273
## [745] 0.72727273 0.72727273 0.73737374 0.73737374 0.73737374 0.73737374
## [751] 0.73737374 0.73737374 0.73737374 0.73737374 0.73737374 0.73737374
## [757] 0.73737374 0.73737374 0.73737374 0.73737374 0.73737374 0.73737374
## [763] 0.74747475 0.74747475 0.74747475 0.75757576 0.75757576 0.75757576
## [769] 0.75757576 0.75757576 0.75757576 0.75757576 0.75757576 0.75757576
## [775] 0.75757576 0.75757576 0.75757576 0.76767677 0.76767677 0.76767677
## [781] 0.76767677 0.76767677 0.77777778 0.77777778 0.77777778 0.77777778
## [787] 0.77777778 0.77777778 0.77777778 0.78787879 0.78787879 0.79797980
## [793] 0.81818182 0.82828283 0.82828283 0.82828283 0.82828283 0.82828283
## [799] 0.82828283 0.83838384 0.83838384 0.83838384 0.83838384 0.83838384
## [805] 0.83838384 0.83838384 0.83838384 0.83838384 0.83838384 0.83838384
## [811] 0.83838384 0.84848485 0.84848485 0.84848485 0.84848485 0.84848485
## [817] 0.84848485 0.84848485 0.84848485 0.84848485 0.84848485 0.84848485
## [823] 0.84848485 0.84848485 0.84848485 0.84848485 0.84848485 0.84848485
## [829] 0.85858586 0.85858586 0.85858586 0.85858586 0.85858586 0.85858586
## [835] 0.85858586 0.85858586 1.00000000
## 
## $specificity
##   [1] 0.9871858 0.9817644 0.9788073 0.9748645 0.9728931 0.9699359 0.9689502
##   [8] 0.9674717 0.9645145 0.9640217 0.9635288 0.9615574 0.9595860 0.9590931
##  [15] 0.9586003 0.9556432 0.9546575 0.9541646 0.9536718 0.9521932 0.9502218
##  [22] 0.9487432 0.9462790 0.9457861 0.9452932 0.9428290 0.9418433 0.9398719
##  [29] 0.9393790 0.9388862 0.9383933 0.9379004 0.9374076 0.9364219 0.9349433
##  [36] 0.9329719 0.9324791 0.9314933 0.9310005 0.9305076 0.9295219 0.9290291
##  [43] 0.9270577 0.9255791 0.9241005 0.9236077 0.9226220 0.9216363 0.9206506
##  [50] 0.9176934 0.9157220 0.9152292 0.9137506 0.9127649 0.9107935 0.9093149
##  [57] 0.9078364 0.9073435 0.9063578 0.9053721 0.9048793 0.9038935 0.9029078
##  [64] 0.9024150 0.9014293 0.8999507 0.8989650 0.8974864 0.8965007 0.8960079
##  [71] 0.8950222 0.8935436 0.8930508 0.8925579 0.8915722 0.8910793 0.8905865
##  [78] 0.8891079 0.8871365 0.8861508 0.8856580 0.8846723 0.8831937 0.8822080
##  [85] 0.8812223 0.8802366 0.8797437 0.8792509 0.8777723 0.8767866 0.8758009
##  [92] 0.8748152 0.8733366 0.8718581 0.8703795 0.8689009 0.8679152 0.8669295
##  [99] 0.8664367 0.8654510 0.8649581 0.8639724 0.8634795 0.8615081 0.8605224
## [106] 0.8590439 0.8585510 0.8580582 0.8575653 0.8560867 0.8546082 0.8541153
## [113] 0.8526368 0.8506654 0.8501725 0.8496796 0.8477082 0.8467225 0.8447511
## [120] 0.8447511 0.8442583 0.8422868 0.8408083 0.8398226 0.8393297 0.8378512
## [127] 0.8368655 0.8334155 0.8324298 0.8319369 0.8299655 0.8279941 0.8265155
## [134] 0.8260227 0.8250370 0.8230655 0.8206013 0.8186299 0.8176442 0.8171513
## [141] 0.8161656 0.8146870 0.8132085 0.8127156 0.8117299 0.8102514 0.8097585
## [148] 0.8097585 0.8092656 0.8087728 0.8077871 0.8068014 0.8063085 0.8058157
## [155] 0.8048300 0.8043371 0.8038443 0.8033514 0.8028586 0.8023657 0.8018728
## [162] 0.8008871 0.8003943 0.7989157 0.7984229 0.7979300 0.7964515 0.7954657
## [169] 0.7944800 0.7934943 0.7930015 0.7925086 0.7920158 0.7905372 0.7900444
## [176] 0.7885658 0.7875801 0.7870872 0.7856087 0.7846230 0.7846230 0.7841301
## [183] 0.7821587 0.7811730 0.7801873 0.7792016 0.7787087 0.7782159 0.7777230
## [190] 0.7772302 0.7767373 0.7752587 0.7742730 0.7732873 0.7713159 0.7703302
## [197] 0.7698374 0.7693445 0.7683588 0.7668802 0.7654017 0.7649088 0.7634303
## [204] 0.7624446 0.7619517 0.7619517 0.7614588 0.7589946 0.7580089 0.7575160
## [211] 0.7565303 0.7550517 0.7545589 0.7535732 0.7520946 0.7511089 0.7506161
## [218] 0.7496304 0.7486447 0.7466732 0.7456875 0.7451947 0.7442090 0.7437161
## [225] 0.7432233 0.7427304 0.7417447 0.7402661 0.7392804 0.7387876 0.7373090
## [232] 0.7353376 0.7333662 0.7328733 0.7323805 0.7313948 0.7309019 0.7279448
## [239] 0.7269591 0.7240020 0.7230163 0.7205520 0.7200591 0.7171020 0.7151306
## [246] 0.7136520 0.7126663 0.7121735 0.7116806 0.7102021 0.7087235 0.7072449
## [253] 0.7062592 0.7052735 0.7047807 0.7042878 0.7033021 0.7023164 0.7013307
## [260] 0.6998521 0.6993593 0.6983736 0.6973879 0.6964022 0.6959093 0.6934450
## [267] 0.6929522 0.6919665 0.6914736 0.6904879 0.6895022 0.6890094 0.6880237
## [274] 0.6865451 0.6850665 0.6845737 0.6835880 0.6826023 0.6816166 0.6806309
## [281] 0.6796451 0.6786594 0.6761952 0.6757023 0.6737309 0.6727452 0.6712666
## [288] 0.6688024 0.6683095 0.6673238 0.6658452 0.6653524 0.6638738 0.6628881
## [295] 0.6619024 0.6609167 0.6594381 0.6589453 0.6579596 0.6564810 0.6550025
## [302] 0.6545096 0.6530310 0.6525382 0.6515525 0.6500739 0.6495811 0.6485954
## [309] 0.6471168 0.6466240 0.6456382 0.6451454 0.6436668 0.6431740 0.6426811
## [316] 0.6421883 0.6416954 0.6407097 0.6402169 0.6387383 0.6382454 0.6362740
## [323] 0.6343026 0.6338098 0.6333169 0.6323312 0.6313455 0.6303598 0.6288812
## [330] 0.6283884 0.6274027 0.6264170 0.6254312 0.6244455 0.6239527 0.6224741
## [337] 0.6209956 0.6205027 0.6200099 0.6195170 0.6195170 0.6190241 0.6175456
## [344] 0.6165599 0.6155742 0.6145885 0.6126171 0.6121242 0.6111385 0.6106456
## [351] 0.6101528 0.6096599 0.6091671 0.6086742 0.6081814 0.6067028 0.6057171
## [358] 0.6042385 0.6027600 0.6022671 0.6017743 0.6012814 0.6002957 0.5988172
## [365] 0.5983243 0.5978314 0.5968457 0.5953672 0.5948743 0.5933958 0.5924101
## [372] 0.5904386 0.5894529 0.5879744 0.5869887 0.5864958 0.5860030 0.5855101
## [379] 0.5850172 0.5845244 0.5840315 0.5830458 0.5825530 0.5810744 0.5800887
## [386] 0.5795959 0.5786102 0.5766387 0.5761459 0.5751602 0.5731888 0.5722031
## [393] 0.5712173 0.5697388 0.5682602 0.5672745 0.5662888 0.5653031 0.5643174
## [400] 0.5628388 0.5623460 0.5618531 0.5613603 0.5598817 0.5593889 0.5588960
## [407] 0.5579103 0.5574174 0.5559389 0.5549532 0.5539675 0.5534746 0.5524889
## [414] 0.5519961 0.5515032 0.5505175 0.5500246 0.5490389 0.5485461 0.5480532
## [421] 0.5470675 0.5450961 0.5441104 0.5426318 0.5416461 0.5401676 0.5396747
## [428] 0.5372104 0.5357319 0.5352390 0.5347462 0.5332676 0.5327748 0.5317891
## [435] 0.5298176 0.5293248 0.5283391 0.5278462 0.5263677 0.5258748 0.5248891
## [442] 0.5239034 0.5229177 0.5219320 0.5209463 0.5204534 0.5194677 0.5189749
## [449] 0.5170034 0.5160177 0.5150320 0.5145392 0.5140463 0.5130606 0.5120749
## [456] 0.5091178 0.5081321 0.5076392 0.5056678 0.5051750 0.5046821 0.5032035
## [463] 0.5012321 0.5002464 0.5002464 0.4997536 0.4977822 0.4967965 0.4958107
## [470] 0.4948250 0.4938393 0.4923608 0.4908822 0.4898965 0.4894036 0.4889108
## [477] 0.4879251 0.4869394 0.4864465 0.4854608 0.4849680 0.4844751 0.4829966
## [484] 0.4825037 0.4810251 0.4805323 0.4800394 0.4795466 0.4785609 0.4775752
## [491] 0.4760966 0.4751109 0.4751109 0.4741252 0.4736323 0.4721538 0.4706752
## [498] 0.4696895 0.4682109 0.4677181 0.4657467 0.4652538 0.4647610 0.4637753
## [505] 0.4632824 0.4618038 0.4603253 0.4598324 0.4588467 0.4578610 0.4563825
## [512] 0.4558896 0.4553967 0.4544110 0.4539182 0.4534253 0.4519468 0.4514539
## [519] 0.4494825 0.4484968 0.4470182 0.4465254 0.4460325 0.4450468 0.4440611
## [526] 0.4435683 0.4420897 0.4415968 0.4411040 0.4401183 0.4391326 0.4376540
## [533] 0.4371612 0.4351897 0.4346969 0.4337112 0.4317398 0.4312469 0.4307541
## [540] 0.4302612 0.4287827 0.4282898 0.4277969 0.4273041 0.4268112 0.4263184
## [547] 0.4253327 0.4248398 0.4233613 0.4223756 0.4213898 0.4208970 0.4199113
## [554] 0.4194184 0.4184327 0.4174470 0.4169542 0.4159685 0.4144899 0.4135042
## [561] 0.4130113 0.4115328 0.4105471 0.4100542 0.4090685 0.4085757 0.4070971
## [568] 0.4066042 0.4061114 0.4051257 0.4046328 0.4031543 0.4026614 0.4016757
## [575] 0.4001971 0.3987186 0.3977329 0.3977329 0.3967472 0.3957615 0.3952686
## [582] 0.3937900 0.3923115 0.3918186 0.3898472 0.3888615 0.3878758 0.3868901
## [589] 0.3859044 0.3844258 0.3839330 0.3829473 0.3819616 0.3814687 0.3804830
## [596] 0.3785116 0.3780187 0.3775259 0.3770330 0.3760473 0.3745688 0.3735830
## [603] 0.3716116 0.3701331 0.3691474 0.3686545 0.3671759 0.3656974 0.3652045
## [610] 0.3622474 0.3602760 0.3592903 0.3587974 0.3568260 0.3558403 0.3528832
## [617] 0.3523903 0.3518975 0.3509118 0.3504189 0.3484475 0.3479547 0.3469690
## [624] 0.3454904 0.3440118 0.3425333 0.3410547 0.3405619 0.3400690 0.3395761
## [631] 0.3380976 0.3361262 0.3351405 0.3341548 0.3336619 0.3321833 0.3316905
## [638] 0.3302119 0.3297191 0.3292262 0.3277477 0.3262691 0.3238048 0.3228191
## [645] 0.3223263 0.3218334 0.3213406 0.3203549 0.3198620 0.3188763 0.3169049
## [652] 0.3164120 0.3159192 0.3149335 0.3144406 0.3139478 0.3129621 0.3119763
## [659] 0.3104978 0.3095121 0.3085264 0.3085264 0.3080335 0.3075407 0.3070478
## [666] 0.3065550 0.3060621 0.3035978 0.3026121 0.3016264 0.3006407 0.2986693
## [673] 0.2981764 0.2976836 0.2971907 0.2957122 0.2947265 0.2942336 0.2942336
## [680] 0.2932479 0.2917693 0.2912765 0.2907836 0.2897979 0.2888122 0.2873337
## [687] 0.2868408 0.2863480 0.2853622 0.2838837 0.2824051 0.2819123 0.2814194
## [694] 0.2804337 0.2789552 0.2779694 0.2774766 0.2769837 0.2755052 0.2745195
## [701] 0.2725481 0.2725481 0.2715623 0.2710695 0.2705766 0.2690981 0.2681124
## [708] 0.2661410 0.2646624 0.2621981 0.2612124 0.2607196 0.2602267 0.2577624
## [715] 0.2562839 0.2552982 0.2548053 0.2543125 0.2523411 0.2518482 0.2503696
## [722] 0.2498768 0.2488911 0.2474125 0.2469197 0.2459340 0.2449483 0.2444554
## [729] 0.2429768 0.2414983 0.2405126 0.2400197 0.2395269 0.2385412 0.2375554
## [736] 0.2365697 0.2350912 0.2336126 0.2326269 0.2321341 0.2316412 0.2301626
## [743] 0.2291769 0.2286841 0.2276984 0.2267127 0.2252341 0.2237555 0.2227698
## [750] 0.2217841 0.2207984 0.2198127 0.2193199 0.2188270 0.2183342 0.2173484
## [757] 0.2168556 0.2163627 0.2148842 0.2138985 0.2129128 0.2124199 0.2114342
## [764] 0.2089699 0.2074914 0.2069985 0.2065057 0.2055200 0.2045343 0.2035485
## [771] 0.2030557 0.2025628 0.2010843 0.2000986 0.1981272 0.1971414 0.1966486
## [778] 0.1956629 0.1951700 0.1941843 0.1927058 0.1917201 0.1902415 0.1892558
## [785] 0.1882701 0.1877772 0.1867915 0.1858058 0.1848201 0.1838344 0.1823558
## [792] 0.1813701 0.1793987 0.1793987 0.1784130 0.1779202 0.1774273 0.1759487
## [799] 0.1749630 0.1749630 0.1739773 0.1734845 0.1729916 0.1705274 0.1695416
## [806] 0.1690488 0.1660917 0.1651060 0.1636274 0.1616560 0.1611631 0.1601774
## [813] 0.1596846 0.1591917 0.1582060 0.1577132 0.1567275 0.1552489 0.1547560
## [820] 0.1542632 0.1527846 0.1498275 0.1473632 0.1458847 0.1448990 0.1424347
## [827] 0.1414490 0.1409561 0.1399704 0.1389847 0.1384919 0.1370133 0.1365205
## [834] 0.1350419 0.1335633 0.1315919 0.0000000
## 
## $AUC
## [1] 0.4059192
## 
## $targets
##  [1] "8877"   "26999"  "10061"  "4050"   "11054"  "79586"  "51171"  "5260"  
##  [9] "127544" "8645"   "3702"   "4621"   "2000"   "7248"   "8518"   "844"   
## [17] "6950"   "845"    "10237"  "11044"  "462"    "54739"  "1856"   "7077"  
## [25] "10645"  "3759"   "3455"   "140465" "80273"  "87178"  "84913"  "5313"  
## [33] "9592"   "6696"   "6426"   "26119"  "9510"   "8993"   "4843"   "6383"  
## [41] "5133"   "1289"   "5395"   "10478"  "4129"   "4897"   "10538"  "1737"  
## [49] "715"    "5739"   "3112"   "5880"   "567"    "3651"   "9971"   "10856" 
## [57] "7096"   "6258"   "1642"   "16"     "7127"   "3764"   "4316"   "6518"  
## [65] "3489"   "1984"   "290"    "1345"   "7475"   "5650"   "7298"   "1823"  
## [73] "633"    "657"    "9456"   "39"     "8476"   "5641"   "3659"   "6881"  
## [81] "10840"  "79109"  "6820"   "597"    "8630"   "10134"  "4778"   "967"   
## [89] "6370"   "4504"   "4482"   "6150"   "26275"  "55109"  "9935"   "4209"  
## [97] "51148"  "134429" "1213"  
## 
## $P.values
##      Cutoff Sig Hypergeometric_P
## [1,]  0e+00   2      0.355903161
## [2,]  1e-05   2      0.355903161
## [3,]  1e-04   2      0.355903161
## [4,]  1e-03   2      0.524176993
## [5,]  1e-02   4      0.460291014
## [6,]  1e-01  11      0.728237649
## [7,]  5e-01  35      0.966856290
## [8,]  1e+00  99      0.009988581
## 
## $Q.values
##      Cutoff Sig Hypergeometric_P
## [1,]  0e+00   2      0.355903161
## [2,]  1e-05   2      0.355903161
## [3,]  1e-04   2      0.355903161
## [4,]  1e-03   2      0.355903161
## [5,]  1e-02   2      0.355903161
## [6,]  1e-01   2      0.602972318
## [7,]  5e-01   7      0.378559005
## [8,]  1e+00  99      0.009988581

Or alternatively you could test for ontological enrichment within the depleted/enriched targets via the sparrow package:

#Create a geneset database using one of the many helper functions
genesetdb <- sparrow::getMSigGeneSetDb(collection = 'h', species = 'human', id.type = 'entrez')

ct.seas(resultsDF, gdb = genesetdb)
## $result
## SparrowResult (max FDR by collection set to 0.20%)
## --------------------------------------------------- 
##   collection method geneset_count sig_count sig_up sig_down
## 1          H  fgsea            50         1      0        1
# If you have a library that targets elements unevenly (e.g., variable numbers of 
# elements/promoters per genes), you can conform it via `sparrow::convertIdentifiers()`

genesetdb.GREAT <- sparrow::convertIdentifiers(genesetdb, 
                                               from = 'geneID', 
                                               to = 'geneSymbol', 
                                               xref = ann)
ct.seas(resultsDF, gdb = genesetdb.GREAT)
## $result
## SparrowResult (max FDR by collection set to 0.20%)
## --------------------------------------------------- 
##   collection method geneset_count sig_count sig_up sig_down
## 1          H  fgsea            50         1      0        1

See the Contrast_Comparisons vignette for more advanced use cases of gCrisprTools and extension to complex experiments and study designs.

Finally, you can make reports in a directory of interest:

path2report <-      #Make a report of the whole experiment
  ct.makeReport(fit = fit, 
                eset = es, 
                sampleKey = sk, 
                annotation = ann, 
                results = resultsDF, 
                aln = aln, 
                outdir = ".") 

path2QC <-          #Or one focusing only on experiment QC
  ct.makeQCReport(es, 
                  trim = 1000, 
                  log2.ratio = 0.05, 
                  sampleKey = sk, 
                  annotation = ann, 
                  aln = aln, 
                  identifier = 'Crispr_QC_Report',
                  lib.size = NULL
                  )                

path2Contrast <-    #Or Contrast-specific one
  ct.makeContrastReport(eset = es, 
                        fit = fit, 
                        sampleKey = sk, 
                        results = resultsDF, 
                        annotation = ann, 
                        comparison.id = NULL, 
                        identifier = 'Crispr_Contrast_Report')            
sessionInfo()
## R Under development (unstable) (2023-10-22 r85388)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] limma_3.59.0        gCrisprTools_2.9.0  Biobase_2.63.0     
## [4] BiocGenerics_0.49.0
## 
## loaded via a namespace (and not attached):
##   [1] DBI_1.1.3                   bitops_1.0-7               
##   [3] gridExtra_2.3               BiasedUrn_2.0.11           
##   [5] GSEABase_1.65.0             BiocSet_1.17.0             
##   [7] rlang_1.1.1                 magrittr_2.0.3             
##   [9] clue_0.3-65                 GetoptLong_1.0.5           
##  [11] msigdbr_7.5.1               sparrow_1.9.0              
##  [13] matrixStats_1.0.0           compiler_4.4.0             
##  [15] RSQLite_2.3.1               png_0.1-8                  
##  [17] vctrs_0.6.4                 pkgconfig_2.0.3            
##  [19] shape_1.4.6                 crayon_1.5.2               
##  [21] fastmap_1.1.1               backports_1.4.1            
##  [23] magick_2.8.1                XVector_0.43.0             
##  [25] labeling_0.4.3              utf8_1.2.4                 
##  [27] rmarkdown_2.25              markdown_1.11              
##  [29] graph_1.81.0                purrr_1.0.2                
##  [31] bit_4.0.5                   xfun_0.40                  
##  [33] zlibbioc_1.49.0             cachem_1.0.8               
##  [35] GenomeInfoDb_1.39.0         jsonlite_1.8.7             
##  [37] blob_1.2.4                  DelayedArray_0.29.0        
##  [39] BiocParallel_1.37.0         irlba_2.3.5.1              
##  [41] parallel_4.4.0              cluster_2.1.4              
##  [43] R6_2.5.1                    RColorBrewer_1.1-3         
##  [45] GenomicRanges_1.55.0        Rcpp_1.0.11                
##  [47] SummarizedExperiment_1.33.0 iterators_1.0.14           
##  [49] knitr_1.44                  IRanges_2.37.0             
##  [51] Matrix_1.6-1.1              tidyselect_1.2.0           
##  [53] yaml_2.3.7                  viridis_0.6.4              
##  [55] abind_1.4-5                 doParallel_1.0.17          
##  [57] codetools_0.2-19            plyr_1.8.9                 
##  [59] lattice_0.22-5              tibble_3.2.1               
##  [61] withr_2.5.1                 KEGGREST_1.43.0            
##  [63] evaluate_0.22               ontologyIndex_2.11         
##  [65] circlize_0.4.15             Biostrings_2.71.0          
##  [67] pillar_1.9.0                MatrixGenerics_1.15.0      
##  [69] checkmate_2.2.0             foreach_1.5.2              
##  [71] stats4_4.4.0                plotly_4.10.3              
##  [73] generics_0.1.3              RCurl_1.98-1.12            
##  [75] S4Vectors_0.41.0            ggplot2_3.4.4              
##  [77] commonmark_1.9.0            munsell_0.5.0              
##  [79] scales_1.2.1                xtable_1.8-4               
##  [81] glue_1.6.2                  lazyeval_0.2.2             
##  [83] tools_4.4.0                 BiocIO_1.13.0              
##  [85] data.table_1.14.8           fgsea_1.29.0               
##  [87] annotate_1.81.0             locfit_1.5-9.8             
##  [89] babelgene_22.9              XML_3.99-0.14              
##  [91] fastmatch_1.1-4             cowplot_1.1.1              
##  [93] grid_4.4.0                  Cairo_1.6-1                
##  [95] tidyr_1.3.0                 AnnotationDbi_1.65.0       
##  [97] edgeR_4.1.0                 colorspace_2.1-0           
##  [99] GenomeInfoDbData_1.2.11     cli_3.6.1                  
## [101] fansi_1.0.5                 viridisLite_0.4.2          
## [103] S4Arrays_1.3.0              ComplexHeatmap_2.19.0      
## [105] dplyr_1.1.3                 gtable_0.3.4               
## [107] digest_0.6.33               SparseArray_1.3.0          
## [109] farver_2.1.1                rjson_0.2.21               
## [111] htmlwidgets_1.6.2           memoise_2.0.1              
## [113] htmltools_0.5.6.1           lifecycle_1.0.3            
## [115] httr_1.4.7                  mime_0.12                  
## [117] GlobalOptions_0.1.2         statmod_1.5.0              
## [119] bit64_4.0.5