Abstract
Genome-wide CRISPR (clustered regularly interspaced short palindrome repeats) coupled with nuclease Cas9 (CRISPR/Cas9) screens represent a promising technology to systematically evaluate gene functions. Data analysis for CRISPR/Cas9 screens is a critical process that includes identifying screen hits and exploring biological functions for these hits in downstream analysis. We have previously developed two algorithms, MAGeCK and MAGeCK-VISPR, to analyze CRISPR/Cas9 screen data in various scenarios. These two algorithms allow users to perform quality control, read count generation and normalization, and calculate beta score to evaluate gene selection performance. In downstream analysis, biological functional analysis is required for understanding biological functions of these identified genes with different screening purposes. Here, We developed MAGeCKFlute for supporting downstream analysis, utilizing the data provided through MAGeCK and MAGeCK-VISPR. MAGeCKFlute provides several strategies to remove potential biases within sgRNA-level read counts and gene-level beta scores. The downstream analysis with the package includes identifying essential, non-essential, and target-associated genes, and performing biological functional category analysis and pathway enrichment analysis of these genes. The package also visualizes genes in the context of pathways to benefit users exploring screening data. Collectively, MAGeCKFlute enables accurate identification of essential, non-essential, and targeted genes, as well as their related biological functions. This vignette explains the use of the package and demonstrates typical workflows. MAGeCKFlute package version: 1.0.0Note: if you use MAGeCKFlute in published research, please cite:
Any and all MAGeCKFlute questions should be posted to the Bioconductor support site, which serves as a searchable knowledge base of questions and answers:
https://support.bioconductor.org
Posting a question and tagging with “MAGeCKFlute” will automatically send an alert to the package authors to respond on the support site. See the first question in the list of Frequently Asked Questions (FAQ) for information about how to construct an informative post.
You should not email your question to the package authors, as we will just reply that the question should be posted to the Bioconductor support site.
As input, the MAGeCKFlute package expects gene summary file as obtained by running commands mageck test
or mageck mle
in MAGeCK (Wei Li and Liu. 2014) and MAGeCK-VISPR (Wei Li and Liu. 2015), which are developed by our lab previously, to analyze CRISPR/Cas9 screen data in different scenarios(Tim Wang 2014, Hiroko Koike-Yusa (2014), Ophir Shalem1 (2014), Luke A.Gilbert (2014), Silvana Konermann (2015)). Both algorithms use negative binomial models to model the variances of sgRNAs, and use Robust Rank Aggregation (for MAGeCK) or maximum likelihood framework (for MAGeCK-VISPR) for a robust identification of selected genes.
MAGeCK-MLE can be used to analyze CRISPR screen data from multi-conditioned experiments. MAGeCK-MLE also normalizes the data across multiple samples, making them comparable to each other. The most important ouput of MAGeCK MLE is “gene_summary” file, which includes the beta scores of multiple conditions and the associated statistics. The ‘beta score’ for each gene describes how the gene is selected: a positive beta score indicates a positive selection, and a negative beta score indicates a negative selection.
MAGeCK RRA allows for the comparison between two experimental conditions. It can identify genes and sgRNAs are significantly selected between the two conditions. The most important output of MAGeCK RRA is the file “gene_summary.txt”. MAGeCK RRA will output both the negative score and positive score for each gene. A smaller score indicates higher gene importance. MAGeCK RRA will also output the statistical value for the scores of each gene. Genes that are significantly positively and negatively selected can be identified based on the p-value or FDR.
Here we show the most basic steps for integrative analysis pipeline using gene summary file
. Before using MAGeCKFlute, analysing CRISPR/Cas9 screen data using MAGeCK RRA (in MAGeCK (Wei Li and Liu. 2014)) or MAGeCK MLE (in MAGeCK-VISPR (Wei Li and Liu. 2015)) is neccessary, which result in the generation of the gene summary file.
To run MAGeCKFlute pipeline, we need gene summary file generated by running MAGeCK RRA or MAGeCK MLE. MAGeCKFlute package provides two example data, one is MLE_Data
and the other is RRA_Data
. We will work with them in this document.
Downstream analysis pipeline for MAGeCK MLE
library(MAGeCKFlute)
##Load gene summary data in MAGeCK MLE results
data("MLE_Data")
##Run the downstream analysis pipeline for MAGeCK MLE
FluteMLE(MLE_Data, ctrlname=c("D7_R1", "D7_R2"), treatname=c("PLX7_R1","PLX7_R2"), prefix="BRAF", organism="hsa")
#All pipeline results are written into local directory
# "./BRAF_Flute_Results/", and all figures are integrated into file "BRAF_Flute.mle_summary.pdf".
Downstream analysis pipeline for MAGeCK RRA
##Load gene summary data in MAGeCK RRA results
data("RRA_Data")
##Run the downstream analysis pipeline for MAGeCK RRA
FluteRRA(RRA_Data, prefix="BRAF", organism="hsa")
#All pipeline results are written into local directory
# "./BRAF_Flute_Results/" too, and all figures are integrated into file
# "BRAF_Flute.rra_summary.pdf".
As the input of MAGeCKFlute package, the gene summary file in MAGeCK-MLE results includes beta scores of all genes in multiple condition samples. Use function ‘data’ to load the dataset, and have a look at the file with a text editor to see how it is formatted.
## Gene D7_R1.beta D7_R2.beta PLX7_R1.beta PLX7_R2.beta
## 1 SMPD2 -0.458770 -0.514690 -0.67770 -0.510360
## 2 MIXL1 0.094571 0.122780 -0.14383 -0.119670
## 3 ENPEP -0.082621 -0.082305 -0.11271 -0.061796
## 4 VAV2 -0.024270 -0.034077 0.14537 -0.015095
## 5 C5orf24 0.452310 0.226940 0.15809 0.191940
## 6 TLCD1 -0.883320 -0.779930 -0.80811 -0.854900
Then, extract beta scores of control and treatment samples from the gene summary table(can be a file path of ‘gene_summary’ or data frame), and have a look at the beta score matrix.
gene_summary=MLE_Data
ctrlname = c("D7_R1", "D7_R2")
treatname = c("PLX7_R1", "PLX7_R2")
#Read beta scores from gene summary table in MAGeCK MLE results
dd=ReadBeta(gene_summary, organism="hsa")
head(dd)
## Gene ENTREZID D7_R1 D7_R2 PLX7_R1 PLX7_R2
## 1 SMPD2 6610 -0.458770 -0.514690 -0.67770 -0.510360
## 2 MIXL1 83881 0.094571 0.122780 -0.14383 -0.119670
## 3 ENPEP 2028 -0.082621 -0.082305 -0.11271 -0.061796
## 4 VAV2 7410 -0.024270 -0.034077 0.14537 -0.015095
## 5 C5orf24 134553 0.452310 0.226940 0.15809 0.191940
## 6 TLCD1 116238 -0.883320 -0.779930 -0.80811 -0.854900
Is there batch effect? This is a common asked question before perform later analysis. In our package, we provide HeatmapView
to ensure whether batch effect exists in data, and use BatchRemove
to remove easily if same batch samples cluster together.
batchMat = data.frame(samples = c(ctrlname, treatname), batch = c(1,2,1,2), cov = c(1,1,2,2))
dd1 = BatchRemove(dd, batchMat)$data
## Standardizing Data across genes
It is difficult to control all samples with a consistent cell cycle in an CRISPR screen experiment with multi conditions. Besides, beta score among different conditions with an inconsistent cell cycle are incomparable. So it is necessary to do the normalization when comparing the beta scores in different conditions. Essential genesare thosegenes 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 cell cycle normalization method to shorten the gap of cell cycle in different conditions. In addition, a previous normalization method called loess normalization is available in this package.(Laurent Gautier 2004)
dd_essential = NormalizeBeta(dd, samples=c(ctrlname, treatname), method="cell_cycle")
head(dd_essential)
## Gene ENTREZID D7_R1 D7_R2 PLX7_R1 PLX7_R2
## 1 SMPD2 6610 -0.6456275 -0.75734255 -1.2049714 -0.90687137
## 2 MIXL1 83881 0.1330899 0.18066510 -0.2557341 -0.21264460
## 3 ENPEP 2028 -0.1162726 -0.12110800 -0.2004018 -0.10980685
## 4 VAV2 7410 -0.0341552 -0.05014273 0.2584723 -0.02682268
## 5 C5orf24 134553 0.6365364 0.33393172 0.2810889 0.34106296
## 6 TLCD1 116238 -1.2430972 -1.14763096 -1.4368444 -1.51909306
## Gene ENTREZID D7_R1 D7_R2 PLX7_R1 PLX7_R2
## 1 SMPD2 6610 -0.43634044 -0.50935627 -0.6946410 -0.52118230
## 2 MIXL1 83881 0.09494536 0.13491059 -0.1511553 -0.12484968
## 3 ENPEP 2028 -0.08082140 -0.07423802 -0.1182929 -0.06607972
## 4 VAV2 7410 -0.02141159 -0.02232013 0.1339356 -0.01827593
## 5 C5orf24 134553 0.45497915 0.25058056 0.1384579 0.18526237
## 6 TLCD1 116238 -0.83451078 -0.77556366 -0.8339421 -0.88224349
After normalization, the distribution of beta scores in different conditions should be similar. We can evaluate the distribution of beta scores using function ‘ViolinView’, ‘DensityView’, and ‘DensityDiffView’.
#we can also use function 'MAView' to evaluate the data quality of normalized
#beta score profile.
MAView(dd_essential, ctrlname, treatname, cex=1, main="Cell cycle normalized")
After normalization, the cell cycle time in different condition should be almost consistent. Here we use linear fitting to estimate the cell cycle time, and use function CellCycleView
to view the cell cycle time of all samples.
##Fitting beta score of all genes
CellCycleView(dd_essential[,-2], ctrlname, main="Cell cycle normalized")
The function ScatterView
can group all genes into three groups, positive selection genes (GroupA), negative selection genes (GroupB), and others, and visualize these three grouped genes in scatter plot. We can also use function RankView
to rank the beta score deviation between control and treatment and mark top selected genes in figure.
## Add column of 'diff'
dd_essential$Control = rowMeans(dd_essential[,ctrlname])
dd_essential$Treatment = rowMeans(dd_essential[,treatname])
dd_essential$diff = dd_essential$Treatment - dd_essential$Control
p2 = RankView(dd_essential[, c("Gene", "diff")], main="Cell cycle normalized")
print(p2)
For gene set enrichment analysis, we provide five methods in this package, including “ORT”(Over-Representing Test (Guangchuang Yu and He. 2012)), “GSEA”(Gene Set Enrichment Analysis (Aravind Subramanian and Mesirov. 2005)), “DAVID” (Huang da W 2009), “GOstats” (S. Falcon 2007), and “HGT”(HyperGemetric test). And gene ontology(GO) term (Consortium. 2014) and Kyoto encyclopedia of genes and genomes (KEGG) pathway (Minoru Kanehisa 2014) enrichment analysis are available in this package. The enrichment analysis can be done easily using function enrichment_analysis
, which return a list containing gridPlot
(ggplot object) and enrichRes
(enrichResult instance). Alternatively, you can do enrichment analysis using function enrich.ORT
for “ORT”, enrich.GSE
for GSEA, enrich.DAVID
for DAVID, enrich.GOstats
for GOstats, and enrich.HGT
for “HGT”, which return an enrichResult instance. Function EnrichedView
and EnrichedGSEView
(for enrich.GSE
) can be used to generate gridPlot
from enrichRes
easily, as shown below.
## Get information of positive and negative selection genes
groupAB = p1$data
## select positive selection genes
idx1=groupAB$group=="up"
genes=as.character(groupAB$ENTREZID[idx1])
geneList=groupAB$diff[idx1]
names(geneList)=genes
universe=as.character(groupAB$ENTREZID)
## Do enrichment analysis using HGT method
keggA = enrichment_analysis(geneList = genes, universe = universe, method = "HGT",
type = "KEGG", organism = "human", plotTitle = "Positive selection")
#same as
kegg_A = enrich.HGT(genes, universe, type = "KEGG", organism = "human")
keggA_grid = EnrichedView(kegg_A@result, plotTitle = "Positive selection")
## look at the results
head(keggA$enrichRes@result)
## ID Description pvalue p.adjust
## 191 hsa04714 Thermogenesis 4.385206e-11 1.403266e-08
## 14 hsa00190 Oxidative phosphorylation 1.431591e-09 2.290546e-07
## 249 hsa05012 Parkinson's disease 2.688353e-07 2.867577e-05
## 251 hsa05016 Huntington's disease 3.940556e-06 3.152445e-04
## 248 hsa05010 Alzheimer's disease 7.483474e-05 4.789423e-03
## 90 hsa03010 Ribosome 1.398850e-04 7.460532e-03
## nLogpvalue GeneRatio BgRatio
## 191 7.852860 69/210 210/29399
## 14 6.640061 43/115 115/29399
## 249 4.542485 40/121 121/29399
## 251 3.501352 49/174 174/29399
## 248 2.319717 42/158 158/29399
## 90 2.127230 32/114 114/29399
## geneID
## 191 9377/6195/60/84699/8110/1327/1337/4716/55744/1329/57521/196883/54539/51287/51103/6389/4708/1340/4719/4702/6602/6392/513/4717/5606/4715/5568/55967/8289/114/4713/51079/6598/4709/57104/498/107/28958/4701/6199/10632/79133/7385/29796/493753/7248/64223/65260/2778/4720/10063/7386/55471/1355/515/6605/1339/4728/4724/4718/90639/506/1353/90993/7381/1268/7249/4722/112
## 14 9377/1327/1337/155066/4716/1329/54539/6389/4708/23545/1340/4719/4702/6392/513/525/4717/4715/55967/4713/51079/4709/498/4701/10632/7385/29796/534/4720/10063/7386/1355/515/9114/1339/4728/4724/4718/526/506/1353/7381/4722
## 249 9377/1327/1337/4716/7317/1329/7332/54539/6389/4708/1340/4719/4702/6392/513/4717/4715/5568/55967/6531/4713/51079/7054/4709/498/135/4701/7385/29796/4720/7386/515/1339/4728/4724/4718/506/7381/7345/4722
## 251 9377/7019/5330/84699/581/160/1327/7157/1337/4716/127602/1329/6648/54539/6389/4708/1340/4719/4702/6392/513/2876/4717/4715/55967/6667/4713/5434/51079/1173/4709/498/4701/1387/7385/29796/6647/4720/7386/4899/515/1339/4728/4724/4718/506/90993/7381/4722
## 248 2903/9377/5330/1327/1337/4716/1329/8883/100506742/3028/54539/6389/4708/1340/4719/4702/6392/513/4717/4311/4715/55967/4713/351/51079/4709/498/4701/9451/7385/29796/2906/4720/7386/515/1339/4728/4724/4718/506/7381/4722
## 90 51021/51116/6230/6208/64979/63931/124995/51318/29093/64968/51264/29088/28998/6183/6191/55052/64981/6193/64969/6234/6181/63875/55168/6182/6150/6146/11222/6209/6169/64963/54948/51073
## geneName
## 191 COX5A/RPS6KA1/ACTB/CREB3L3/DPF3/COX4I1/COX6A1/NDUFB10/COA1/COX5B/RPTOR/ADCY4/NDUFB11/COA4/NDUFAF1/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SMARCD1/SDHD/ATP5F1D/NDUFC1/MAP2K3/NDUFB9/PRKACG/NDUFA12/ARID1A/ADCY8/NDUFB7/NDUFA13/SMARCB1/NDUFB3/PNPLA2/ATP5F1A/ADCY1/COA3/NDUFA7/RPS6KB2/ATP5MG/NDUFAF5/UQCRC2/UQCR10/COA5/TSC1/MLST8/COA7/GNAS/NDUFS2/COX17/UQCRFS1/NDUFAF7/COX15/ATP5PB/SMARCE1/COX6A2/NDUFS8/NDUFS4/NDUFC2/COX19/ATP5F1B/COX11/CREB3L1/UQCRB/CNR1/TSC2/NDUFS3/ADCY6
## 14 COX5A/COX4I1/COX6A1/ATP6V0E2/NDUFB10/COX5B/NDUFB11/SDHA/NDUFB2/ATP6V0A2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/ATP6V1B1/NDUFC1/NDUFB9/NDUFA12/NDUFB7/NDUFA13/NDUFB3/ATP5F1A/NDUFA7/ATP5MG/UQCRC2/UQCR10/ATP6V1G2/NDUFS2/COX17/UQCRFS1/COX15/ATP5PB/ATP6V0D1/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP6V1B2/ATP5F1B/COX11/UQCRB/NDUFS3
## 249 COX5A/COX4I1/COX6A1/NDUFB10/UBA1/COX5B/UBE2L3/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/NDUFC1/NDUFB9/PRKACG/NDUFA12/SLC6A3/NDUFB7/NDUFA13/TH/NDUFB3/ATP5F1A/ADORA2A/NDUFA7/UQCRC2/UQCR10/NDUFS2/UQCRFS1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/UQCRB/UCHL1/NDUFS3
## 251 COX5A/TFAM/PLCB2/CREB3L3/BAX/AP2A1/COX4I1/TP53/COX6A1/NDUFB10/DNAH14/COX5B/SOD2/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/GPX1/NDUFC1/NDUFB9/NDUFA12/SP1/NDUFB7/POLR2E/NDUFA13/AP2M1/NDUFB3/ATP5F1A/NDUFA7/CREBBP/UQCRC2/UQCR10/SOD1/NDUFS2/UQCRFS1/NRF1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/CREB3L1/UQCRB/NDUFS3
## 248 GRIN2A/COX5A/PLCB2/COX4I1/COX6A1/NDUFB10/COX5B/NAE1/CASP12/HSD17B10/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/NDUFC1/MME/NDUFB9/NDUFA12/NDUFB7/APP/NDUFA13/NDUFB3/ATP5F1A/NDUFA7/EIF2AK3/UQCRC2/UQCR10/GRIN2D/NDUFS2/UQCRFS1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/UQCRB/NDUFS3
## 90 MRPS16/MRPS2/RPS25/RPS14/MRPL36/MRPS14/MRPL10/MRPL35/MRPL22/MRPS6/MRPL27/MRPL15/MRPL13/MRPS12/RPS4X/MRPL20/MRPL34/RPS5/MRPS5/RPS28/RPLP2/MRPL17/MRPS18A/MRPL12/MRPL23/RPL22/MRPL3/RPS15/RPL38/MRPS11/MRPL16/MRPL4
## Count
## 191 69
## 14 43
## 249 40
## 251 49
## 248 42
## 90 32
## ID Description pvalue p.adjust
## 191 hsa04714 Thermogenesis 4.385206e-11 1.403266e-08
## 14 hsa00190 Oxidative phosphorylation 1.431591e-09 2.290546e-07
## 249 hsa05012 Parkinson's disease 2.688353e-07 2.867577e-05
## 251 hsa05016 Huntington's disease 3.940556e-06 3.152445e-04
## 248 hsa05010 Alzheimer's disease 7.483474e-05 4.789423e-03
## 90 hsa03010 Ribosome 1.398850e-04 7.460532e-03
## nLogpvalue GeneRatio BgRatio
## 191 7.852860 69/210 210/29399
## 14 6.640061 43/115 115/29399
## 249 4.542485 40/121 121/29399
## 251 3.501352 49/174 174/29399
## 248 2.319717 42/158 158/29399
## 90 2.127230 32/114 114/29399
## geneID
## 191 9377/6195/60/84699/8110/1327/1337/4716/55744/1329/57521/196883/54539/51287/51103/6389/4708/1340/4719/4702/6602/6392/513/4717/5606/4715/5568/55967/8289/114/4713/51079/6598/4709/57104/498/107/28958/4701/6199/10632/79133/7385/29796/493753/7248/64223/65260/2778/4720/10063/7386/55471/1355/515/6605/1339/4728/4724/4718/90639/506/1353/90993/7381/1268/7249/4722/112
## 14 9377/1327/1337/155066/4716/1329/54539/6389/4708/23545/1340/4719/4702/6392/513/525/4717/4715/55967/4713/51079/4709/498/4701/10632/7385/29796/534/4720/10063/7386/1355/515/9114/1339/4728/4724/4718/526/506/1353/7381/4722
## 249 9377/1327/1337/4716/7317/1329/7332/54539/6389/4708/1340/4719/4702/6392/513/4717/4715/5568/55967/6531/4713/51079/7054/4709/498/135/4701/7385/29796/4720/7386/515/1339/4728/4724/4718/506/7381/7345/4722
## 251 9377/7019/5330/84699/581/160/1327/7157/1337/4716/127602/1329/6648/54539/6389/4708/1340/4719/4702/6392/513/2876/4717/4715/55967/6667/4713/5434/51079/1173/4709/498/4701/1387/7385/29796/6647/4720/7386/4899/515/1339/4728/4724/4718/506/90993/7381/4722
## 248 2903/9377/5330/1327/1337/4716/1329/8883/100506742/3028/54539/6389/4708/1340/4719/4702/6392/513/4717/4311/4715/55967/4713/351/51079/4709/498/4701/9451/7385/29796/2906/4720/7386/515/1339/4728/4724/4718/506/7381/4722
## 90 51021/51116/6230/6208/64979/63931/124995/51318/29093/64968/51264/29088/28998/6183/6191/55052/64981/6193/64969/6234/6181/63875/55168/6182/6150/6146/11222/6209/6169/64963/54948/51073
## geneName
## 191 COX5A/RPS6KA1/ACTB/CREB3L3/DPF3/COX4I1/COX6A1/NDUFB10/COA1/COX5B/RPTOR/ADCY4/NDUFB11/COA4/NDUFAF1/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SMARCD1/SDHD/ATP5F1D/NDUFC1/MAP2K3/NDUFB9/PRKACG/NDUFA12/ARID1A/ADCY8/NDUFB7/NDUFA13/SMARCB1/NDUFB3/PNPLA2/ATP5F1A/ADCY1/COA3/NDUFA7/RPS6KB2/ATP5MG/NDUFAF5/UQCRC2/UQCR10/COA5/TSC1/MLST8/COA7/GNAS/NDUFS2/COX17/UQCRFS1/NDUFAF7/COX15/ATP5PB/SMARCE1/COX6A2/NDUFS8/NDUFS4/NDUFC2/COX19/ATP5F1B/COX11/CREB3L1/UQCRB/CNR1/TSC2/NDUFS3/ADCY6
## 14 COX5A/COX4I1/COX6A1/ATP6V0E2/NDUFB10/COX5B/NDUFB11/SDHA/NDUFB2/ATP6V0A2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/ATP6V1B1/NDUFC1/NDUFB9/NDUFA12/NDUFB7/NDUFA13/NDUFB3/ATP5F1A/NDUFA7/ATP5MG/UQCRC2/UQCR10/ATP6V1G2/NDUFS2/COX17/UQCRFS1/COX15/ATP5PB/ATP6V0D1/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP6V1B2/ATP5F1B/COX11/UQCRB/NDUFS3
## 249 COX5A/COX4I1/COX6A1/NDUFB10/UBA1/COX5B/UBE2L3/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/NDUFC1/NDUFB9/PRKACG/NDUFA12/SLC6A3/NDUFB7/NDUFA13/TH/NDUFB3/ATP5F1A/ADORA2A/NDUFA7/UQCRC2/UQCR10/NDUFS2/UQCRFS1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/UQCRB/UCHL1/NDUFS3
## 251 COX5A/TFAM/PLCB2/CREB3L3/BAX/AP2A1/COX4I1/TP53/COX6A1/NDUFB10/DNAH14/COX5B/SOD2/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/GPX1/NDUFC1/NDUFB9/NDUFA12/SP1/NDUFB7/POLR2E/NDUFA13/AP2M1/NDUFB3/ATP5F1A/NDUFA7/CREBBP/UQCRC2/UQCR10/SOD1/NDUFS2/UQCRFS1/NRF1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/CREB3L1/UQCRB/NDUFS3
## 248 GRIN2A/COX5A/PLCB2/COX4I1/COX6A1/NDUFB10/COX5B/NAE1/CASP12/HSD17B10/NDUFB11/SDHA/NDUFB2/COX6B1/NDUFS1/NDUFA8/SDHD/ATP5F1D/NDUFC1/MME/NDUFB9/NDUFA12/NDUFB7/APP/NDUFA13/NDUFB3/ATP5F1A/NDUFA7/EIF2AK3/UQCRC2/UQCR10/GRIN2D/NDUFS2/UQCRFS1/ATP5PB/COX6A2/NDUFS8/NDUFS4/NDUFC2/ATP5F1B/UQCRB/NDUFS3
## 90 MRPS16/MRPS2/RPS25/RPS14/MRPL36/MRPS14/MRPL10/MRPL35/MRPL22/MRPS6/MRPL27/MRPL15/MRPL13/MRPS12/RPS4X/MRPL20/MRPL34/RPS5/MRPS5/RPS28/RPLP2/MRPL17/MRPS18A/MRPL12/MRPL23/RPL22/MRPL3/RPS15/RPL38/MRPS11/MRPL16/MRPL4
## Count
## 191 69
## 14 43
## 249 40
## 251 49
## 248 42
## 90 32
## Do enrichment analysis using GSEA method
gseA = enrichment_analysis(geneList = geneList, method = "GSEA", type = "KEGG",
organism="human", plotTitle="Positive selection")
#same as
gse_A = enrich.GSE(geneList, type = "KEGG", organism = "human")
gseA_grid = EnrichedGSEView(gse_A@result, plotTitle = "Positive selection")
head(gseA$enrichRes@result)
## ID Description setSize enrichmentScore
## hsa05016 hsa05016 Huntington's disease 49 0.4350534
## hsa00190 hsa00190 Oxidative phosphorylation 43 0.4691172
## hsa05010 hsa05010 Alzheimer's disease 42 0.4802060
## hsa05012 hsa05012 Parkinson's disease 40 0.4806812
## hsa04120 hsa04120 Ubiquitin mediated proteolysis 30 0.4902955
## hsa05169 hsa05169 Epstein-Barr virus infection 28 0.4848621
## NES pvalue p.adjust qvalues rank
## hsa05016 1.919289 0.001019368 0.05402055 0.05145702 625
## hsa00190 2.019095 0.001024590 0.05402055 0.05145702 625
## hsa05010 2.057331 0.001025641 0.05402055 0.05145702 625
## hsa05012 2.027985 0.001030928 0.05402055 0.05145702 625
## hsa04120 1.948030 0.001048218 0.05402055 0.05145702 186
## hsa05169 1.894961 0.001062699 0.05402055 0.05145702 523
## leading_edge
## hsa05016 tags=49%, list=24%, signal=38%
## hsa00190 tags=53%, list=24%, signal=42%
## hsa05010 tags=55%, list=24%, signal=42%
## hsa05012 tags=55%, list=24%, signal=43%
## hsa04120 tags=27%, list=7%, signal=25%
## hsa05169 tags=50%, list=20%, signal=41%
## core_enrichment
## hsa05016 7019/1387/1327/7381/9377/4715/6392/160/4702/4708/506/4701/5434/498/1329/513/29796/6389/7385/1337/7386/4709/4719/1340
## hsa00190 10063/1327/7381/9377/4715/6392/4702/10632/4708/1355/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/4719/1340
## hsa05010 3028/1327/8883/7381/9377/4715/6392/4702/4708/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/2906/4719/1340
## hsa05012 1327/7381/9377/4715/6392/4702/7332/4708/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/135/4719/1340
## hsa04120 8452/9040/6923/9039/9616/8065/7332/9817
## hsa05169 1460/1387/7979/5608/3305/5704/5434/3312/7186/3106/171568/6693/51728/10621
## geneName
## hsa05016 TFAM/CREBBP/COX4I1/UQCRB/COX5A/NDUFB9/SDHD/AP2A1/NDUFA8/NDUFB2/ATP5F1B/NDUFA7/POLR2E/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/NDUFS1/COX6B1
## hsa00190 COX17/COX4I1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/ATP5MG/NDUFB2/COX15/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/NDUFS1/COX6B1
## hsa05010 HSD17B10/COX4I1/NAE1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/NDUFB2/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/GRIN2D/NDUFS1/COX6B1
## hsa05012 COX4I1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/UBE2L3/NDUFB2/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/ADORA2A/NDUFS1/COX6B1
## hsa04120 CUL3/UBE2M/ELOB/UBA3/RNF7/CUL5/UBE2L3/KEAP1
## hsa05169 CSNK2B/CREBBP/SEM1/MAP2K6/HSPA1L/PSMC4/POLR2E/HSPA8/TRAF2/HLA-B/POLR3H/SPN/POLR3K/POLR3F
## ID Description setSize enrichmentScore
## hsa00190 hsa00190 Oxidative phosphorylation 43 0.4691172
## hsa05010 hsa05010 Alzheimer's disease 42 0.4802060
## hsa04714 hsa04714 Thermogenesis 69 0.4110571
## hsa05012 hsa05012 Parkinson's disease 40 0.4806812
## hsa05016 hsa05016 Huntington's disease 49 0.4350534
## hsa04120 hsa04120 Ubiquitin mediated proteolysis 30 0.4902955
## NES pvalue p.adjust qvalues rank
## hsa00190 1.986966 0.001018330 0.07790224 0.06860328 625
## hsa05010 2.029919 0.001018330 0.07790224 0.06860328 625
## hsa04714 1.877031 0.002012072 0.07830092 0.06895437 625
## hsa05012 2.008690 0.002047083 0.07830092 0.06895437 625
## hsa05016 1.878569 0.003039514 0.09211318 0.08111791 625
## hsa04120 1.934879 0.004214963 0.09211318 0.08111791 186
## leading_edge
## hsa00190 tags=53%, list=24%, signal=42%
## hsa05010 tags=55%, list=24%, signal=42%
## hsa04714 tags=46%, list=24%, signal=36%
## hsa05012 tags=55%, list=24%, signal=43%
## hsa05016 tags=49%, list=24%, signal=38%
## hsa04120 tags=27%, list=7%, signal=25%
## core_enrichment
## hsa00190 10063/1327/7381/9377/4715/6392/4702/10632/4708/1355/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/4719/1340
## hsa05010 3028/1327/8883/7381/9377/4715/6392/4702/4708/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/2906/4719/1340
## hsa04714 6598/10063/65260/1327/7381/9377/4715/6392/79133/4702/10632/4708/1355/506/4701/8110/498/1329/513/29796/6389/7385/6199/1337/7386/60/55744/8289/4709/51103/4719/1340
## hsa05012 1327/7381/9377/4715/6392/4702/7332/4708/506/4701/498/1329/513/29796/6389/7385/1337/7386/4709/135/4719/1340
## hsa05016 7019/1387/1327/7381/9377/4715/6392/160/4702/4708/506/4701/5434/498/1329/513/29796/6389/7385/1337/7386/4709/4719/1340
## hsa04120 8452/9040/6923/9039/9616/8065/7332/9817
## geneName
## hsa00190 COX17/COX4I1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/ATP5MG/NDUFB2/COX15/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/NDUFS1/COX6B1
## hsa05010 HSD17B10/COX4I1/NAE1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/NDUFB2/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/GRIN2D/NDUFS1/COX6B1
## hsa04714 SMARCB1/COX17/COA7/COX4I1/UQCRB/COX5A/NDUFB9/SDHD/NDUFAF5/NDUFA8/ATP5MG/NDUFB2/COX15/ATP5F1B/NDUFA7/DPF3/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/RPS6KB2/COX6A1/UQCRFS1/ACTB/COA1/ARID1A/NDUFB3/NDUFAF1/NDUFS1/COX6B1
## hsa05012 COX4I1/UQCRB/COX5A/NDUFB9/SDHD/NDUFA8/UBE2L3/NDUFB2/ATP5F1B/NDUFA7/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/ADORA2A/NDUFS1/COX6B1
## hsa05016 TFAM/CREBBP/COX4I1/UQCRB/COX5A/NDUFB9/SDHD/AP2A1/NDUFA8/NDUFB2/ATP5F1B/NDUFA7/POLR2E/ATP5F1A/COX5B/ATP5F1D/UQCR10/SDHA/UQCRC2/COX6A1/UQCRFS1/NDUFB3/NDUFS1/COX6B1
## hsa04120 CUL3/UBE2M/ELOB/UBA3/RNF7/CUL5/UBE2L3/KEAP1
For enriched pathways, we can use function KeggPathwayView
to visualize the beta score level in control and treatment on pathway map.(Weijun Luo 2013)
genedata = dd_essential[,c("Control","Treatment")]
rownames(genedata)=dd_essential$ENTREZID
keggID = keggA$enrichRes@result$ID[1]
#The pathway map will be located on current workspace
KeggPathwayView(gene.data = genedata, pathway.id = keggID, species="hsa")
##Read the figure into R
pngname=paste0(keggID, ".pathview.multi.png")
grid.arrange(grid::rasterGrob(png::readPNG(pngname)))
## [1] TRUE
Considering the difference of beta scores in control and treatment sample, we developed a 9-square model, which group all genes into several subgroups. Among these subgroups, four subgroup genes are treatment-associated, which correspond to specific functions. Group1 and Group3 genes are not selected in control sample, while they are significantly selected in treatment sample, so they may related to drug resistance. Group2 and Group4 genes are selected in control, but they are not selected in treatment, so maybe these genes are associated with drug targets.
In this package, use function SquareView
can select these four treatment-associated subgroup genes and view them in 9-Square scatter plot.
Same as section above. We can do enrichment analysis for treatment-associated genes.
##Get information of treatment-associated genes
Square9 = p3$data
##==select group1 genes in 9-Square
idx=Square9$group=="Group1"
genes=as.character(Square9$ENTREZID[idx])
universe=as.character(Square9$ENTREZID)
#====KEGG_enrichment=====
kegg1=enrichment_analysis(geneList = genes, universe = universe,
type = "KEGG",plotTitle = "KEGG: Group1")
## look at the results
head(kegg1$enrichRes@result)
## ID Description GeneRatio BgRatio
## hsa00190 hsa00190 Oxidative phosphorylation 19/324 107/6243
## hsa03010 hsa03010 Ribosome 15/324 71/6243
## hsa04714 hsa04714 Thermogenesis 24/324 204/6243
## hsa05016 hsa05016 Huntington's disease 20/324 160/6243
## hsa05010 hsa05010 Alzheimer's disease 19/324 152/6243
## hsa05012 hsa05012 Parkinson's disease 14/324 116/6243
## pvalue p.adjust qvalue
## hsa00190 1.899049e-06 0.000330253 0.000330253
## hsa03010 2.492476e-06 0.000330253 0.000330253
## hsa04714 1.275775e-04 0.011269343 0.011269343
## hsa05016 2.066394e-04 0.013689860 0.013689860
## hsa05010 2.971804e-04 0.015750560 0.015750560
## hsa05012 2.559705e-03 0.113053640 0.113053640
## geneID
## hsa00190 10063/7381/4708/1340/4701/1355/526/515/23545/4722/6390/4695/529/51079/498/6392/29796/126328/55967
## hsa03010 63875/64979/6182/51116/55052/124995/6230/6183/54948/6155/51021/65008/51069/6150/2197
## hsa04714 10063/7381/4708/6598/1340/65260/4701/1355/7249/515/4722/493753/6390/4695/90993/51079/498/6392/29796/126328/51287/55967/57521/79133
## hsa05016 7381/4708/1340/6648/4701/515/5436/4722/160/6390/4695/90993/51079/498/6392/29796/1767/5330/126328/55967
## hsa05010 3028/7381/4708/8883/1340/4701/515/4722/6390/4695/51079/498/6392/29796/2906/5330/126328/55967/4311
## hsa05012 7381/4708/1340/4701/515/4722/6390/4695/51079/498/6392/29796/126328/55967
## Count
## hsa00190 19
## hsa03010 15
## hsa04714 24
## hsa05016 20
## hsa05010 19
## hsa05012 14
## geneName
## hsa00190 COX17/UQCRB/NDUFB2/COX6B1/NDUFA7/COX15/ATP6V1B2/ATP5PB/ATP6V0A2/NDUFS3/SDHB/NDUFA2/ATP6V1E1/NDUFA13/ATP5F1A/SDHD/UQCR10/NDUFA11/NDUFA12
## hsa03010 MRPL17/MRPL36/MRPL12/MRPS2/MRPL20/MRPL10/RPS25/MRPS12/MRPL16/RPL27/MRPS16/MRPL1/MRPL2/MRPL23/FAU
## hsa04714 COX17/UQCRB/NDUFB2/SMARCB1/COX6B1/COA7/NDUFA7/COX15/TSC2/ATP5PB/NDUFS3/COA5/SDHB/NDUFA2/CREB3L1/NDUFA13/ATP5F1A/SDHD/UQCR10/NDUFA11/COA4/NDUFA12/RPTOR/NDUFAF5
## hsa05016 UQCRB/NDUFB2/COX6B1/SOD2/NDUFA7/ATP5PB/POLR2G/NDUFS3/AP2A1/SDHB/NDUFA2/CREB3L1/NDUFA13/ATP5F1A/SDHD/UQCR10/DNAH5/PLCB2/NDUFA11/NDUFA12
## hsa05010 HSD17B10/UQCRB/NDUFB2/NAE1/COX6B1/NDUFA7/ATP5PB/NDUFS3/SDHB/NDUFA2/NDUFA13/ATP5F1A/SDHD/UQCR10/GRIN2D/PLCB2/NDUFA11/NDUFA12/MME
## hsa05012 UQCRB/NDUFB2/COX6B1/NDUFA7/ATP5PB/NDUFS3/SDHB/NDUFA2/NDUFA13/ATP5F1A/SDHD/UQCR10/NDUFA11/NDUFA12
Also, pathway visualization can be done using function KeggPathwayView
, the same as section above.
genedata = dd_essential[, c("Control","Treatment")]
rownames(genedata) = dd_essential$ENTREZID
keggID = kegg1$enrichRes@result$ID[1]
KeggPathwayView(gene.data = genedata, pathway.id = keggID, species="hsa")
##Read the figure into R
pngname=paste0(keggID, ".pathview.multi.png")
grid.arrange(grid::rasterGrob(png::readPNG(pngname)))
## [1] TRUE
For experiments with two experimental conditions, we recommend to use MAGeCK-RRA to identify essential genes from CRISPR/Cas9 knockout screens and tests the statistical significance of each observed change between two conditions. Gene summary file in MAGeCK-RRA results summarize the statistical significance of positive selection and negative selection. Use function ‘data’ to load the dataset, and have a look at the file with a text editor to see how it is formatted.
## id neg.fdr pos.fdr
## 1 CMIP 0.001650 0.977480
## 2 MCL1 0.001650 1.000000
## 3 ITGB2 0.001650 0.999967
## 4 SCN2A 0.003713 0.886655
## 5 GRB2 0.048680 0.984303
## 6 EGFR 0.048680 0.993803
Then, extract “neg.fdr” and “pos.fdr” from the gene summary table.
## Official neg.fdr pos.fdr ENTREZID
## 1 CMIP 0.001650 0.977480 80790
## 2 MCL1 0.001650 1.000000 4170
## 3 ITGB2 0.001650 0.999967 3689
## 4 SCN2A 0.003713 0.886655 6326
## 5 GRB2 0.048680 0.984303 2885
## 6 EGFR 0.048680 0.993803 1956
Take 0.05 as cutoff, get negative selection and positive selection genes and do enrichment analysis on KEGG pathway and GO BP terms.
##Negative selection
universe=dd.rra$ENTREZID
genes = dd.rra[dd.rra$neg.fdr<0.05, "ENTREZID"]
kegg.neg=enrichment_analysis(geneList = genes, universe=universe,
type = "KEGG", plotTitle="KEGG: neg")
bp.neg=enrichment_analysis(geneList = genes, universe=universe,
type = "BP", plotTitle="BP: neg")
print(kegg.neg$gridPlot)
print(bp.neg$gridPlot)
##Positive selection
universe=dd.rra$ENTREZID
genes = dd.rra[dd.rra$pos.fdr<0.05, "ENTREZID"]
kegg.pos=enrichment_analysis(geneList = genes, universe=universe,
type = "KEGG", plotTitle="KEGG: pos")
bp.pos=enrichment_analysis(geneList = genes, universe=universe,
type = "BP", plotTitle="BP: pos")
print(kegg.pos$gridPlot)
print(bp.pos$gridPlot)
## R version 3.5.0 (2018-04-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.7-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.7-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] MAGeCKFlute_1.0.0 gridExtra_2.3 pathview_1.20.0
## [4] org.Hs.eg.db_3.6.0 AnnotationDbi_1.42.0 IRanges_2.14.0
## [7] S4Vectors_0.18.0 Biobase_2.40.0 BiocGenerics_0.26.0
## [10] ggplot2_2.2.1
##
## loaded via a namespace (and not attached):
## [1] fgsea_1.6.0 assertive.base_0.0-7
## [3] colorspace_1.3-2 ggridges_0.5.0
## [5] rprojroot_1.3-2 qvalue_2.12.0
## [7] XVector_0.20.0 ggrepel_0.7.0
## [9] bit64_0.9-7 codetools_0.2-15
## [11] splines_3.5.0 GOSemSim_2.6.0
## [13] knitr_1.20 annotate_1.58.0
## [15] GO.db_3.6.0 png_0.1-7
## [17] pheatmap_1.0.8 pathological_0.1-2
## [19] graph_1.58.0 ggforce_0.1.1
## [21] shiny_1.0.5 compiler_3.5.0
## [23] httr_1.3.1 rvcheck_0.0.9
## [25] GOstats_2.46.0 backports_1.1.2
## [27] assertthat_0.2.0 Matrix_1.2-14
## [29] lazyeval_0.2.1 limma_3.36.0
## [31] later_0.7.1 tweenr_0.1.5
## [33] htmltools_0.3.6 tools_3.5.0
## [35] bindrcpp_0.2.2 igraph_1.2.1
## [37] gtable_0.2.0 glue_1.2.0
## [39] Category_2.46.0 reshape2_1.4.3
## [41] DO.db_2.9 dplyr_0.7.4
## [43] fastmatch_1.1-0 Rcpp_0.12.16
## [45] enrichplot_1.0.0 Biostrings_2.48.0
## [47] nlme_3.1-137 udunits2_0.13
## [49] assertive.files_0.0-2 ggraph_1.0.1
## [51] stringr_1.3.0 mime_0.5
## [53] miniUI_0.1.1 clusterProfiler_3.8.0
## [55] XML_3.98-1.11 DOSE_3.6.0
## [57] zlibbioc_1.26.0 MASS_7.3-50
## [59] scales_0.5.0 promises_1.0.1
## [61] RBGL_1.56.0 KEGGgraph_1.40.0
## [63] RColorBrewer_1.1-2 assertive.strings_0.0-3
## [65] yaml_2.1.18 memoise_1.1.0
## [67] UpSetR_1.3.3 reshape_0.8.7
## [69] ggExtra_0.8 stringi_1.1.7
## [71] RSQLite_2.1.0 genefilter_1.62.0
## [73] BiocParallel_1.14.0 matrixStats_0.53.1
## [75] rlang_0.2.0 pkgconfig_2.0.1
## [77] bitops_1.0-6 evaluate_0.10.1
## [79] lattice_0.20-35 purrr_0.2.4
## [81] bindr_0.1.1 labeling_0.3
## [83] assertive.properties_0.0-4 cowplot_0.9.2
## [85] bit_1.1-12 GSEABase_1.42.0
## [87] AnnotationForge_1.22.0 ggsci_2.8
## [89] plyr_1.8.4 magrittr_1.5
## [91] R6_2.2.2 DBI_0.8
## [93] mgcv_1.8-23 pillar_1.2.2
## [95] assertive.numbers_0.0-2 units_0.5-1
## [97] survival_2.42-3 KEGGREST_1.20.0
## [99] RCurl_1.95-4.10 tibble_1.4.2
## [101] assertive.types_0.0-3 rmarkdown_1.9
## [103] viridis_0.5.1 grid_3.5.0
## [105] sva_3.28.0 data.table_1.10.4-3
## [107] blob_1.1.1 Rgraphviz_2.24.0
## [109] digest_0.6.15 xtable_1.8-2
## [111] tidyr_0.8.0 httpuv_1.4.1
## [113] munsell_0.4.3 viridisLite_0.3.0
## [115] assertive.reflection_0.0-4
Aravind Subramanian, Vamsi K. Moothaa, Pablo Tamayo, and Jill P. Mesirov. 2005. “Gene set enrichment analysis: A knowledge-based approach for interpreting genome-wide expression profiles.” http://www.pnas.org/content/102/43/15545.full.
Consortium., The Gene Ontology. 2014. “Gene Ontology Consortium: going forward.” https://academic.oup.com/nar/article/43/D1/D1049/2439067.
Guangchuang Yu, Yanyan Han, Li-Gen Wang, and Qing-Yu He. 2012. “clusterProfiler: an R Package for Comparing Biological Themes Among Gene Clusters.” http://online.liebertpub.com/doi/abs/10.1089/omi.2011.0118.
Hiroko Koike-Yusa, E-Pien Tan, Yilong Li. 2014. “Genome-wide recessive genetic screening in mammalian cells with a lentiviral CRISPR-guide RNA library.” http://science.sciencemag.org/content/343/6166/80.long.
Huang da W, Lempicki RA., Sherman BT. 2009. “Bioinformatics enrichment tools: paths toward the comprehensive functional analysis of large gene lists.” https://academic.oup.com/nar/article-lookup/doi/10.1093/nar/gkn923.
Laurent Gautier, Benjamin M. Bolstad, Leslie Cope. 2004. “affy—analysis of Affymetrix GeneChip data at the probe level.” https://academic.oup.com/bioinformatics/article/20/3/307/185980.
Luke A.Gilbert, BrittAdamson, Max A.Horlbeck. 2014. “Genome-Scale CRISPR-Mediated Control of Gene Repression and Activation.” https://linkinghub.elsevier.com/retrieve/pii/S0092-8674(14)01178-7.
Minoru Kanehisa, Yoko Sato, Susumu Goto. 2014. “Data, information, knowledge and principle: back to metabolism in KEGG.” https://academic.oup.com/nar/article-lookup/doi/10.1093/nar/gkt1076.
Ophir Shalem1, *, 2. 2014. “Genome-scale CRISPR-Cas9 knockout screening in human cells.” http://science.sciencemag.org/content/343/6166/84.long.
S. Falcon, R. Gentleman. 2007. “Using GOstats to test gene lists for GO term association.” https://academic.oup.com/bioinformatics/article/23/2/257/204776.
Silvana Konermann, Alexandro E. Trevino, Mark D. Brigham. 2015. “Genome-scale transcriptional activation by an engineered CRISPR-Cas9 complex.” https://www.nature.com/nature/journal/vnfv/ncurrent/full/nature14136.html.
Tim Wang, David M. Sabatini, Jenny J. Wei1. 2014. “Genetic Screens in Human Cells Using the CRISPR-Cas9 System.” http://science.sciencemag.org/content/343/6166/80.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.
Weijun Luo, Cory Brouwer. 2013. “Pathview: an R/Bioconductor package for pathway-based data integration and visualization.” https://academic.oup.com/bioinformatics/article-lookup/doi/10.1093/bioinformatics/btt285.