This function computes PSD values of a given k-mer of interest in a set of input sequences. It also iterates the PSD calculation process over shuffled sequences, if n_shuffling is used.

getPeriodicityWithIterations(x, ...)

# S3 method for DNAStringSet
getPeriodicityWithIterations(
  x,
  motif,
  n_shuffling = 10,
  cores_shuffling = 1,
  cores_computing = 1,
  order = 1,
  verbose = 1,
  ...
)

# S3 method for GRanges
getPeriodicityWithIterations(x, genome, ...)

Arguments

x

DNAStringSet, sequences of interest

...

Arguments passed to S3 methods

motif

character, k-mer of interest

n_shuffling

integer, Number of shuffling

cores_shuffling

integer, Number of cores used for shuffling

cores_computing

integer, split the workload over several processors using BiocParallel

order

Integer, which order to take into consideration for shuffling (ushuffle python library must be installed for orders > 1)

verbose

integer, Should the function be verbose?

genome

genome ID, BSgenome or DNAStringSet object (optional, if x is a GRanges)

Value

Several metrics

Methods (by class)

  • DNAStringSet: S3 method for DNAString

  • GRanges: S3 method for GRanges

Examples

data(ce11_proms_seqs) res <- getPeriodicityWithIterations( ce11_proms_seqs[1:10], genome = 'BSgenome.Celegans.UCSC.ce11', motif = 'TT', cores_shuffling = 1 )
#> - Calculating observed PSD
#> - Mapping k-mers.
#> - 12945 pairwise distances measured.
#> - Calculating pairwise distance distribution.
#> - Normalizing distogram vector.
#> - Applying Fast Fourier Transform to the normalized distogram.
#> - Shuffling 1/10
#> - Shuffling 2/10
#> - Shuffling 3/10
#> - Shuffling 4/10
#> - Shuffling 5/10
#> - Shuffling 6/10
#> - Shuffling 7/10
#> - Shuffling 8/10
#> - Shuffling 9/10
#> - Shuffling 10/10
#> Only 10 shufflings. Cannot compute accurate empirical p-values. To compute empirical p-values, set up n_shuffling to at least 100. Only l2FC values are returned
res$observed_PSD
#> NULL
res$shuffled_PSD
#> NULL