@@ -440,6 +440,8 @@ msScope_annotation <- function(libs) {
440
440
# ' annotateCTSS(exampleCAGEexp, exampleZv9_annot)
441
441
# ' colData(exampleCAGEexp)
442
442
# '
443
+ # ' @importFrom GenomicFeatures genes
444
+ # '
443
445
# ' @export
444
446
445
447
setGeneric ("annotateCTSS ", function(object, ranges, upstream=500, downstream=500)
@@ -459,6 +461,22 @@ setMethod("annotateCTSS", c("CAGEexp", "GRanges"), function (object, ranges, ups
459
461
object
460
462
})
461
463
464
+ # ' @rdname annotateCTSS
465
+
466
+ setMethod ("annotateCTSS ", c("CAGEexp", "TxDb"), function (object, ranges){
467
+ g <- genes(ranges )
468
+ g $ gene_name <- g $ gene_id
469
+ annotateCTSS(object , g )
470
+ CTSScoordinatesGR(object )$ genes <- ranges2genes(CTSScoordinatesGR(object ), g )
471
+ CTSScoordinatesGR(object )$ annotation <- txdb2annot(CTSScoordinatesGR(object ), ranges )
472
+
473
+ annot <- sapply( CTSStagCountDF(object )
474
+ , function (X ) tapply(X , CTSScoordinatesGR(object )$ annotation , sum ))
475
+ colData(object )[levels(CTSScoordinatesGR(object )$ annotation )] <- DataFrame(t(annot ))
476
+ validObject(object )
477
+ object
478
+ })
479
+
462
480
# ' `annotateTagClusters` annotates the _tag clusters_ of a _CAGEr_ object.
463
481
# '
464
482
# ' @return `annotateTagClusters` returns the input object with the same
@@ -619,6 +637,33 @@ ranges2annot <- function(ranges, annot, upstream=500, downstream=500) {
619
637
Rle(annot )
620
638
}
621
639
640
+ # ' @name txdb2annot
641
+ # '
642
+ # ' @rdname ranges2annot
643
+ # '
644
+ # ' @importFrom GenomicFeatures promoters exons transcripts
645
+
646
+ txdb2annot <- function (ranges , annot ) {
647
+ findOverlapsBool <- function (A , B ) {
648
+ overlap <- findOverlaps(A , B )
649
+ overlap <- as(overlap , " List" )
650
+ any(overlap )
651
+ }
652
+
653
+ classes <- c(" promoter" , " exon" , " intron" , " unknown" )
654
+ p <- findOverlapsBool(ranges , trim(suppressWarnings(promoters(annot , 500 , 500 ))))
655
+ e <- findOverlapsBool(ranges , exons(annot ))
656
+ t <- findOverlapsBool(ranges , transcripts(annot ))
657
+ annot <- sapply( 1 : length(ranges ), function (i ) {
658
+ if (p [i ]) {classes [1 ]}
659
+ else if (e [i ]) {classes [2 ]}
660
+ else if (t [i ]) {classes [3 ]}
661
+ else {classes [4 ]}
662
+ })
663
+
664
+ annot <- factor (annot , levels = classes )
665
+ Rle(annot )
666
+ }
622
667
623
668
# ' ranges2genes
624
669
# '
0 commit comments