The cleanUpdTSeq user’s guide

Introduction

3’ ends of transcripts have generally been poorly annotated. With the advent of deep sequencing, many methods have been developed to identify 3’ ends. The majority of these methods use an oligo-dT primer, which can bind to internal adenine-rich sequences, and lead to artifactual identification of polyadenylation sites. Heuristic filtering methods rely on a certain number of adenines in the genomic sequence downstream of a putative polyadenylation site to remove internal priming events. We introduce a package to provide a robust method to classify putative polyadenylation sites. cleanUpdTSeq uses a naïve Bayes classifier, implemented through the e1071 [1], and sequence features surrounding the putative polyadenylation sites for classification.

The package includes a training dataset constructed from 6 different Zebrafish sequencing dataset, and functions for fetching surrounding sequences using BSgenome [2], building feature vectors and classifying whether the putative polyadenylations site is a true polyadenylation site or a mis-primed false site.

Citation

If you use cleanUpdTSeq, please cite:

Sheppard, S., Lawson, N.D. and Zhu, L.J., 2013. Accurate identification of polyadenylation sites from 3’ end deep sequencing using a naive Bayes classifier. Bioinformatics, 29(20), pp.2564-2571 and 2014, 30(4), pp.596, https://doi.org/10.1093/bioinformatics/btt714

step-by-step guide

Here is a step-by-step guide on how to use the cleanUpdTSeq package to classify a list of putative polyadenylation sites

Step 1. Load the cleanUpdTSeq package, and then use the function BED6WithSeq2GRangesSeq to convert the test dataset in a extended BED6 file to a GRanges object with or without sequence information.

suppressPackageStartupMessages(library(cleanUpdTSeq))
testFile <- system.file("extdata", "test.bed", 
                        package = "cleanUpdTSeq")
peaks <- BED6WithSeq2GRangesSeq(file = testFile, 
                               skip = 1L, withSeq = FALSE)

If test dataset contains sequence information already, then use the following command instead.

peaks <- BED6WithSeq2GRangesSeq(file = testFile, 
                                skip = 1L, withSeq = TRUE)

To work with your own test dataset, please set testFile to the file path that contains the putative pA site information.

Here is how the test dataset look like.

head(read.delim(testFile, header = TRUE, skip = 0))
##     chr    start     stop        name score strand
## 1 chr10  2965327  2965328 6hpas-22249     1      -
## 2 chr10  2966558  2966559 6hpas-22250     1      -
## 3 chr10  2974251  2974252 6hpas-22251     2      -
## 4 chr10  2978441  2978442 6hpas-22252     1      -
## 5 chr11 16772291 16772292 6hpas-33204     1      -
## 6 chr11 16777848 16777849 6hpas-33205     1      -
##                                   upstream                     downstream
## 1 TCTTCATCATGGTCATCTCGCACCAGAGAGTGTGCCAGGG CAGGAAGTTTTACCTGTCTGTCATTATCGT
## 2 ACCCTGGTGAGGGTATAGAGCTGGTCCAGTGTGCCACGGC AAAGAGGAAAACAGCATTGTTCCTCCTGGA
## 3 TGATTTGTTTGTAACTGATTTTATCTTTTAATAAAAAAGA AAAAAGAAAGTCAAGCCAAGAGGCAAATAC
## 4 GGAGCGCGACCGCATCAACAAAATCTTGCAGGATTATCAG AAGAAAAAGATGGTGAGTTATTATCATTCA
## 5 AGGGAAATAAATACAAAAGAATAAAAATATGATTCATTGT AAGAAAAACACTTTAGCTACAAAAGTCCTT
## 6 ATTTAGTTGGGTATTATTTCAAATAAAGAGAGAGAGAGAC ACAAAAACTACATCAAATTTGAGGACAAAA

Step2. Build feature vectors for the classifier using the function buildFeatureVector.

The zebrafish genome from BSgenome is used in this example for obtaining flanking sequences. For a list of other genomes available through BSgenome, please refer to the BSgenome package documentation [2].

suppressPackageStartupMessages(library(BSgenome.Drerio.UCSC.danRer7))
testSet.NaiveBayes <- buildFeatureVector(peaks, 
                                         genome = Drerio,
                                         upstream = 40L, 
                                         downstream = 30L, 
                                         wordSize = 6L,
                                         alphabet = c("ACGT"),
                                         sampleType = "unknown", 
                                         replaceNAdistance = 30, 
                                         method = "NaiveBayes",
                                         fetchSeq = TRUE,
                                         return_sequences = TRUE)

If sequences are present in the test dataset already, then set fetchSeq = FALSE.

Step 3. Load the training dataset and classify putative polyadenylation sites.

The output file is a tab-delimited file containing the name of the putative polyadenylation sites, the probability that the putative polyadenylation site is false/oligodT internally primed, the probability the putative polyadenylation site is true, the predicted class based on the assignment cutoff and the sequence surrounding the putative polyadenylation site.

data(data.NaiveBayes)
if(interactive()){
   out <- predictTestSet(data.NaiveBayes$Negative,
                         data.NaiveBayes$Positive, 
                         testSet.NaiveBayes = testSet.NaiveBayes, 
                         outputFile = file.path(tempdir(), 
                                          "test-predNaiveBayes.tsv"), 
                        assignmentCutoff = 0.5)
}

Alternatively, instead of passing in a positive and a negative training dataset, set the parameter classifier to a pre-built PASclassifier to speed up the process. To built a PASclassifier using the training dataset, please use function buildClassifier. A PASclassifier named as classifier is included in the package which is generated using the included training dataset with upstream = 40, downstream = 30, and wordSize = 6. Please note that in order to use this pre-built classier, you need to build feature vector using buildFeatureVector from your test dataset with the same setting, i.e., upstream = 40, downstream = 30, and wordSize = 6.

