Last data update: 2014.03.03

R: Automatic merging of overlapped nucleosome calls
mergeCallsR Documentation

Automatic merging of overlapped nucleosome calls

Description

This function joints close nucleosome calls into one larger, fuzzy nucleosome.

Usage

mergeCalls(calls, min.overlap=50, discard.low=0.2, mc.cores=1, verbose=TRUE)

Arguments

calls

RangedData with scored and ranged nucleosome calls from peakScoring or peakDetection(..., score=TRUE).

min.overlap

Minimum overlap between two reads for merge them

discard.low

Discard low covered calls (i.e. calls with score_h < discard.low will be discarded)

mc.cores

Number of cores available to parallel data processing.

verbose

Show progress info?

Details

This functions looks for overlapped calls and join those with more than min.overlap bases overlapped. More than two reads can be joined in one single call if all of them are overlapped at least that distance with almost another read in the range.

Joining is performed in chain, so if nucleosome call A is close to B and B is close to C, the final call will comprise the range A-B-C. The resulting scores (mixed, width, height) of the final joined call will be the average value of the individual scores.

The parameter discard.low allows to ignore the small peaks that could be merged with larger ones, originating large calls. In the case that all of the overlapped reads in a given position have score_h less than discard.low, all of them will be selected instead of deleting that call.

Value

RangedData with merged calls and the additional data column nmerge, with the count of how many original ranges are merged in the resulting range.

Author(s)

Oscar Flores oflores@mmb.pcb.ub.es

See Also

peakScoring

Examples

	#Generate a synthetic coverage map (assuming reads of 40bp and fragments of 130)
	map = syntheticNucMap(wp.num=20, fuz.num=20,  nuc.len=40, lin.len=130, rnd.seed=1)
	cover = filterFFT(coverage(map$syn.reads))

	#Find peaks over FFT filtered coverage
	calls = peakDetection(filterFFT(cover, pcKeepComp=0.02), width=130, score=TRUE)

	#Merge overlapped calls
	merged_calls = mergeCalls(calls)

	plotPeaks(merged_calls, cover)

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(nucleR)
Loading required package: ShortRead
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Loading required package: BiocParallel
Loading required package: Biostrings
Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: XVector
Loading required package: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: GenomicAlignments
Loading required package: SummarizedExperiment
Loading required package: Biobase
Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/nucleR/mergeCalls.Rd_%03d_medium.png", width=480, height=480)
> ### Name: mergeCalls
> ### Title: Automatic merging of overlapped nucleosome calls
> ### Aliases: mergeCalls
> ### Keywords: manip
> 
> ### ** Examples
> 
> 	#Generate a synthetic coverage map (assuming reads of 40bp and fragments of 130)
> 	map = syntheticNucMap(wp.num=20, fuz.num=20,  nuc.len=40, lin.len=130, rnd.seed=1)
> 	cover = filterFFT(coverage(map$syn.reads))
> 
> 	#Find peaks over FFT filtered coverage
> 	calls = peakDetection(filterFFT(cover, pcKeepComp=0.02), width=130, score=TRUE)
> 
> 	#Merge overlapped calls
> 	merged_calls = mergeCalls(calls)
* Starting space: 1
 - Finding overlapped reads
 - Constructing merge list
 - Merging calls
 - Formatting results
 - Done (6 non-overlapped | 6 merged calls)
> 
> 	plotPeaks(merged_calls, cover)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>