Contents

knitr::opts_chunk$set(tidy=FALSE, cache=TRUE,
                      dev="png", message=FALSE, error=FALSE, warning=TRUE)

Note: if you use MAGeCKFlute in published research, please cite: Binbin Wang, Mei Wang, Wubing Zhang. “Integrative analysis of pooled CRISPR genetic screens using MAGeCKFlute.” Nature Protocols (2019), doi: 10.1038/s41596-018-0113-7.

0.1 Install and load the required packages

if(!"MAGeCKFlute" %in% installed.packages()) BiocManager::install("MAGeCKFlute")
if(!"clusterProfiler" %in% installed.packages()) BiocManager::install("clusterProfiler")
if(!"ggplot2" %in% installed.packages()) BiocManager::install("ggplot2")

library(MAGeCKFlute)
library(clusterProfiler)
library(ggplot2)

1 Quick start

1.1 Downstream analysis of MAGeCK RRA

The MAGeCK (mageck test) uses Robust Rank Aggregation (RRA) for robust identification of CRISPR-screen hits, and outputs the summary results at both sgRNA and gene level. Before performing the downstream analysis, please make sure you have got the gene summary and sgRNA summary results from mageck test. MAGeCKFlute incorporates an example datasets (Deng Pan 2018), which study the gene functions in T cell mediate tumor killing by performing CRISPR screens in B16F10 cell lines co-cultured with Pmel-1 T cells. There are three samples in the data, including Pmel1_Input (B16F10 cells without T cell co-culture), Pmel1_Ctrl (B16F10 cells co-cultured with control T cells), and Pmel1 (B16F10 cells co-cultured with antigen specific T cells). We compared Pmel1 with Pmel1_Ctrl using mageck test, which identifies genes associated with T cell mediated tumor killing.

1.1.1 Run the FluteRRA pipeline

MAGeCKFlute processes MAGeCK RRA results (“gene_summary” and “sgrna_summary”) with the function FluteRRA, which identifies positively and negatively selected genes and performs functional analysis.

## path to the gene summary file (required)
file1 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/rra.gene_summary.txt")
## path to the sgRNA summary file (optional)
file2 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/rra.sgrna_summary.txt")
# Run FluteRRA with only gene summary file
FluteRRA(file1, proj="Pmel1", organism="mmu", outdir = "./")

# Run FluteRRA with both gene summary file and sgRNA summary file
FluteRRA(file1, file2, proj="Pmel1", organism="mmu", outdir = "./")

1.1.2 Important parameters

  • incorporateDepmap will allow users to compare the selected genes in the input data and Depmap data, which could clean the selected genes by eliminating lethal genes (false positive hits in some CRISPR screens) and then perform downstream enrichment analysis. You can set the incorporateDepmap to be TRUE if you want to include this analysis. Users could specify Depmap data of specific cell lines by setting the cell_lines and lineages.
FluteRRA(file1, proj="Pmel1", organism="mmu", incorporateDepmap = TRUE, outdir = "./")
  • omitEssential allows users to eliminate all common essential genes selected in the Depmap from all of the downstream analysis.
FluteRRA(gdata, proj="Pmel1", organism="mmu", omitEssential = TRUE, outdir = "./")

Hints: all figures and intermediate data are saved into local directory “./MAGeCKFlute_Test/”, and all figures are integrated into file “FluteRRA_Test.pdf”.

For more available parameters in FluteRRA, please read the documentation using the command ?FluteRRA.

1.2 Downstream analysis of MAGeCK MLE

The MAGeCK-VISPR (mageck mle) utilizes a maximum likelihood estimation (MLE) for robust identification of CRISPR screen hits. It outputs beta scores and the associated statistics in multiple conditions. The beta score describes how a gene is selected: a positive beta score indicates positive selection, and a negative beta score indicates negative selection. Before using FluteMLE, you should first get gene summary result from MAGeCK-VISPR (mageck mle). MAGeCKFlute uses the same dataset as before for demonstration. Using mageck mle, we removed the baseline effect (plasmid sample) from all the three samples, including Pmel1_Input (B16F10 cells without T cell co-culture), Pmel1_Ctrl (B16F10 cells co-cultured with control T cells), and Pmel1 (B16F10 cells co-cultured with antigen specific T cells).

