getPeriodicityTrack | R Documentation |
This function takes a set of GRanges in a genome, recover the corresponding sequences and divides them using a sliding window. For each sub-sequence, it then computes the PSD value of a k-mer of interest at a chosen period, and generates a linear .bigWig track from these values.
getPeriodicityTrack( genome = NULL, granges, motif = "WW", period = 10, BPPARAM = setUpBPPARAM(1), extension = 1000, window_size = 100, step_size = 2, range_spectrum = seq(5, 50), smooth_track = 20, bw_file = NULL )
genome |
DNAStringSet, BSgenome or genome ID |
granges |
GRanges object |
motif |
character, k-mer of interest. |
period |
Integer, the period of the k-mer to study (default=10). |
BPPARAM |
split the workload over several processors using BiocParallel |
extension |
Integer, the width the GRanges are going to be extended to (default 1000). |
window_size |
Integer, the width of the bins to split the GRanges objects in (default 100). |
step_size |
Integer, the increment between bins over GRanges (default 2). |
range_spectrum |
Numeric vector, the distances between nucleotides to take into consideration when performing Fast Fourier Transform (default seq_len(50)). |
smooth_track |
Integer, smooth the resulting track |
bw_file |
character, the name of the output bigWig track |
Rlelist and a bigWig track in the working directory.
data(ce11_proms) track <- getPeriodicityTrack( genome = 'BSgenome.Celegans.UCSC.ce11', ce11_proms[1], extension = 200, window_size = 100, step_size = 10, smooth_track = 1, motif = 'WW', period = 10, BPPARAM = setUpBPPARAM(1) ) track unlink( 'BSgenome.Celegans.UCSC.ce11_WW_10-bp-periodicity_g-100^10_smooth-1.bw' )
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.