generateControl: generateControl

generateControlR Documentation

generateControl

Description

function to generate random control region for given regions.

Usage

generateControl(regions, ...)

## S4 method for signature 'GenomicRanges'
generateControl(
  regions,
  genome,
  gene.txdb = NULL,
  n.random = 5,
  random.seed = NULL
)

Arguments

regions

GenomicRanges for generating control regions

genome

BSgenome object, or short string signifying genome build recognized by getBSgenome.

gene.txdb

TxDb gene annotation for control regions, See Details.(default=NULL)

n.random

number for generating random control.For each region, generating n.random control regions.(default=5)

random.seed

seed for random programs

Details

If the gene.txdb is specified, control regions will be generated from the same TSS distance to random gene of given region. Otherwise, the control region will be generated by random location with the same chromosome and region width of given regions.

Value

GenomicRanges random control regions with n.random times the length of the input regions

Methods (by class)

  • generateControl(GenomicRanges): GenomicRanges

Examples

library(BSgenome.Hsapiens.UCSC.hg19)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)

example_motifs <- getJasparMotifs(species = "Homo sapiens",
                                  collection = "CORE")
# Make a set of input regions
Input <- GenomicRanges::GRanges(seqnames = c("chr1","chr2","chr2"),
                                ranges = IRanges::IRanges(start = c(76585873,42772928,
                                                                    100183786),
                                                          width = 500))
# Make a set of control regions by generateControl
Control <- generateControl(Input, genome=BSgenome.Hsapiens.UCSC.hg19,
                           gene.txdb=TxDb.Hsapiens.UCSC.hg19.knownGene)

# Scan motif for example motifs
motif_ix_input <- motifScan(example_motifs, Input, genome = "BSgenome.Hsapiens.UCSC.hg19")
motif_ix_control <- motifScan(example_motifs, Control, genome = "BSgenome.Hsapiens.UCSC.hg19")

# Find Enrichment motif of input by control
Enrichment_result <- motifEnrichment(motif_ix_input, motif_ix_control)


LouisKwok-PICB/motifscanr documentation built on July 22, 2024, 6:36 a.m.