1.2.1 Run the FluteMLE pipeline

MAGeCKFlute processes MAGeCK MLE results (“gene_summary”) with the function FluteMLE, which performs QC and data normalization based on the beta score, identifies essential genes by comparing control and treatment samples, and finally performs functional enrichment analysis.

file3 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/mle.gene_summary.txt")
FluteMLE(file3, treatname="Pmel1", ctrlname="Pmel1_Ctrl", proj="Pmel1", organism="mmu")

1.2.2 Important parameters

  • incorporateDepmap will allow users to compare the selected genes in the input data and Depmap data, which could clean the selected genes by eliminating lethal genes (false positive hits in some CRISPR screens) and then perform downstream enrichment analysis. You can set the incorporateDepmap to be TRUE if you want to include this analysis. Users could specify Depmap data of specific cell lines by setting the cell_lines and lineages.
## Take Depmap screen as control
FluteMLE(gdata, treatname="Pmel1_Ctrl", ctrlname="Depmap", proj="Pmel1", organism="mmu", incorporateDepmap = TRUE)
  • omitEssential allows users to eliminate all common essential genes selected in the Depmap from all of the downstream analysis.
FluteMLE(gdata, treatname="Pmel1", ctrlname="Pmel1_Ctrl", proj="Pmel1", organism="mmu", omitEssential = TRUE)

Hint: All pipeline results are written into local directory “./MAGeCKFlute_Test/”, and all figures are integrated into file “FluteMLE_Test.pdf”.

For more available parameters in FluteMLE, please read the documentation using the command ?FluteMLE.

2 Step by step analysis

2.1 Section I: Quality control

2.1.1 Input data

MAGeCK/MAGeCK-VISPR outputs a count summary file, which summarizes some basic QC scores at raw count level, including map ratio, Gini index, and NegSelQC. MAGeCKFlute incorporates an example datasets (Ophir Shalem* 2014) for demonstration.

file4 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/countsummary.txt")
countsummary = read.delim(file4, check.names = FALSE)
head(countsummary)
##                                   File    Label    Reads   Mapped Percentage
## 1 ../data/GSC_0131_Day23_Rep1.fastq.gz day23_r1 62818064 39992777     0.6366
## 2  ../data/GSC_0131_Day0_Rep2.fastq.gz  day0_r2 47289074 31709075     0.6705
## 3  ../data/GSC_0131_Day0_Rep1.fastq.gz  day0_r1 51190401 34729858     0.6784
## 4 ../data/GSC_0131_Day23_Rep2.fastq.gz day23_r2 58686580 37836392     0.6447
##   TotalsgRNAs Zerocounts GiniIndex NegSelQC NegSelQCPval
## 1       64076         57   0.08510        0            1
## 2       64076         17   0.07496        0            1
## 3       64076         14   0.07335        0            1
## 4       64076         51   0.08587        0            1
##   NegSelQCPvalPermutation NegSelQCPvalPermutationFDR NegSelQCGene
## 1                       1                          1            0
## 2                       1                          1            0
## 3                       1                          1            0
## 4                       1                          1            0

2.1.2 Visualize the QC results

# Gini index
BarView(countsummary, x = "Label", y = "GiniIndex",
        ylab = "Gini index", main = "Evenness of sgRNA reads")

# Missed sgRNAs
countsummary$Missed = log10(countsummary$Zerocounts)
BarView(countsummary, x = "Label", y = "Missed", fill = "#394E80",
        ylab = "Log10 missed gRNAs", main = "Missed sgRNAs")

# Read mapping
MapRatesView(countsummary)

