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.
library(MAGeCKFlute)
library(ggplot2)
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 (Ophir Shalem1 2014) for demonstration, shown as below.
file1 = file.path(system.file("extdata", package = "MAGeCKFlute"),
"testdata/rra.gene_summary.txt")
# Read and visualize the file format
gdata = read.delim(file1, check.names = FALSE)
head(gdata)
## id num neg|score neg|p-value neg|fdr neg|rank neg|goodsgrna neg|lfc
## 1 CREBBP 4 1.00000 1.00000 1 8269 0 0.96608
## 2 EP300 4 1.00000 1.00000 1 8270 0 1.02780
## 3 CHD 4 0.99999 0.99999 1 8268 0 0.59265
## 4 C16orf72 4 0.99998 0.99998 1 8267 0 0.82307
## 5 CACNB2 5 0.99863 0.99864 1 8243 0 0.39268
## 6 FGF12 6 0.54609 0.71111 1 4862 1 0.40822
## pos|score pos|p-value pos|fdr pos|rank pos|goodsgrna pos|lfc
## 1 3.17e-12 2.9700e-07 0.001650 1 4 0.96608
## 2 3.32e-10 2.9700e-07 0.001650 2 4 1.02780
## 3 1.22e-05 4.4300e-05 0.092203 3 4 0.59265
## 4 2.06e-05 8.0000e-05 0.147965 4 4 0.82307
## 5 4.45e-05 2.1084e-04 0.319082 5 5 0.39268
## 6 5.46e-05 2.7627e-04 0.336280 6 5 0.40822
You can also read the file using ReadRRA
in MAGeCKFlute
gdata = ReadRRA(file1)
head(gdata)
## id Score FDR
## 1 CREBBP 0.96608 0.001650
## 2 EP300 1.02780 0.001650
## 3 CHD 0.59265 0.092203
## 4 C16orf72 0.82307 0.147965
## 5 CACNB2 0.39268 0.319082
## 6 FGF12 0.40822 0.336280
Hints: you can also use a data from other analysis, just make sure the three columns (id, Score, and FDR) are avaible in the data.
file2 = file.path(system.file("extdata", package = "MAGeCKFlute"),
"testdata/rra.sgrna_summary.txt")
sdata = read.delim(file2)
head(sdata)
## sgrna Gene control_count treatment_count control_mean treat_mean LFC
## 1 s_36798 NF2 8917/21204 5020.7/5127.9 15061.0 5074.30 -1.5693
## 2 s_45763 RAB6A 3375.8/3667.7 372.88/357.79 3521.8 365.33 -3.2655
## 3 s_50164 SF1 3657.8/3352.6 453.62/628.28 3505.2 540.95 -2.6937
## 4 s_20780 FDXR 3444.8/3191.1 552.51/555.38 3317.9 553.95 -2.5803
## 5 s_36796 NF2 5492.9/11396 3832.2/3794.6 8444.6 3813.40 -1.1467
## 6 s_36799 NF2 1805.2/5541.7 695.86/882.47 3673.5 789.16 -2.2173
## control_var adj_var score p.low p.high p.twosided FDR
## 1 75491000 78871 35.559 2.9804e-277 1 5.9609e-277 1.2732e-272
## 2 42617 15519 25.338 6.1638e-142 1 1.2328e-141 1.9748e-137
## 3 46575 15438 23.857 4.2365e-126 1 8.4731e-126 9.0487e-122
## 4 32179 14520 22.938 9.7946e-117 1 1.9589e-116 1.5690e-112
## 5 17425000 41247 22.803 2.1321e-115 1 4.2641e-115 3.0359e-111
## 6 6980700 16267 22.615 1.5592e-113 1 3.1183e-113 1.9981e-109
## high_in_treatment
## 1 False
## 2 False
## 3 False
## 4 False
## 5 False
## 6 False
You can also read the file using ReadsgRRA
in MAGeCKFlute
sdata = ReadsgRRA(file2)
head(sdata)
## sgrna Gene LFC FDR
## 1 s_36798 NF2 -1.5693 1.2732e-272
## 2 s_45763 RAB6A -3.2655 1.9748e-137
## 3 s_50164 SF1 -2.6937 9.0487e-122
## 4 s_20780 FDXR -2.5803 1.5690e-112
## 5 s_36796 NF2 -1.1467 3.0359e-111
## 6 s_36799 NF2 -2.2173 1.9981e-109
Run the downstream analysis pipeline with both gene summary and sgrna summary
FluteRRA(file1, file2, proj="Test", organism="hsa", scale_cutoff = 1, outdir = "./")
# Or
FluteRRA(gdata, sdata, proj="Test", organism="hsa", scale_cutoff = 1, outdir = "./")
Run the downstream analysis pipeline with only gene summary file
FluteRRA(file1, proj="Test", organism="hsa", outdir = "./")
# Or
FluteRRA(gdata, proj="Test", organism="hsa", outdir = "./")
Incorporate Depmap data into analysis
FluteRRA(gdata, proj="Test", organism="hsa", incorporateDepmap = TRUE,
outdir = "./")
Omit common essential genes in the analysis
FluteRRA(gdata, proj="Test", organism="hsa", incorporateDepmap = TRUE,
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
?FluteRRA
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.
file3 = file.path(system.file("extdata", package = "MAGeCKFlute"),
"testdata/mle.gene_summary.txt")
# Read and visualize the file format
gdata = read.delim(file3, check.names = FALSE)
head(gdata)
## Gene sgRNA dmso|beta dmso|z dmso|p-value dmso|fdr dmso|wald-p-value
## 1 FEZ1 6 -0.045088 -0.66798 0.79649 0.97939 5.0415e-01
## 2 TNN 6 0.094325 1.36120 0.34176 0.89452 1.7344e-01
## 3 OAS2 8 -0.271210 -4.76860 0.46995 0.93572 1.8555e-06
## 4 JOSD1 4 0.092888 1.08610 0.34538 0.89543 2.7745e-01
## 5 CFH 6 -0.113230 -1.71070 0.95069 0.99513 8.7140e-02
## 6 C11orf68 4 -0.163690 -2.26570 0.77504 0.97694 2.3472e-02
## dmso|wald-fdr plx|beta plx|z plx|p-value plx|fdr plx|wald-p-value
## 1 6.3060e-01 -0.036721 -0.54346 0.81604 0.98345 5.8681e-01
## 2 2.8578e-01 0.065533 0.94344 0.47309 0.93207 3.4546e-01
## 3 1.4126e-05 -0.289010 -5.07170 0.40411 0.90933 3.9431e-07
## 4 4.0754e-01 0.094913 1.10840 0.39207 0.90504 2.6769e-01
## 5 1.6606e-01 -0.060018 -0.90751 0.90090 0.98876 3.6414e-01
## 6 5.7351e-02 -0.094403 -1.30760 0.97492 0.99833 1.9103e-01
## plx|wald-fdr
## 1 6.9940e-01
## 2 4.7400e-01
## 3 3.5296e-06
## 4 3.9013e-01
## 5 4.9438e-01
## 6 2.9934e-01
You can also read beta scores from the data using ReadBeta
in MAGeCKFlute
gdata = ReadBeta(file3)
head(gdata)
## Gene dmso plx
## 1 FEZ1 -0.045088 -0.036721
## 2 TNN 0.094325 0.065533
## 3 OAS2 -0.271210 -0.289010
## 4 JOSD1 0.092888 0.094913
## 5 CFH -0.113230 -0.060018
## 6 C11orf68 -0.163690 -0.094403
Hints: you can also run FluteMLE using other data, in which the first column is “Gene”, and other columns represent samples.
FluteMLE(file3, treatname="plx", ctrlname="dmso", proj="Test", organism="hsa")
# Or
FluteMLE(gdata, treatname="plx", ctrlname="dmso", proj="Test", organism="hsa")
If your data only include one condition, you can take Depmap screens as control.
## Take Depmap screen as control
FluteMLE(gdata, treatname="plx", ctrlname="Depmap", proj="PLX", organism="hsa", incorporateDepmap = TRUE)
If you are not interested in common essential genes, you can omit them in the analysis by setting a parameter “omitEssential”
FluteMLE(gdata, treatname="plx", ctrlname="Depmap", proj="PLX", organism="hsa", incorporateDepmap = TRUE, 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
?FluteMLE
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 Shalem1 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
# 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"))
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"))
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.
file1 = file.path(system.file("extdata", package = "MAGeCKFlute"),
"testdata/rra.gene_summary.txt")
gdata = ReadRRA(file1)
head(gdata)
## id Score FDR
## 1 CREBBP 0.96608 0.001650
## 2 EP300 1.02780 0.001650
## 3 CHD 0.59265 0.092203
## 4 C16orf72 0.82307 0.147965
## 5 CACNB2 0.39268 0.319082
## 6 FGF12 0.40822 0.336280
file2 = file.path(system.file("extdata", package = "MAGeCKFlute"),
"testdata/rra.sgrna_summary.txt")
sdata = ReadsgRRA(file2)
head(sdata)
## sgrna Gene LFC FDR
## 1 s_36798 NF2 -1.5693 1.2732e-272
## 2 s_45763 RAB6A -3.2655 1.9748e-137
## 3 s_50164 SF1 -2.6937 9.0487e-122
## 4 s_20780 FDXR -2.5803 1.5690e-112
## 5 s_36796 NF2 -1.1467 3.0359e-111
## 6 s_36799 NF2 -2.2173 1.9981e-109
## The first run must be time-consuming for downloading Depmap data automatically.
depmap_similarity = ResembleDepmap(gdata, symbol = "id", score = "Score")
gdata = OmitCommonEssential(gdata)
sdata = OmitCommonEssential(sdata, symbol = "Gene")
# Compute the similarity with Depmap screens based on subset genes
depmap_similarity = ResembleDepmap(gdata, symbol = "id", score = "Score")
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)
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)
Label interested hits using parameter toplabels
(in ScatterView) and genelist
(in RankView).
ScatterView(gdata, x = "Rank", y = "Score", label = "id",
auto_cut_y = TRUE, groups = c("top", "bottom"),
ylab = "Log2FC", toplabels = c("EP300", "NF2"))
## Warning: ggrepel: 1 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
Plot Log2FC at x-axis
ScatterView(gdata, x = "Score", y = "Rank", label = "id",
auto_cut_x = TRUE, groups = c("left", "right"),
xlab = "Log2FC", top = 3)
Or
geneList= gdata$Score
names(geneList) = gdata$id
p2 = RankView(geneList, top = 5, bottom = 10) + xlab("Log2FC")
print(p2)
## Warning: ggrepel: 9 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
RankView(geneList, top = 0, bottom = 0, genelist = c("EP300", "NF2")) + xlab("Log2FC")
Only plot positive selection
gdata$Rank = rank(-gdata$Score)
ScatterView(gdata[gdata$Score>0,], x = "Rank", y = "Score", label = "id",
auto_cut_y = TRUE, groups = c("top", "bottom"),
ylab = "Log2FC", top = 5)
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
p2 = sgRankView(sdata, top = 4, bottom = 4)
print(p2)
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 = EnrichAnalyzer(geneList = geneList[geneList>0.5],
method = "HGT", type = "KEGG")
EnrichedView(enrich, mode = 1, top = 5)
EnrichedView(enrich, mode = 2, top = 5)
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.
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 dmso plx
## 1 FEZ1 -0.045088 -0.036721
## 2 TNN 0.094325 0.065533
## 3 OAS2 -0.271210 -0.289010
## 4 JOSD1 0.092888 0.094913
## 5 CFH -0.113230 -0.060018
## 6 C11orf68 -0.163690 -0.094403
Is there batch effects? This is a commonly asked question before perform later analysis. In our package, we provide HeatmapView
to ensure whether the batch effect exists in data and use BatchRemove
to remove easily if same batch samples cluster together.
##Before batch removal
edata = matrix(c(rnorm(2000, 5), rnorm(2000, 8)), 1000)
colnames(edata) = paste0("s", 1:4)
HeatmapView(cor(edata))
## After batch removal
batchMat = data.frame(sample = colnames(edata), batch = rep(1:2, each = 2))
edata1 = BatchRemove(edata, batchMat)
head(edata1$data)
## s1 s2 s3 s4
## [1,] 7.191277 6.793563 7.272144 6.712056
## [2,] 7.006768 6.340287 7.453082 5.892356
## [3,] 6.051565 6.539662 6.568575 6.022173
## [4,] 7.324603 7.482559 8.112867 6.692783
## [5,] 7.879187 6.444309 7.369581 6.951306
## [6,] 5.883239 6.617949 7.065960 5.432753
print(edata1$p)
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 = "dmso"
treatname = "plx"
gdata_cc = NormalizeBeta(gdata, samples=c(ctrlname, treatname), method="cell_cycle")
head(gdata_cc)
## Gene dmso plx
## 1 FEZ1 -0.05913258 -0.05702550
## 2 TNN 0.12370654 0.10176880
## 3 OAS2 -0.35568991 -0.44881511
## 4 JOSD1 0.12182193 0.14739417
## 5 CFH -0.14850031 -0.09320434
## 6 C11orf68 -0.21467823 -0.14660217
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)
gdata_cc$Control = rowMeans(gdata_cc[,ctrlname, drop = FALSE])
gdata_cc$Treatment = rowMeans(gdata_cc[,treatname, drop = FALSE])
p1 = ScatterView(gdata_cc, "Control", "Treatment", groups = c("top", "bottom"), auto_cut_diag = TRUE, display_cut = TRUE, toplabels = c("NF1", "NF2", "EP300"))
print(p1)
gdata_cc$Diff = gdata_cc$Treatment - gdata_cc$Control
gdata_cc$Rank = rank(gdata_cc$Diff)
p1 = ScatterView(gdata_cc, x = "Diff", y = "Rank", label = "Gene",
top = 5, model = "rank")
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)
## Warning: ggrepel: 10 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
p1 = ScatterView(gdata_cc, x = "dmso", y = "plx", label = "Gene",
model = "ninesquare", top = 5, display_cut = TRUE, force = 2)
print(p1)
## Warning: ggrepel: 9 unlabeled data points (too many overlaps). Consider
## increasing max.overlaps
Customize the cutoff
p1 = ScatterView(gdata_cc, x = "dmso", y = "plx", label = "Gene",
model = "ninesquare", top = 5, display_cut = TRUE,
x_cut = c(-1,1), y_cut = c(-1,1))
print(p1)
## Warning: ggrepel: 7 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)
# 9-square groups
Square9 = p1$data
idx=Square9$group=="topcenter"
geneList = Square9$Diff
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).
genedata = gdata_cc[, c("Control","Treatment")]
arrangePathview(genedata, pathways = "hsa01521", organism = "hsa", sub = NULL)
sessionInfo()
## R version 4.1.0 (2021-05-18)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.13-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] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] ggplot2_3.3.3 MAGeCKFlute_1.12.0 BiocStyle_2.20.0
##
## loaded via a namespace (and not attached):
## [1] fgsea_1.18.0 colorspace_2.0-1 ggtree_3.0.0
## [4] ellipsis_0.3.2 qvalue_2.24.0 XVector_0.32.0
## [7] aplot_0.0.6 farver_2.1.0 graphlayouts_0.7.1
## [10] ggrepel_0.9.1 bit64_4.0.5 AnnotationDbi_1.54.0
## [13] fansi_0.4.2 scatterpie_0.1.6 codetools_0.2-18
## [16] splines_4.1.0 cachem_1.0.5 GOSemSim_2.18.0
## [19] knitr_1.33 polyclip_1.10-0 jsonlite_1.7.2
## [22] annotate_1.70.0 GO.db_3.13.0 png_0.1-7
## [25] pheatmap_1.0.12 ggforce_0.3.3 msigdbr_7.4.1
## [28] BiocManager_1.30.15 compiler_4.1.0 httr_1.4.2
## [31] rvcheck_0.1.8 assertthat_0.2.1 Matrix_1.3-3
## [34] fastmap_1.1.0 lazyeval_0.2.2 limma_3.48.0
## [37] tweenr_1.0.2 htmltools_0.5.1.1 tools_4.1.0
## [40] igraph_1.2.6 gtable_0.3.0 glue_1.4.2
## [43] GenomeInfoDbData_1.2.6 reshape2_1.4.4 DO.db_2.9
## [46] dplyr_1.0.6 fastmatch_1.1-0 Rcpp_1.0.6
## [49] enrichplot_1.12.0 Biobase_2.52.0 jquerylib_0.1.4
## [52] vctrs_0.3.8 Biostrings_2.60.0 babelgene_21.4
## [55] ape_5.5 nlme_3.1-152 ggraph_2.0.5
## [58] xfun_0.23 stringr_1.4.0 lifecycle_1.0.0
## [61] clusterProfiler_4.0.0 XML_3.99-0.6 DOSE_3.18.0
## [64] edgeR_3.34.0 zlibbioc_1.38.0 MASS_7.3-54
## [67] scales_1.1.1 tidygraph_1.2.0 parallel_4.1.0
## [70] RColorBrewer_1.1-2 yaml_2.2.1 memoise_2.0.0
## [73] gridExtra_2.3 downloader_0.4 sass_0.4.0
## [76] stringi_1.6.2 RSQLite_2.2.7 genefilter_1.74.0
## [79] highr_0.9 S4Vectors_0.30.0 tidytree_0.3.3
## [82] BiocGenerics_0.38.0 BiocParallel_1.26.0 GenomeInfoDb_1.28.0
## [85] matrixStats_0.58.0 rlang_0.4.11 pkgconfig_2.0.3
## [88] bitops_1.0-7 evaluate_0.14 lattice_0.20-44
## [91] purrr_0.3.4 treeio_1.16.0 patchwork_1.1.1
## [94] labeling_0.4.2 cowplot_1.1.1 shadowtext_0.0.8
## [97] bit_4.0.4 tidyselect_1.1.1 plyr_1.8.6
## [100] magrittr_2.0.1 bookdown_0.22 R6_2.5.0
## [103] IRanges_2.26.0 magick_2.7.2 generics_0.1.0
## [106] DBI_1.1.1 mgcv_1.8-35 pillar_1.6.1
## [109] withr_2.4.2 survival_3.2-11 KEGGREST_1.32.0
## [112] RCurl_1.98-1.3 tibble_3.1.2 crayon_1.4.1
## [115] utf8_1.2.1 rmarkdown_2.8 viridis_0.6.1
## [118] locfit_1.5-9.4 grid_4.1.0 sva_3.40.0
## [121] data.table_1.14.0 blob_1.2.1 digest_0.6.27
## [124] xtable_1.8-4 tidyr_1.1.3 stats4_4.1.0
## [127] munsell_0.5.0 viridisLite_0.4.0 bslib_0.2.5.1
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 Shalem1, *, 2. 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.