analyze_fragmentation: Provides fragment ends analysis

Description Usage Arguments Details Value See Also Examples

View source: R/analyze_fragmentation.R

Description

Calculates the number of fragment ends and the Windowed Protection Score (WPS) in genomic tiles within targets

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
analyze_fragmentation(
  bam,
  targets,
  tag = "",
  window_size = 120,
  step_size = 5,
  min_size = 120,
  max_size = 180,
  ...
)

Arguments

bam

the input bam file

targets

The targets to restrict the windows within. Must have the columns chr, start and end. In case of whole-genome, specify full chromosomes targets.

tag

the RG tag if the bam has more than one sample.

window_size

The window (bin) size to use within the targets

step_size

The step size to use in case of overlapping bins.

min_size

Restrict fragments to this minimum size.

max_size

Restrict fragments to this maximum size.

...

Other parameters passed to get_fragment_size

Details

Fragment length will extracted from the bam file according to the parameters passed to get_fragment_size, and the number of fragment ends, and the Windowed Protection Score (WPS) will be computed in the binned input targets. Binning is done according to the window_size and step_size parameters.

WPS is defined as the number of fragments completely spanning a window (bin) minus the number of fragments with an endpoint within the same window as reported by Snyder et al., Cell 2016.

The output include both the fragment end counts and the WPS in their raw format as well as after adjustment by coverage in the bin.

Minimum and maximum bounds of the fragment size are applied before computing WPS and fragment ends counts.

Value

a data frame with the first three columns having the bins coordinates and other columns having the WPS (raw and adjusted by coverage) and number of fragment ends (raw and adjusted by coverage).

See Also

get_fragment_size bin_fragment_size summarize_fragment_size

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
data("targets", package = "ctDNAtools")
bamN1 <- system.file("extdata", "N1.bam", package = "ctDNAtools")

## basic usage
analyze_fragmentation(bam = bamN1, targets = targets)

## more options
analyze_fragmentation(
  bam = bamN1, targets = targets,
  step_size = 10, window_size = 50
)

alkodsi/ctDNAtools documentation built on Feb. 22, 2022, 9:40 a.m.