R/RNAmodR.R
, R/AllGenerics.R
, R/Modifier-utils.R
, and 1 more
RNAmodR-development.Rd
These functions are not intended for general use, but are used for additional package development.
getData
is used to load data into a
SequenceData
object and must be
implented for all SequenceData
classes. The results must match the
requirements outlined in the value section.
In addition the following functions should be implemented for complete functionality:
aggregateData
for each SequenceData
and Modifier
class.
See also aggregateData
findMod
for each Modifier
class. See also
findMod
.
plotData
/plotDataByCoord
for each Modifier
and ModifierSet
class. See also
plotData
.
The following helper function can be called from within findMod
to
construct a coordinate for each modification found:
constructModRanges
constructs a GRanges
object describing the
location, type and associated scores of a modification.
constructModRanges
is typically called from the modify
function, which must be implemented for all
Modifier
classes.
constructModRanges(range, data, modType, scoreFun, source, type)
getData(x, bamfiles, grl, sequences, param, args)
# S4 method for GRanges,DataFrame
constructModRanges(range, data, modType, scoreFun, source, type)
for constructModRanges
: a GRanges
object
for constructModRanges
: a DataFrame
object
for constructModRanges
: a valid shortName for the
modification found. Must be present in shortName(ModRNAString())
.
for constructModRanges
: a custom function for
extracting scores from data
. The result must be a list
.
for constructModRanges
: a single character vector for
populating the source column of the result.
for constructModRanges
: a single character vector for
populating the source column of the result.
for getData
:a SequenceData
object.
for getData
:a BamFileList
object.
for getData
:a GRangesList
object.
for getData
:a XStringSet
object.
for getData
:a ScanBamParam
object.
for getData
: a list
with optional arguments.
getData
: returns a list with elements per BamFile in
bamfiles
. Elements can be
IntegerList
,
NumericList
or a
CompressedSplitDataFrameList
. The
data in the elements must be order by increasing positions numbers. However,
names and rownames will be discarded.
constructModRanges
: returns a GRanges
object with
genomic coordinates of modified nucleotides in the associated transcripts.
# new SequenceData class
setClass(Class = "ExampleSequenceData",
contains = "SequenceData",
prototype = list(minQuality = 5L))
ExampleSequenceData <- function(bamfiles, annotation, sequences, seqinfo, ...){
RNAmodR:::SequenceData("Example", bamfiles = bamfiles,
annotation = annotation, sequences = sequences,
seqinfo = seqinfo, ...)
}
setMethod("getData",
signature = c(x = "ExampleSequenceData",
bamfiles = "BamFileList",
grl = "GRangesList",
sequences = "XStringSet",
param = "ScanBamParam"),
definition = function(x, bamfiles, grl, sequences, param, args){
###
}
)
setMethod("aggregateData",
signature = c(x = "ExampleSequenceData"),
function(x, condition = c("Both","Treated","Control")){
###
}
)
setMethod(
f = "getDataTrack",
signature = c(x = "ExampleSequenceData"),
definition = function(x, name, ...) {
###
}
)
# new Modifier class
setClass("ModExample",
contains = "Modifier",
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, ...)
}
setMethod(f = "aggregateData",
signature = c(x = "ModExample"),
definition =
function(x, force = FALSE){
# Some data with element per transcript
}
)
setMethod("findMod",
signature = c(x = "ModExample"),
function(x){
# an element per modification found.
}
)
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", ...) {
}
)
# new ModifierSet 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, ...)
}
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", ...) {
}
)