# Or
countsummary$Unmapped = countsummary$Reads - countsummary$Mapped
gg = reshape2::melt(countsummary[, c("Label", "Mapped", "Unmapped")], id.vars = "Label")
gg$variable = factor(gg$variable, levels = c("Unmapped", "Mapped"))
gg = gg[order(gg$Label, gg$variable), ]
p = BarView(gg, x = "Label", y = "value", fill = "variable", 
            position = "stack", xlab = NULL, ylab = "Reads", main = "Map ratio")
p + scale_fill_manual(values = c("#9BC7E9", "#1C6DAB"))

2.2 Section II: Downstream analysis of MAGeCK RRA

For CRISPR/Cas9 screens with two experimental conditions, MAGeCK-RRA is available for identification of essential genes. In MAGeCK-RRA results, the sgRNA summary and gene summary file summarizes the statistical significance of positive selections and negative selections at sgRNA level and gene level.

2.2.1 Read the required data

file1 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/rra.gene_summary.txt")
gdata = ReadRRA(file1)
head(gdata)
##       id   Score      FDR
## 1  Cd274 -3.4698 0.002475
## 2  Psmb8 -3.7260 0.002475
## 3   Rela -3.2413 0.008251
## 4 Otulin -3.6018 0.008663
## 5  Ikbkb -3.4256 0.010891
## 6  Tcea1 -3.4236 0.042079
file2 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/rra.sgrna_summary.txt")
sdata = ReadsgRRA(file2)
head(sdata)
##                  sgrna  Gene     LFC FDR
## 1 AAACATATAGTGTACCTCTA  Jak1 10.8910   0
## 2 TCCGAACCGAATCATCACTG  Jak1 10.9170   0
## 3 TGAATAAATCCATCAGACAG  Jak1 10.5970   0
## 4 GGATAGACGCCCAGCCACTG Stat1  9.9921   0
## 5 TTAATGACGAGCTCGTGGAG Stat1  9.2728   0
## 6 GAAAAGCAAGCGTAATCTCC Stat1  8.7931   0

2.3 To incorporate depmap data that are profiled in human cell lines, we will convert mouse gene names to homologous human genes for this dataset.

gdata$HumanGene = TransGeneID(gdata$id, fromType = "symbol", toType = "symbol",
                              fromOrg = "mmu", toOrg = "hsa")
sdata$HumanGene = TransGeneID(sdata$Gene, fromType = "symbol", toType = "symbol",
                              fromOrg = "mmu", toOrg = "hsa")

2.3.1 Compute the similarity between the CRISPR screen with Depmap screens

## Remove missing or duplicate human genes
idx = duplicated(gdata$HumanGene)|is.na(gdata$HumanGene)
gdata = gdata[!idx, ]
depmap_similarity = ResembleDepmap(gdata, symbol = "HumanGene", score = "Score")
head(depmap_similarity)

2.3.2 Omit common essential genes from the data

gdata = OmitCommonEssential(gdata, symbol = "HumanGene")
sdata = OmitCommonEssential(sdata, symbol = "HumanGene")

2.3.3 Visualization of negative selections and positive selections

2.3.3.1 Volcano plot

gdata$LogFDR = -log10(gdata$FDR)
p1 = ScatterView(gdata, x = "Score", y = "LogFDR", label = "id", 
                 model = "volcano", top = 5)
print(p1)

# Or
p2 = VolcanoView(gdata, x = "Score", y = "FDR", Label = "id")
print(p2)

2.3.3.2 Rank plot

Rank all the genes based on their scores and label genes in the rank plot.

gdata$Rank = rank(gdata$Score)
p1 = ScatterView(gdata, x = "Rank", y = "Score", label = "id", 
                 top = 5, auto_cut_y = TRUE, ylab = "Log2FC", 
                 groups = c("top", "bottom"))
print(p1)
## Warning: ggrepel: 5 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

