count5ends: Prepare counts matrix for enrichment analysis

View source: R/count5ends.R

count5endsR Documentation

Prepare counts matrix for enrichment analysis

Description

Prepare the counts matrix by 5'end of reads.

Usage

count5ends(
  bam,
  index = bam,
  yieldSize = 1e+05,
  positive = 4L,
  negative = 5L,
  bindingSites,
  bindingSitesWithGap,
  bindingSitesWithProximal,
  bindingSitesWithProximalAndGap,
  bindingSitesWithDistal
)

Arguments

bam

A character vector indicates the file names of the bams or an object of BamFile.

index

The names of the index file of the 'BAM' file being processed; This is given without the '.bai' extension.

yieldSize

Number of records to yield each time the file is read. See BamFile for details.

positive, negative

integer(1). the size to be shift for positive/negative strand. If the bam file is 5'end shifed files, please set the parameter to 0.

bindingSites

A object of GenomicRanges indicates candidate binding sites. The prepareBindingSites function is a helper function to generate the binding sites. Users can also use other software for example fimo to generate the list.

bindingSitesWithGap

bindingSites with gaps and in both ends,

bindingSitesWithProximal

bindingSites with gaps and proximal region in both ends,

bindingSitesWithProximalAndGap

bindingSites with gaps, and then proximal and gaps in both ends,

bindingSitesWithDistal

bindingSites with gap, proximal, gap and distal regions.

Value

A RangedSummarizedExperiment object with assays of count matrix with bindingSites, proximalRegion and distalRegion as column names and bindingSites GRanges object as rowRanges.

Author(s)

Jianhong Ou

Examples

bam <- system.file("extdata",
                   "KD.shift.rep1.bam",
                   package="ATACseqTFEA")
bsl <- system.file("extdata", "bindingSites.rds",
                   package="ATACseqTFEA")
bindingSites <- readRDS(bsl)
## get the count regions
bsEx <- expandBindingSites(bindingSites)
res <- count5ends(bam, positive=0L, negative=0L,
                  bindingSites=bindingSites,
                  bindingSitesWithGap=bsEx$bindingSitesWithGap,
                  bindingSitesWithProximal=bsEx$bindingSitesWithProximal,
                  bindingSitesWithProximalAndGap=
                      bsEx$bindingSitesWithProximalAndGap,
                  bindingSitesWithDistal=bsEx$bindingSitesWithDistal)
head(res)

jianhong/ATACseqTFEA documentation built on Nov. 5, 2024, 9:46 a.m.