find_spacers: Find crispr spacers in targetranges

Description Usage Arguments Value See Also Examples

View source: R/06_find_spacers.R

Description

Find crispr spacers in targetranges

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
find_spacers(
  gr,
  bsgenome,
  spacer = strrep("N", 20),
  pam = "NGG",
  complement = TRUE,
  ontargetmethod = c("Doench2014", "Doench2016")[1],
  offtargetmethod = c("bowtie", "pdict")[1],
  offtargetfilterby = character(0),
  subtract_targets = FALSE,
  mismatches = 2,
  indexedgenomesdir = INDEXEDGENOMESDIR,
  outdir = OUTDIR,
  verbose = TRUE,
  plot = TRUE,
  ...
)

Arguments

gr

GRanges-class

bsgenome

BSgenome-class

spacer

string: spacer pattern in extended IUPAC alphabet

pam

string: pam pattern in extended IUPAC alphabet

complement

TRUE (default) or FALSE: also search in compl ranges?

ontargetmethod

'Doench2016','Doench2016' or NULL (no on-target score)

offtargetmethod

'bowtie', 'pdict', or NULL (no offtarget analysis)

offtargetfilterby

filter for best off-target counts by this variable

subtract_targets

TRUE or FALSE (default): whether to subtract target (mis)matches from offtarget counts

mismatches

0-3: allowed mismatches in offtargetanalysis (choose mismatch=-1 to suppress offtarget analysis)

indexedgenomesdir

directory with Bowtie-indexed genomes (as produced with index_genome)

outdir

directory where bowtie analysis results are written to

verbose

TRUE (default) or FALSE

plot

TRUE (default) or FALSE

...

passed to plot_intervals

Value

GRanges-class

See Also

find_primespacers to find prime editing spacers

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
# PE example
#-----------
    require(magrittr)
    bsgenome <- BSgenome.Hsapiens.UCSC.hg38::BSgenome.Hsapiens.UCSC.hg38  
    gr <- char_to_granges(c(PRNP = 'chr20:4699600:+',             # snp
                            HBB  = 'chr11:5227002:-',             # snp
                            HEXA = 'chr15:72346580-72346583:-',   # del
                            CFTR = 'chr7:117559593-117559595:+'), # ins
                          bsgenome)
    plot_intervals(gr)
    find_primespacers(gr, bsgenome)
    find_spacers(extend_for_pe(gr), bsgenome, complement=FALSE, mismatches=0)
          # complement = FALSE because extend_for_pe  already 
          # adds  reverse complements and does so in a strand-specific 
          # manner
    
# TFBS example
#-------------
    bsgenome <- BSgenome.Mmusculus.UCSC.mm10::BSgenome.Mmusculus.UCSC.mm10
    bedfile  <- system.file('extdata/SRF.bed', package='multicrispr')
    gr <- bed_to_granges(bedfile, 'mm10') %>% extend()
    gr %<>% extract(1:100)
    find_spacers(gr, bsgenome, subtract_targets = TRUE)

multicrispr documentation built on Nov. 8, 2020, 5:10 p.m.