find_primespacers: Find prime editing spacers

Description Usage Arguments Details Value See Also Examples

Description

Find prime editing spacers around target ranges

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
find_primespacers(
  gr,
  bsgenome,
  edits = get_plus_seq(bsgenome, gr),
  nprimer = 13,
  nrt = 16,
  ontargetmethod = c("Doench2014", "Doench2016")[1],
  offtargetmethod = c("bowtie", "pdict")[1],
  mismatches = 0,
  nickmatches = 2,
  indexedgenomesdir = INDEXEDGENOMESDIR,
  outdir = OUTDIR,
  verbose = TRUE,
  plot = TRUE,
  ...
)

Arguments

gr

GRanges-class

bsgenome

BSgenome-class

edits

character vector: desired edits on '+' strand. If named, names should be identical to those of gr

nprimer

n primer nucleotides (default 13, max 17)

nrt

n rev transcr nucleotides (default 16, recomm. 10-16)

ontargetmethod

'Doench2014' or 'Doench2016': on-target scoring method

offtargetmethod

'bowtie' or 'pdict'

mismatches

no of primespacer mismatches (default 0, to suppress offtarget analysis: -1)

nickmatches

no of nickspacer offtarget mismatches (default 2, to suppresses offtarget analysis: -1)

indexedgenomesdir

directory with indexed genomes (as created by index_genome)

outdir

directory whre offtarget analysis output is written

verbose

TRUE (default) or FALSE

plot

TRUE (default) or FALSE

...

passed to plot_intervals

Details

Below the architecture of a prime editing site. Edits can be performed anywhere in the revtranscript area.

spacer pam ——————–=== primer revtranscript ————-================ 1..............17....GG.......... .....................CC.......... ———-extension———-

Value

GRanges-class with prime editing spacer ranges and following mcols: * crisprspacer: N20 spacers * crisprpam: NGG PAMs * crisprprimer: primer (on PAM strand) * crisprtranscript: reverse transcript (on PAM strand) * crisprextension: 3' extension of gRNA contains: reverse transcription template + primer binding site sequence can be found on non-PAM strand * crisprextrange: genomic range of crispr extension * Doench2016|4: on-target efficiency scores * off0, off1, off2: number of offtargets with 0, 1, 2 mismatches * off: total number of offtargets: off = off0 + off1 + ... * nickrange: nickspacer range * nickspacer: nickspacer sequence * nickDoench2016|4: nickspacer Doench scores * nickoff: nickspacer offtarget counts

See Also

find_spacers to find standard crispr sites

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
# Find PE spacers for 4 clinically relevant loci (Anzalone et al, 2019)
    bsgenome <- BSgenome.Hsapiens.UCSC.hg38::BSgenome.Hsapiens.UCSC.hg38  
    gr <- char_to_granges(c(
        PRNP = 'chr20:4699600:+',             # snp: prion disease
        HBB  = 'chr11:5227002:-',             # snp: sickle cell anemia
        HEXA = 'chr15:72346580-72346583:-',   # del: tay sachs disease
        CFTR = 'chr7:117559593-117559595:+'), # ins: cystic fibrosis
        bsgenome)
    spacers <- find_primespacers(gr, bsgenome)
    spacers <- find_spacers(extend_for_pe(gr), bsgenome, complement = FALSE)
    
# Edit PRNP locus for resistance against prion disease (Anzalone et al, 2019)
    bsgenome <- BSgenome.Hsapiens.UCSC.hg38::BSgenome.Hsapiens.UCSC.hg38  
    gr <- char_to_granges(c(PRNP = 'chr20:4699600:+'), bsgenome)
    find_primespacers(gr, bsgenome)
    find_primespacers(gr, bsgenome, edits = 'T')

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