RNAmodR 1.21.0
For users interested in the general aspect of any RNAmodR
based package please
have a look at the main vignette of the package.
This vignette is aimed at developers and researchers, who want to use the
functionality of the RNAmodR
package to develop a new modification strategy
based on high throughput sequencing data.
library(RNAmodR)
Two classes have to be considered to establish a new analysis pipeline using
RNAmodR
. These are the SequenceData
and the Modifier
class.
SequenceData
classFirst, the SequenceData
class has to be considered. Several classes are
already implemented, which are:
End5SequenceData
End3SequenceData
EndSequenceData
ProtectedEndSequenceData
CoverageSequenceData
PileupSequenceData
NormEnd5SequenceData
NormEnd3SequenceData
If these cannot be reused, a new class can be implemented quite easily. First
the DataFrame class, the Data class and a constructor has to defined. The only
value, which has to be provided, is a default minQuality
integer value and
some basic information.
setClass(Class = "ExampleSequenceDataFrame",
contains = "SequenceDFrame")
ExampleSequenceDataFrame <- function(df, ranges, sequence, replicate,
condition, bamfiles, seqinfo){
RNAmodR:::.SequenceDataFrame("Example",df, ranges, sequence, replicate,
condition, bamfiles, seqinfo)
}
setClass(Class = "ExampleSequenceData",
contains = "SequenceData",
slots = c(unlistData = "ExampleSequenceDataFrame"),
prototype = list(unlistData = ExampleSequenceDataFrame(),
unlistType = "ExampleSequenceDataFrame",
minQuality = 5L,
dataDescription = "Example data"))
ExampleSequenceData <- function(bamfiles, annotation, sequences, seqinfo, ...){
RNAmodR:::SequenceData("Example", bamfiles = bamfiles,
annotation = annotation, sequences = sequences,
seqinfo = seqinfo, ...)
}
Second, the getData
function has to be implemented. This is used to load
the data from a bam file and must return a named list IntegerList
,
NumericList
or CompressedSplitDataFrameList
per file.
setMethod("getData",
signature = c(x = "ExampleSequenceData",
bamfiles = "BamFileList",
grl = "GRangesList",
sequences = "XStringSet",
param = "ScanBamParam"),
definition = function(x, bamfiles, grl, sequences, param, args){
###
}
)
Third, the aggregate
function has to be implemented. This function is used to
aggregate data over replicates for all or one of the conditions. The resulting
data is passed on to the Modifier
class.
setMethod("aggregateData",
signature = c(x = "ExampleSequenceData"),
function(x, condition = c("Both","Treated","Control")){
###
}
)
Modifier
classA new Modifier
class is probably the main class, which needs to be
implemented. Three variable have to be set. mod
must be a single element from
the Modstrings::shortName(Modstrings::ModRNAString())
. score
is the default
score, which is used for several function. A column with this name should be
returned from the aggregate
function. dataType
defines the SequenceData
class to be used. dataType
can contain multiple names of a SequenceData
class, which are then combined to form a SequenceDataSet
.
setClass("ModExample",
contains = c("RNAModifier"),
prototype = list(mod = "X",
score = "score",
dataType = "ExampleSequenceData"))
ModExample <- function(x, annotation, sequences, seqinfo, ...){
RNAmodR:::Modifier("ModExample", x = x, annotation = annotation,
sequences = sequences, seqinfo = seqinfo, ...)
}
dataType
can also be a list
of character
vectors, which leads then to the
creation of SequenceDataList
. However, for now this is a hypothetical case and
should only be used, if the detection of a modification requires bam files from
two or more different methods to be used to detect one modification.
The settings<-
function can be amended to save specifc settings (
.norm_example_args
must be defined seperatly to normalize input arguments in
any way one sees fit).
setReplaceMethod(f = "settings",
signature = signature(x = "ModExample"),
definition = function(x, value){
x <- callNextMethod()
# validate special setting here
x@settings[names(value)] <- unname(.norm_example_args(value))
x
})
The aggregateData
function is used to take the aggregated data from the
SequenceData
object and to calculate the specific scores, which are then
stored in the aggregate
slot.
setMethod(f = "aggregateData",
signature = signature(x = "ModExample"),
definition =
function(x, force = FALSE){
# Some data with element per transcript
}
)
The findMod
function takes the aggregate data and searches for modifications,
which are then returned as a GRanges object and stored in the modifications
slot.
setMethod("findMod",
signature = c(x = "ModExample"),
function(x){
# an element per modification found.
}
)
ModifierSet
classThe ModifierSet
class is implemented very easily by defining the class and
the constructor. The functionality is defined by the Modifier
class.
setClass("ModSetExample",
contains = "ModifierSet",
prototype = list(elementType = "ModExample"))
ModSetExample <- function(x, annotation, sequences, seqinfo, ...){
RNAmodR:::ModifierSet("ModExample", x = x, annotation = annotation,
sequences = sequences, seqinfo = seqinfo, ...)
}
Additional functions, which need to be implemented, are getDataTrack
for the
new SequenceData
and new Modifier
classes and
plotData
/plotDataByCoord
for the new Modifier
and ModifierSet
classes. name
defines a transcript name found in names(ranges(x))
and
type
is the data type typically found as a column in the aggregate
slot.
setMethod(
f = "getDataTrack",
signature = signature(x = "ExampleSequenceData"),
definition = function(x, name, ...) {
###
}
)
setMethod(
f = "getDataTrack",
signature = signature(x = "ModExample"),
definition = function(x, name, type, ...) {
}
)
setMethod(
f = "plotDataByCoord",
signature = signature(x = "ModExample", coord = "GRanges"),
definition = function(x, coord, type = "score", window.size = 15L, ...) {
}
)
setMethod(
f = "plotData",
signature = signature(x = "ModExample"),
definition = function(x, name, from, to, type = "score", ...) {
}
)
setMethod(
f = "plotDataByCoord",
signature = signature(x = "ModSetExample", coord = "GRanges"),
definition = function(x, coord, type = "score", window.size = 15L, ...) {
}
)
setMethod(
f = "plotData",
signature = signature(x = "ModSetExample"),
definition = function(x, name, from, to, type = "score", ...) {
}
)
If unsure, how to modify these functions, have a look a the code in the
Modifier-Inosine-viz.R
file of this package.
As suggested directly above, for a more detailed example have a look at the
ModInosine
class source code found in the Modifier-Inosine-class.R
and
Modifier-Inosine-viz.R
files of this package.
sessionInfo()
## R Under development (unstable) (2024-10-21 r87258)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB 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: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] RNAmodR_1.21.0 Modstrings_1.23.0 RNAmodR.Data_1.19.0
## [4] ExperimentHubData_1.33.0 AnnotationHubData_1.37.0 futile.logger_1.4.3
## [7] ExperimentHub_2.15.0 AnnotationHub_3.15.0 BiocFileCache_2.15.0
## [10] dbplyr_2.5.0 txdbmaker_1.3.0 GenomicFeatures_1.59.0
## [13] AnnotationDbi_1.69.0 Biobase_2.67.0 Rsamtools_2.23.0
## [16] Biostrings_2.75.0 XVector_0.47.0 rtracklayer_1.67.0
## [19] GenomicRanges_1.59.0 GenomeInfoDb_1.43.0 IRanges_2.41.0
## [22] S4Vectors_0.45.0 BiocGenerics_0.53.0 BiocStyle_2.35.0
##
## loaded via a namespace (and not attached):
## [1] BiocIO_1.17.0 bitops_1.0-9
## [3] filelock_1.0.3 tibble_3.2.1
## [5] graph_1.85.0 XML_3.99-0.17
## [7] rpart_4.1.23 lifecycle_1.0.4
## [9] httr2_1.0.5 lattice_0.22-6
## [11] ensembldb_2.31.0 OrganismDbi_1.49.0
## [13] backports_1.5.0 magrittr_2.0.3
## [15] Hmisc_5.2-0 sass_0.4.9
## [17] rmarkdown_2.28 jquerylib_0.1.4
## [19] yaml_2.3.10 RUnit_0.4.33
## [21] Gviz_1.51.0 DBI_1.2.3
## [23] RColorBrewer_1.1-3 abind_1.4-8
## [25] zlibbioc_1.53.0 purrr_1.0.2
## [27] AnnotationFilter_1.31.0 biovizBase_1.55.0
## [29] RCurl_1.98-1.16 nnet_7.3-19
## [31] VariantAnnotation_1.53.0 rappdirs_0.3.3
## [33] GenomeInfoDbData_1.2.13 AnnotationForge_1.49.0
## [35] codetools_0.2-20 DelayedArray_0.33.0
## [37] xml2_1.3.6 tidyselect_1.2.1
## [39] UCSC.utils_1.3.0 farver_2.1.2
## [41] matrixStats_1.4.1 base64enc_0.1-3
## [43] GenomicAlignments_1.43.0 jsonlite_1.8.9
## [45] Formula_1.2-5 tools_4.5.0
## [47] progress_1.2.3 stringdist_0.9.12
## [49] Rcpp_1.0.13 glue_1.8.0
## [51] gridExtra_2.3 SparseArray_1.7.0
## [53] BiocBaseUtils_1.9.0 xfun_0.48
## [55] MatrixGenerics_1.19.0 dplyr_1.1.4
## [57] withr_3.0.2 formatR_1.14
## [59] BiocManager_1.30.25 fastmap_1.2.0
## [61] latticeExtra_0.6-30 fansi_1.0.6
## [63] digest_0.6.37 R6_2.5.1
## [65] mime_0.12 colorspace_2.1-1
## [67] jpeg_0.1-10 dichromat_2.0-0.1
## [69] biomaRt_2.63.0 RSQLite_2.3.7
## [71] utf8_1.2.4 generics_0.1.3
## [73] data.table_1.16.2 prettyunits_1.2.0
## [75] httr_1.4.7 htmlwidgets_1.6.4
## [77] S4Arrays_1.7.0 pkgconfig_2.0.3
## [79] gtable_0.3.6 blob_1.2.4
## [81] htmltools_0.5.8.1 bookdown_0.41
## [83] RBGL_1.83.0 ProtGenerics_1.39.0
## [85] scales_1.3.0 png_0.1-8
## [87] colorRamps_2.3.4 knitr_1.48
## [89] lambda.r_1.2.4 rstudioapi_0.17.1
## [91] reshape2_1.4.4 rjson_0.2.23
## [93] checkmate_2.3.2 curl_5.2.3
## [95] biocViews_1.75.0 cachem_1.1.0
## [97] stringr_1.5.1 BiocVersion_3.21.1
## [99] parallel_4.5.0 foreign_0.8-87
## [101] restfulr_0.0.15 pillar_1.9.0
## [103] grid_4.5.0 vctrs_0.6.5
## [105] cluster_2.1.6 htmlTable_2.4.3
## [107] evaluate_1.0.1 tinytex_0.53
## [109] magick_2.8.5 cli_3.6.3
## [111] compiler_4.5.0 futile.options_1.0.1
## [113] rlang_1.1.4 crayon_1.5.3
## [115] labeling_0.4.3 interp_1.1-6
## [117] plyr_1.8.9 stringi_1.8.4
## [119] deldir_2.0-4 BiocParallel_1.41.0
## [121] BiocCheck_1.43.0 munsell_0.5.1
## [123] lazyeval_0.2.2 Matrix_1.7-1
## [125] BSgenome_1.75.0 hms_1.1.3
## [127] bit64_4.5.2 ggplot2_3.5.1
## [129] KEGGREST_1.47.0 SummarizedExperiment_1.37.0
## [131] highr_0.11 ROCR_1.0-11
## [133] memoise_2.0.1 bslib_0.8.0
## [135] bit_4.5.0