Label interested hits using parameter toplabels (in ScatterView) and genelist (in RankView).

ScatterView(gdata, x = "Rank", y = "Score", label = "id", top = 5, 
            auto_cut_y = TRUE, groups = c("top", "bottom"), 
            ylab = "Log2FC", toplabels = c("Pbrm1", "Arid2", "Brd7"))
## Warning: ggrepel: 8 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

You can also use the function RankView to draw the figure.

geneList= gdata$Score
names(geneList) = gdata$id
p2 = RankView(geneList, top = 5, bottom = 10)
print(p2)

RankView(geneList, top = 0, bottom = 0, genelist = c("Pbrm1", "Arid2", "Brd7"))

2.3.3.3 Dot plot

Visualize negative and positive selected genes separately.

gdata$RandomIndex = sample(1:nrow(gdata), nrow(gdata))
gdata = gdata[order(-gdata$Score), ]
gg = gdata[gdata$Score>0, ]
p1 = ScatterView(gg, x = "RandomIndex", y = "Score", label = "id",
                 y_cut = CutoffCalling(gdata$Score,2), 
                 groups = "top", top = 5, ylab = "Log2FC")
p1

gg = gdata[gdata$Score<0, ]
p2 = ScatterView(gg, x = "RandomIndex", y = "Score", label = "id",
                 y_cut = CutoffCalling(gdata$Score,2), 
                 groups = "bottom", top = 5, ylab = "Log2FC")
p2

2.3.3.4 sgRankView - visualize the rank of sgRNAs targeting top selected genes.

p2 = sgRankView(sdata, top = 4, bottom = 4)
print(p2)

2.3.4 Enrichment analysis

For more information about functional enrichment analysis in MAGeCKFlute, please read the MAGeCKFlute_enrichment document, in which we introduce all the available options and methods.

geneList= gdata$Score
names(geneList) = gdata$id
enrich_pos = EnrichAnalyzer(geneList = geneList[geneList>0.5], 
                            method = "HGT", type = "KEGG")
enrich_neg = EnrichAnalyzer(geneList = geneList[geneList< -0.5], 
                            method = "HGT", type = "KEGG")

2.3.4.1 Visualization of enrichment results

EnrichedView(enrich_pos, mode = 1, top = 5, bottom = 0)

EnrichedView(enrich_pos, mode = 2, top = 5, bottom = 0)

EnrichedView(enrich_neg, mode = 2, top = 0, bottom = 5)

2.4 Section III: Downstream analysis of MAGeCK MLE

The MAGeCK-VISPR (mageck mle) computes beta scores and the associated statistics for all genes in multiple conditions. The beta score describes how the gene is selected: a positive beta score indicates a positive selection, and a negative beta score indicates a negative selection. Before using FluteMLE, you should first get gene summary result from MAGeCK-VISPR (mageck mle). MAGeCKFlute incorporates an example datasets for demonstration.

2.4.2 read gene summary file

file3 = file.path(system.file("extdata", package = "MAGeCKFlute"),
                  "testdata/mle.gene_summary.txt")
# Read and visualize the file format
gdata = ReadBeta(file3)
head(gdata)
##      Gene Pmel1_Input Pmel1_Ctrl     Pmel1
## 1  Defb34    0.087884  -0.010034 -0.068015
## 2   Mndal    0.282860   0.033850 -0.161790
## 3   Cox8c    0.212140  -0.031618 -0.491360
## 4 Poldip3    0.125260  -0.123200 -0.450270
## 5   Bcas2   -0.196670   0.031254 -0.457530
## 6  Klk1b9   -0.220270  -0.019697 -0.443300

2.4.3 Normalization of beta scores

It is difficult to control all samples with a consistent cell cycle in a CRISPR screen experiment with multi conditions. Besides, beta score among different states with an inconsistent cell cycle is incomparable. So it is necessary to do the normalization when comparing the beta scores in different conditions. Essential genes are those genes that are indispensable for its survival. The effect generated by knocking out these genes in different cell types is consistent. Based on this, we developed the cell cycle normalization method to shorten the gap of the cell cycle in different conditions.