data(classifier)
testResults <- predictTestSet(testSet.NaiveBayes = testSet.NaiveBayes,
                              classifier = classifier,
                              outputFile = NULL, 
                              assignmentCutoff = 0.5,
                              return_sequences = TRUE)
head(testResults)
##     peak_name prob_fake_pA prob_true_pA pred_class
## 1 6hpas-22249   0.99778351 2.216486e-03          0
## 2 6hpas-22250   1.00000000 6.699365e-10          0
## 3 6hpas-22251   0.99444054 5.559459e-03          0
## 4 6hpas-22252   0.99999995 4.729351e-08          0
## 5 6hpas-33204   0.02810869 9.718913e-01          1
## 6 6hpas-33205   0.99997144 2.855673e-05          0
##                               upstream_seq                 downstream_seq
## 1 TCTTCATCATGGTCATCTCGCACCAGAGAGTGTGCCAGGG CAGGAAGTTTTACCTGTCTGTCATTATCGT
## 2 ACCCTGGTGAGGGTATAGAGCTGGTCCAGTGTGCCACGGC AAAGAGGAAAACAGCATTGTTCCTCCTGGA
## 3 TGATTTGTTTGTAACTGATTTTATCTTTTAATAAAAAAGA AAAAAGAAAGTCAAGCCAAGAGGCAAATAC
## 4 GGAGCGCGACCGCATCAACAAAATCTTGCAGGATTATCAG AAGAAAAAGATGGTGAGTTATTATCATTCA
## 5 AGGGAAATAAATACAAAAGAATAAAAATATGATTCATTGT AAGAAAAACACTTTAGCTACAAAAGTCCTT
## 6 ATTTAGTTGGGTATTATTTCAAATAAAGAGAGAGAGAGAC ACAAAAACTACATCAAATTTGAGGACAAAA

References

  1. Meyer, D., et al., e1071: Misc Functions of the Department of Statistics (e1071), TU Wien. 2012.

  2. Pages, H., BSgenome: Infrastructure for Biostrings-based genome data packages.

  3. Sheppard, S., Lawson, N.D. and Zhu, L.J., 2013. Accurate identification of polyadenylation sites from 3’ end deep sequencing using a naive Bayes classifier. Bioinformatics, 29(20), pp.2564-2571.

Session Info

sessionInfo()

R version 4.4.2 (2024-10-31) Platform: x86_64-pc-linux-gnu Running under: Ubuntu 24.04.1 LTS

Matrix products: default BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0

locale: [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
[3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
[5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
[7] LC_PAPER=en_US.UTF-8 LC_NAME=C
[9] LC_ADDRESS=C LC_TELEPHONE=C
[11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C

time zone: Etc/UTC tzcode source: system (glibc)

attached base packages: [1] stats4 stats graphics grDevices utils datasets methods
[8] base

other attached packages: [1] cleanUpdTSeq_1.45.0 BSgenome.Drerio.UCSC.danRer7_1.4.0 [3] BSgenome_1.75.0 rtracklayer_1.67.0
[5] BiocIO_1.17.1 Biostrings_2.75.1
[7] XVector_0.47.0 GenomicRanges_1.59.1
[9] GenomeInfoDb_1.43.2 IRanges_2.41.1
[11] S4Vectors_0.45.2 BiocGenerics_0.53.3
[13] generics_0.1.3 BiocStyle_2.35.0

loaded via a namespace (and not attached): [1] SummarizedExperiment_1.37.0 rjson_0.2.23
[3] xfun_0.49 bslib_0.8.0
[5] Biobase_2.67.0 lattice_0.22-6
[7] vctrs_0.6.5 tools_4.4.2
[9] bitops_1.0-9 curl_6.0.1
[11] parallel_4.4.2 proxy_0.4-27
[13] Matrix_1.7-1 lifecycle_1.0.4
[15] GenomeInfoDbData_1.2.13 stringr_1.5.1
[17] compiler_4.4.2 Rsamtools_2.23.1
[19] codetools_0.2-20 htmltools_0.5.8.1
[21] sys_3.4.3 buildtools_1.0.0
[23] class_7.3-22 sass_0.4.9
[25] RCurl_1.98-1.16 yaml_2.3.10
[27] crayon_1.5.3 jquerylib_0.1.4
[29] seqinr_4.2-36 MASS_7.3-61
[31] BiocParallel_1.41.0 cachem_1.1.0
[33] DelayedArray_0.33.2 abind_1.4-8
[35] digest_0.6.37 stringi_1.8.4
[37] restfulr_0.0.15 maketools_1.3.1
[39] ade4_1.7-22 fastmap_1.2.0
[41] grid_4.4.2 cli_3.6.3
[43] SparseArray_1.7.2 magrittr_2.0.3
[45] S4Arrays_1.7.1 XML_3.99-0.17
[47] e1071_1.7-16 UCSC.utils_1.3.0
[49] rmarkdown_2.29 httr_1.4.7
[51] matrixStats_1.4.1 evaluate_1.0.1
[53] knitr_1.49 rlang_1.1.4
[55] Rcpp_1.0.13-1 glue_1.8.0
[57] BiocManager_1.30.25 jsonlite_1.8.9
[59] R6_2.5.1 MatrixGenerics_1.19.0
[61] GenomicAlignments_1.43.0 zlibbioc_1.52.0