This is an example workflow for processing a pooled CRISPR-based screen using the provided sample data. See the various manpages for additional visualization options and algorithmic details.
library(Biobase)
library(limma)
library(gCrisprTools)
data("es", package = "gCrisprTools")
data("ann", package = "gCrisprTools")
data("aln", package = "gCrisprTools")sk <- relevel(as.factor(pData(es)$TREATMENT_NAME), "ControlReference")
names(sk) <- row.names(pData(es))design <- model.matrix(~ 0 + REPLICATE_POOL + TREATMENT_NAME, pData(es))
colnames(design) <- gsub('TREATMENT_NAME', '', colnames(design))
contrasts <-makeContrasts(DeathExpansion - ControlExpansion, levels = design)es <- ct.filterReads(es, trim = 1000, sampleKey = sk)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)ct.filterReads aboveann <- ct.prepareAnnotation(ann, fit, controls = "NoTarget")resultsDF <-
  ct.generateResults(
    fit,
    annotation = ann,
    RRAalphaCutoff = 0.1,
    permutations = 1000,
    scoring = "combined", 
    permutation.seed = 2
  )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)),]ct.alignmentChart(aln, sk)
ct.rawCountDensities(es, sk)ct.gRNARankByReplicate(es, sk) 
ct.gRNARankByReplicate(es, sk, annotation = ann, geneSymb = "NoTarget")  #Show locations of NTC gRNAsct.viewControls(es, ann, sk, normalize = FALSE)
ct.viewControls(es, ann, sk, normalize = TRUE)ct.GCbias(es, ann, sk)
ct.GCbias(fit, ann, sk)ct.stackGuides(es,
               sk,
               plotType = "gRNA",
               annotation = ann,
               nguides = 40)ct.stackGuides(es, 
               sk, 
               plotType = "Target", 
               annotation = ann)ct.stackGuides(es,
               sk,
               plotType = "Target",
               annotation = ann,
               subset = names(sk)[grep('Expansion', sk)])ct.guideCDF(es, sk, plotType = "gRNA")
ct.guideCDF(es, sk, plotType = "Target", annotation = ann)ct.topTargets(fit,
              resultsDF,
              ann,
              targets = 10,
              enrich = TRUE)
ct.topTargets(fit,
              resultsDF,
              ann,
              targets = 10,
              enrich = FALSE)ct.viewGuides("Target1633", fit, ann)
ct.gRNARankByReplicate(es, sk, annotation = ann, geneSymb = "Target1633")enrichmentResults <-
  ct.PantherPathwayEnrichment(
    resultsDF,
    pvalue.cutoff = 0.01,
    enrich = TRUE,
    organism = 'mouse'
  )data("essential.genes", package = "gCrisprTools")
ROCs <- ct.ROC(resultsDF, essential.genes, stat = "deplete.p")
PRCs <- ct.PRC(resultsDF, essential.genes, stat = "deplete.p")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')