ctrlname = "Pmel1_Ctrl"
treatname = "Pmel1"
gdata$HumanGene = TransGeneID(gdata$Gene, fromType = "symbol", toType = "symbol",
                              fromOrg = "mmu", toOrg = "hsa")

gdata_cc = NormalizeBeta(gdata, id = "HumanGene", samples=c(ctrlname, treatname), 
                         method="cell_cycle")
head(gdata_cc)
##      Gene Pmel1_Input  Pmel1_Ctrl      Pmel1 HumanGene
## 1  Defb34    0.087884 -0.06263323 -0.2425685  DEFB106A
## 2   Mndal    0.282860  0.21129508 -0.5770074     IFI16
## 3   Cox8c    0.212140 -0.19736271 -1.7523850     COX8C
## 4 Poldip3    0.125260 -0.76902670 -1.6058418   POLDIP3
## 5   Bcas2   -0.196670  0.19509059 -1.6317338     BCAS2
## 6  Klk1b9   -0.220270 -0.12295064 -1.5809840      KLK3

2.4.4 Distribution of all gene beta scores

After normalization, the distribution of beta scores in different conditions should be similar. We can evaluate the distribution of beta scores using the function ‘DensityView’, and ‘ConsistencyView’.

DensityView(gdata_cc, samples=c(ctrlname, treatname))

ConsistencyView(gdata_cc, ctrlname, treatname)

# Another option MAView
MAView(gdata_cc, ctrlname, treatname)

2.4.5 Positive selection and negative selection

gdata_cc$Control = rowMeans(gdata_cc[,ctrlname, drop = FALSE])
gdata_cc$Treatment = rowMeans(gdata_cc[,treatname, drop = FALSE])

p1 = ScatterView(gdata_cc, "Control", "Treatment", label = "Gene", 
                 auto_cut_diag = TRUE, display_cut = TRUE, 
                 groups = c("top", "bottom"),
                 toplabels = c("Pbrm1", "Brd7", "Arid2", "Jak1", "Stat1", "B2m"))
print(p1)
## Warning: ggrepel: 3 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

2.4.5.1 Rank plot - label top hits

gdata_cc$Diff = gdata_cc$Treatment - gdata_cc$Control
gdata_cc$Rank = rank(gdata_cc$Diff)
p1 = ScatterView(gdata_cc, x = "Rank", y =  "Diff", label = "Gene", 
                 top = 5, auto_cut_y = TRUE, groups = c("top", "bottom"))
print(p1)
## Warning: ggrepel: 5 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

# Or
rankdata = gdata_cc$Treatment - gdata_cc$Control
names(rankdata) = gdata_cc$Gene
RankView(rankdata)

2.4.5.2 Nine-square scatter plot - identify treatment-associated genes

p1 = ScatterView(gdata_cc, x="Pmel1_Ctrl", y="Pmel1", label = "Gene", 
                 model = "ninesquare", top = 5, display_cut = TRUE, force = 2)
print(p1)
## Warning: ggrepel: 10 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

Customize the cutoff

p1 = ScatterView(gdata_cc, x="Pmel1_Ctrl", y="Pmel1", label = "Gene", 
                 model = "ninesquare", top = 5, display_cut = TRUE, 
                 x_cut = c(-1,1), y_cut = c(-1,1))
print(p1)
## Warning: ggrepel: 5 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps

Or

p2 = SquareView(gdata_cc, label = "Gene", 
                x_cutoff = CutoffCalling(gdata_cc$Control, 2), 
                y_cutoff = CutoffCalling(gdata_cc$Treatment, 2))
print(p2)

2.4.6 Functional analysis for treatment-associated genes

# 9-square groups
Square9 = p1$data
idx=Square9$group=="topcenter"
geneList = Square9$Diff[idx]
names(geneList) = Square9$Gene[idx]
universe = Square9$Gene

