makeEpiTxDb {EpiTxDb} | R Documentation |
Creating a EpiTxDb
from user supplied annotations as
data.frame
s
Description
makeEpiTxDb
is a low-level constructor for creating a
EpiTxDb
object from user supplied annotations.
This functions typically will not be used by regular users.
Usage
makeEpiTxDb(
modifications,
reactions = NULL,
specifiers = NULL,
references = NULL,
metadata = NULL,
reassign.ids = FALSE
)
Arguments
modifications |
A data.frame containg the following columns:
mod_id : a unique integer value per modification.
mod_type : the modification type as a character or
factor value. Must be a value from
shortName(ModRNAString()) .
mod_name : a character or factor name for the
specific modification
mod_start : the start position for the modification as
integer value. Usually mod_start = mod_end
mod_end : the end position for the modification as
integer value. Usually mod_start = mod_end
mod_strand : the strand information for the modificaion as a
character or factor .
sn_id : a integer value per unique sequence
sn_name : a character or factor as sequence
name, e.g a chromosome or a transcript identifier like chr1 .
The first six are mandatory, whereas one of the last two has to be set.
sn_id will be generated from sn_name , if sn_id is not
set.
|
reactions |
An optional data.frame containg the following
columns:
mod_id : a integer value per modification and the
link to the modification data.frame .
rx_genename : a character or factor referencing
a genename for the enzyme incorporating the modification.
rx_rank : a integer for sorting enzyme reactions, if
multiple enzymes are involved in the modification's
incorporation/maintenance.
rx_ensembl : a character or factor with an
ensembl identifier for the genename of the enzyme.
rx_ensembltrans : a character or factor with an
ensembl identifier for the transcript being translated into the enzyme.
rx_entrezid : a character or factor with an
entrezid for the genename of the enzyme.
(default: reactions = NULL )
|
specifiers |
An optional data.frame containg the following
columns:
mod_id : a integer value per modification and the
link to the modification data.frame .
spec_type : a character or factor
referencing a type of specifier, e.g. snoRNA . Not checked for
validity.
spec_genename : a character or factor
referencing a genename for the specifier directing an enzyme to the
specific location for the modification to be incorporated.
spec_ensembl : a character or factor with an
ensembl identifier for the genename of the specifier.
spec_ensembltrans : a character or factor with an
ensembl identifier for the transcript being translated into the specifier.
spec_entrezid : a character or factor with an
entrezid for the genename of the specifier.
(default: specifiers = NULL )
|
references |
An optional data.frame containg the following
columns:
mod_id : a integer value per modification and the
link to the modification data.frame .
ref_type : a character or factor with a
reference type, e.g. PMID . Is not checked for validity.
ref : a character or factor with a reference
value, e.g. a specific pubmed id or an journal article. Is not checked
for validity.
(default: references = NULL )
|
metadata |
An optional data.frame containg the following columns:
This dataframe will be returned
by metadata() (default: metadata = NULL )
|
reassign.ids |
TRUE or FALSE Controls how internal
mod_id s should be assigned. If reassign.ids is FALSE
(the default) and if the ids are supplied, then they are used as the
internal ids, otherwise the internal ids are assigned in a way that is
compatible with the order defined by ordering the modifications first by
chromosome, then by strand, then by start, and finally by end.
|
Value
a EpiTxDb
object.
See Also
Examples
mod <- data.frame("mod_id" = 1L,
"mod_type" = "m1A",
"mod_name" = "m1A_1",
"mod_start" = 1L,
"mod_end" = 1L,
"mod_strand" = "+",
"sn_id" = 1L,
"sn_name" = "test")
rx <- data.frame(mod_id = 1L,
rx_genename = "test",
rx_rank = 1L,
rx_ensembl = "test",
rx_ensembltrans = "test",
rx_entrezid = "test")
spec <- data.frame(mod_id = 1L,
spec_type = "test",
spec_genename = "test",
spec_ensembl = "test",
spec_ensembltrans = "test",
spec_entrezid = "test")
ref <- data.frame(mod_id = 1L,
ref_type = "test",
ref = "test")
etdb <- makeEpiTxDb(mod,rx,spec,ref)
[Package
EpiTxDb version 1.5.0
Index]