View source: R/SummitHeatmap.R
SummitHeatmap | R Documentation |
Generate heatmaps of read counts around positions of interest in a genome.
SummitHeatmap(
peaks,
bamFiles,
bamNames = "myreads",
span = 2025,
step = 50,
minOverlap = 1,
useCPM = TRUE,
PairedEnd = FALSE,
minMQS = 255,
strand = 0,
splitOnly = FALSE,
nonSplitOnly = FALSE,
readExtension3 = 0,
readShiftSize = 0,
requireBothEndsMapped = FALSE,
read2pos = 5,
mode = "F",
genome = "BSgenome.Mmusculus.UCSC.mm10"
)
peaks |
A GRanges object containing your regions of interest. Must include seqnames (chromosomes), start, end, strand, and name. |
bamFiles |
Character vector containing the filenames (including the full path) of read alignment files in bam format. |
bamNames |
Character vector containing the names to describe the |
span |
Integer scalar specifyig the distance from the peak center to the left and right that you want your heatmap to extend to. Default is 2025. |
step |
Integer scalar specifyig the window size in which reads are counted. Default is 50. |
minOverlap |
Integer scalar giving the minimum number of overlapping bases required for assigning a read to a heatmap window. For assignment of read pairs, number of overlapping bases from each read in the same pair will be summed. If a negative value is provided, then a gap of up to specified size will be allowed between read and the feature that the read is assigned to. 1 by default. |
useCPM |
Logical scalar indicationg wether to normalize the number of reads per window to the total number of reads in the sample (bam file). Default is TRUE. |
PairedEnd |
Logical scalar, indicating wether reads in bam files were generated with paired-end or single-end sequencing. Default is FALSE (=single-end). |
minMQS |
Integer scalar, specifying the minimum mapping quality that a read must have to be included. Default is 255, which eliminates multimapping
reads in case the STAR aligner was used to generate the |
strand |
Integer vector indicating if strand-specific read counting should be performed. Length of the vector should be either 1 (meaning that the value is applied to all input files), or equal to the total number of input files provided. Each vector element should have one of the following three values: 0 (unstranded), 1 (stranded) and 2 (reversely stranded). Default value of this parameter is 0 (ie. unstranded read counting is performed for all input files). |
splitOnly |
Logical scalar indicating whether only split alignments (their CIGAR strings contain letter 'N') should be included. FALSE by default. |
nonSplitOnly |
Logical scalar indicating whether only non-split alignments (their CIGAR strings do not contain letter 'N') should be included. FALSE by default. |
readExtension3 |
Integer scalar giving the number of bases extended downstream from 3' end of each read. 0 by default. Negative value is not allowed. |
readShiftSize |
Integer scalar specifying the number of bases the reads will be shifted downstream by. 0 by default. Negative value is not allowed. In case of mode = "Q" and PairedEnd = TRUE, it should be set to 'halfInsert' to shift each read by half the fragment size. |
requireBothEndsMapped |
Logical scalar indicating if both ends from the same fragment are required to be successfully aligned before the fragment can be assigned to a feature or meta-feature. This parameter is only appliable when PairedEnd is TRUE. |
read2pos |
Specifying whether each read should be reduced to its 5' most base or 3' most base. It has three possible values: NULL, 5 (denoting 5' most base) and 3 (denoting 3' most base). Default value is 5, ie. only the 5' end of the read will be counted. If a read is reduced to a single base, only that base will be considered for the read assignment. Read reduction is performed after read shifting and extension. |
mode |
Specify if you want to use QuasR (mode = "Q") or FeatureCounts (mode = "F", default) to count the number of rads per window. |
genome |
The reference genome, either a string specifying a BSgenome object or the name of a genome fasta file. Default is "BSgenome.Mmusculus.UCSC.mm10". |
This function generates heatmaps of different experiments (read counts) around positions of interest in the genome, for example the summits of ChIP peaks, or transcription start sites. The center of your GRanges object defines the coordinates of these positions of interest (eg peak summits, TSSs). Features that lie on the minus strand will be reversed in the final output.
A list of matrices of the same length as the number of samples supplied in bamFiles
.
Each matrix contains the number of reads (or cpm) in each window (column headers indicate the middle of the window)
around the center of each region provided in peaks
.
peaks <- SimulatePeaks(1000,rep(100,1000),chromosomeSizes=
system.file("extdata", "chrNameLength_mm10_chr11.txt", package = "MiniChip"))
#peaks <- GenomicRanges::GRanges(
#seqnames = Rle(c("chr1", "chr2", "chr1", "chr3"), c(1, 3, 2, 4)),
#ranges = IRanges(50101:50110, end = 51111:51120),
#strand = Rle(strand(c("-", "+", "*", "+", "-")), c(1, 2, 2, 3, 2)),
#summit = 1:10, name = head(letters, 10))
bamFiles <- list.files(system.file("extdata", package = "MiniChip"),
full.names=TRUE,pattern="*bam$")
SummitHeatmap(peaks=peaks,bamFiles=bamFiles)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.