# Enrichment analysis
kegg1 = EnrichAnalyzer(geneList = geneList, universe = universe)
EnrichedView(kegg1, top = 6, bottom = 0)

Also, pathway visualization can be done using function KeggPathwayView (Luo et al. 2013).

gdata_cc$Entrez = TransGeneID(gdata_cc$Gene, "symbol", "entrez", organism = "mmu")
genedata = gdata_cc[, c("Entrez", "Control","Treatment")]
arrangePathview(genedata, pathways = "mmu04630", organism = "mmu", sub = NULL)
## Warning in data(bods): data set 'bods' not found

3 Session info

sessionInfo()
## R version 4.3.1 (2023-06-16)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.18-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_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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] ggplot2_3.4.4          clusterProfiler_4.10.0 MAGeCKFlute_2.6.0     
## [4] BiocStyle_2.30.0      
## 
## loaded via a namespace (and not attached):
##   [1] splines_4.3.1                 later_1.3.1                  
##   [3] bitops_1.0-7                  ggplotify_0.1.2              
##   [5] filelock_1.0.2                tibble_3.2.1                 
##   [7] polyclip_1.10-6               graph_1.80.0                 
##   [9] XML_3.99-0.14                 lifecycle_1.0.3              
##  [11] edgeR_4.0.0                   lattice_0.22-5               
##  [13] MASS_7.3-60                   magrittr_2.0.3               
##  [15] limma_3.58.0                  sass_0.4.7                   
##  [17] rmarkdown_2.25                jquerylib_0.1.4              
##  [19] yaml_2.3.7                    httpuv_1.6.12                
##  [21] depmap_1.15.0                 cowplot_1.1.1                
##  [23] DBI_1.1.3                     RColorBrewer_1.1-3           
##  [25] zlibbioc_1.48.0               purrr_1.0.2                  
##  [27] ggraph_2.1.0                  BiocGenerics_0.48.0          
##  [29] msigdbr_7.5.1                 RCurl_1.98-1.12              
##  [31] yulab.utils_0.1.0             tweenr_2.0.2                 
##  [33] rappdirs_0.3.3                sva_3.50.0                   
##  [35] GenomeInfoDbData_1.2.11       IRanges_2.36.0               
##  [37] S4Vectors_0.40.0              enrichplot_1.22.0            
##  [39] ggrepel_0.9.4                 tidytree_0.4.5               
##  [41] genefilter_1.84.0             pheatmap_1.0.12              
##  [43] annotate_1.80.0               codetools_0.2-19             
##  [45] DOSE_3.28.0                   ggforce_0.4.1                
##  [47] tidyselect_1.2.0              aplot_0.2.2                  
##  [49] farver_2.1.1                  viridis_0.6.4                
##  [51] matrixStats_1.0.0             stats4_4.3.1                 
##  [53] BiocFileCache_2.10.0          pathview_1.42.0              
##  [55] jsonlite_1.8.7                ellipsis_0.3.2               
##  [57] tidygraph_1.2.3               survival_3.5-7               
##  [59] tools_4.3.1                   treeio_1.26.0                
##  [61] HPO.db_0.99.2                 Rcpp_1.0.11                  
##  [63] glue_1.6.2                    gridExtra_2.3                
##  [65] xfun_0.40                     mgcv_1.9-0                   
##  [67] qvalue_2.34.0                 MatrixGenerics_1.14.0        
##  [69] GenomeInfoDb_1.38.0           dplyr_1.1.3                  
##  [71] withr_2.5.1                   BiocManager_1.30.22          
##  [73] fastmap_1.1.1                 fansi_1.0.5                  
##  [75] digest_0.6.33                 R6_2.5.1                     
##  [77] mime_0.12                     gridGraphics_0.5-1           
##  [79] colorspace_2.1-0              GO.db_3.18.0                 
##  [81] RSQLite_2.3.1                 utf8_1.2.4                   
##  [83] tidyr_1.3.0                   generics_0.1.3               
##  [85] data.table_1.14.8             graphlayouts_1.0.1           
##  [87] httr_1.4.7                    scatterpie_0.2.1             
##  [89] pkgconfig_2.0.3               gtable_0.3.4                 
##  [91] blob_1.2.4                    XVector_0.42.0               
##  [93] shadowtext_0.1.2              htmltools_0.5.6.1            
##  [95] bookdown_0.36                 fgsea_1.28.0                 
##  [97] scales_1.2.1                  Biobase_2.62.0               
##  [99] png_0.1-8                     ggfun_0.1.3                  
## [101] knitr_1.44                    reshape2_1.4.4               
## [103] nlme_3.1-163                  curl_5.1.0                   
## [105] org.Hs.eg.db_3.18.0           cachem_1.0.8                 
## [107] stringr_1.5.0                 BiocVersion_3.18.0           
## [109] parallel_4.3.1                HDO.db_0.99.1                
## [111] AnnotationDbi_1.64.0          pillar_1.9.0                 
## [113] grid_4.3.1                    vctrs_0.6.4                  
## [115] promises_1.2.1                dbplyr_2.3.4                 
## [117] xtable_1.8-4                  Rgraphviz_2.46.0             
## [119] evaluate_0.22                 KEGGgraph_1.62.0             
## [121] magick_2.8.1                  cli_3.6.1                    
## [123] locfit_1.5-9.8                compiler_4.3.1               
## [125] rlang_1.1.1                   crayon_1.5.2                 
## [127] labeling_0.4.3                plyr_1.8.9                   
## [129] fs_1.6.3                      stringi_1.7.12               
## [131] viridisLite_0.4.2             BiocParallel_1.36.0          
## [133] babelgene_22.9                MPO.db_0.99.7                
## [135] munsell_0.5.0                 Biostrings_2.70.0            
## [137] lazyeval_0.2.2                GOSemSim_2.28.0              
## [139] Matrix_1.6-1.1                ExperimentHub_2.10.0         
## [141] patchwork_1.1.3               bit64_4.0.5                  
## [143] KEGGREST_1.42.0               statmod_1.5.0                
## [145] shiny_1.7.5.1                 interactiveDisplayBase_1.40.0
## [147] AnnotationHub_3.10.0          igraph_1.5.1                 
## [149] memoise_2.0.1                 bslib_0.5.1                  
## [151] ggtree_3.10.0                 fastmatch_1.1-4              
## [153] bit_4.0.5                     ape_5.7-1                    
## [155] gson_0.1.0

References

Deng Pan, Peng Jiang, Aya Kobayashi. 2018. “A major chromatin regulator determines resistance of tumor cells to T cell-mediated killing.” https://www.science.org/lookup/doi/10.1126/science.aao1710.

Luo, Weijun, Brouwer, and Cory. 2013. “Pathview: An R/Bioconductor Package for Pathway-Based Data Integration and Visualization.” Bioinformatics 29 (14): 1830–1. https://doi.org/10.1093/bioinformatics/btt285.

Ophir Shalem*, Ella Hartenian, Neville E. Sanjana*. 2014. “Genome-scale CRISPR-Cas9 knockout screening in human cells.” http://science.sciencemag.org/content/343/6166/84.long.

Wei Li, Han Xu, Johannes Köster, and X. Shirley Liu. 2015. “Quality control, modeling, and visualization of CRISPR screens with MAGeCK-VISPR.” https://genomebiology.biomedcentral.com/articles/10.1186/s13059-015-0843-6.

Wei Li, Tengfei Xiao, Han Xu, and X Shirley Liu. 2014. “MAGeCK enables robust identification of essential genes from genome-scale CRISPR/Cas9 knockout screens.” https://genomebiology.biomedcentral.com/articles/10.1186/s13059-014-0554-4.