runBreakpointr: Find breakpoints in Strand-seq data

View source: R/runBreakpointr.R

runBreakpointrR Documentation

Find breakpoints in Strand-seq data

Description

Find breakpoints in Strand-seq data. See section Details on how breakpoints are located.

Usage

runBreakpointr(
  bamfile,
  ID = basename(bamfile),
  pairedEndReads = TRUE,
  chromosomes = NULL,
  windowsize = 1e+06,
  binMethod = "size",
  multi.sizes = NULL,
  trim = 10,
  peakTh = 0.33,
  zlim = 3.291,
  background = 0.05,
  min.mapq = 10,
  pair2frgm = FALSE,
  filtAlt = FALSE,
  genoT = "fisher",
  minReads = 20,
  maskRegions = NULL,
  conf = 0.99
)

Arguments

bamfile

A file with aligned reads in BAM format.

ID

A character string that will serve as identifier in downstream functions.

pairedEndReads

Set to TRUE if you have paired-end reads in your file.

chromosomes

If only a subset of the chromosomes should be binned, specify them here.

windowsize

The window size used to calculate deltaWs, either number of reads or genomic size depending on binMethod.

binMethod

Method used to calculate optimal number of reads in the window ("size", "reads"). By default binMethod='size'.

multi.sizes

User defined multiplications of the original window size.

trim

The amount of outliers in deltaWs removed to calculate the stdev (10 will remove top 10% and bottom 10% of deltaWs).

peakTh

The treshold that the peak deltaWs must pass to be considered a breakpoint (e.g. 0.33 is 1/3 of max(deltaW)).

zlim

The number of stdev that the deltaW must pass the peakTh (ensures only significantly higher peaks are considered).

background

The percent (e.g. 0.05 = 5%) of background reads allowed for WW or CC genotype calls.

min.mapq

Minimum mapping quality when importing from BAM files.

pair2frgm

Set to TRUE if every paired-end read should be merged into a single fragment.

filtAlt

Set to TRUE if you want to filter out alternative alignments defined in 'XA' tag.

genoT

A method ('fisher' or 'binom') to genotype regions defined by a set of breakpoints.

minReads

The minimal number of reads between two breaks required for genotyping.

maskRegions

List of regions to be excluded from the analysis in GRanges-class object.

conf

Desired confidence interval of localized breakpoints.

Details

Breakpoints are located in the following way:

  1. calculate deltaWs chromosome-by-chromsome

  2. localize breaks that pass zlim above the threshold

  3. genotype both sides of breaks to confirm whether strand state changes

  4. write a file of _reads, _deltaWs and _breaks in a chr fold -> can upload on to UCSC Genome browser

  5. write a file for each index with all chromosomes included -> can upload on to UCSC Genome browser

Value

A BreakPoint object.

Author(s)

David Porubsky, Ashley Sanders, Aaron Taudt

Examples

## Get an example file
exampleFolder <- system.file("extdata", "example_bams", package="breakpointRdata")
exampleFile <- list.files(exampleFolder, full.names=TRUE)[1]
## Run breakpointR
brkpts <- runBreakpointr(exampleFile, chromosomes='chr22', pairedEndReads=FALSE)


daewoooo/breakpointR documentation built on May 9, 2024, 5:03 p.m.