Mercurial > repos > rnateam > chipseeker
comparison chipseeker.R @ 1:2019b4dd86a8 draft
planemo upload for repository https://github.com/galaxyproject/tools-iuc/tree/master/tools/chipseeker commit 3419a5a5e19a93369c8c20a39babe5636a309292
| author | rnateam |
|---|---|
| date | Tue, 29 May 2018 15:07:04 -0400 |
| parents | |
| children | 195cba35110e |
comparison
equal
deleted
inserted
replaced
| 0:313332cb7745 | 1:2019b4dd86a8 |
|---|---|
| 1 options( show.error.messages=F, error = function () { cat( geterrmessage(), file=stderr() ); q( "no", 1, F ) } ) | |
| 2 | |
| 3 # we need that to not crash galaxy with an UTF8 error on German LC settings. | |
| 4 loc <- Sys.setlocale("LC_MESSAGES", "en_US.UTF-8") | |
| 5 | |
| 6 suppressPackageStartupMessages({ | |
| 7 library(ChIPseeker) | |
| 8 library(GenomicFeatures) | |
| 9 library(optparse) | |
| 10 }) | |
| 11 | |
| 12 option_list <- list( | |
| 13 make_option(c("-i","--infile"), type="character", help="Peaks file to be annotated"), | |
| 14 make_option(c("-G","--gtf"), type="character", help="GTF to create TxDb."), | |
| 15 make_option(c("-u","--upstream"), type="integer", help="TSS upstream region"), | |
| 16 make_option(c("-d","--downstream"), type="integer", help="TSS downstream region"), | |
| 17 make_option(c("-F","--flankgeneinfo"), type="logical", help="Add flanking gene info"), | |
| 18 make_option(c("-D","--flankgenedist"), type="integer", help="Flanking gene distance"), | |
| 19 make_option(c("-f","--format"), type="character", help="Output format (interval or tabular)."), | |
| 20 make_option(c("-p","--plots"), type="character", help="PDF of plots.") | |
| 21 ) | |
| 22 | |
| 23 parser <- OptionParser(usage = "%prog [options] file", option_list=option_list) | |
| 24 args = parse_args(parser) | |
| 25 | |
| 26 peaks = args$infile | |
| 27 gtf = args$gtf | |
| 28 up = args$upstream | |
| 29 down = args$downstream | |
| 30 format = args$format | |
| 31 plots = args$plots | |
| 32 | |
| 33 peaks <- readPeakFile(peaks) | |
| 34 | |
| 35 # Make TxDb from GTF | |
| 36 txdb <- makeTxDbFromGFF(gtf, format="gtf") | |
| 37 if (!is.null(args$flankgeneinfo)) { | |
| 38 peakAnno <- annotatePeak(peaks, TxDb=txdb, tssRegion=c(-up, down), addFlankGeneInfo=args$flankgeneinfo, flankDistance=args$flankgenedist) | |
| 39 } else { | |
| 40 peakAnno <- annotatePeak(peaks, TxDb=txdb, tssRegion=c(-up, down)) | |
| 41 } | |
| 42 | |
| 43 # Convert from 1-based to 0-based format | |
| 44 res <- as.GRanges(peakAnno) | |
| 45 metacols <- mcols(res) | |
| 46 if (format == "interval") { | |
| 47 metacols <- apply(as.data.frame(metacols), 1, function(col) paste(col, collapse="|")) | |
| 48 resout <- data.frame(Chrom=seqnames(res), | |
| 49 Start=start(res) - 1, | |
| 50 End=end(res), | |
| 51 Comment=metacols) | |
| 52 } else { | |
| 53 resout <- data.frame(Chrom=seqnames(res), | |
| 54 Start=start(res) - 1, | |
| 55 End=end(res), | |
| 56 metacols) | |
| 57 } | |
| 58 | |
| 59 write.table(resout, file="out.tab", sep="\t", row.names=FALSE, quote=FALSE) | |
| 60 | |
| 61 if (!is.null(plots)) { | |
| 62 pdf("out.pdf", width=14) | |
| 63 plotAnnoPie(peakAnno) | |
| 64 plotAnnoBar(peakAnno) | |
| 65 vennpie(peakAnno) | |
| 66 upsetplot(peakAnno) | |
| 67 plotDistToTSS(peakAnno, title="Distribution of transcription factor-binding loci\nrelative to TSS") | |
| 68 dev.off() | |
| 69 } |
