With combineIntoModstrings and separate the construction and deconstruction of ModString Objects from an interacive session avoiding problematic encoding issues. In addition, modification information can be transfered from/to tabular data with these functions.

combineIntoModstrings expects seqnames(gr) or names(gr) to match the available names(x). Only information with strand information * and + are used.

separate when used with a GRanges/GRangesList object will return an object of the same type, but with modifications seperated. For example an element with mod = "m1Am" will be returned as two elements with mod = c("Am","m1A"). The reverse operation is available via combineModifications().

removeIncompatibleModifications filters incompatible modification from a GRanges or GRangesList. incompatibleModifications() returns the logical vector used for this operation.

separate(x, nc.type = "short")

combineIntoModstrings(
  x,
  gr,
  with.qualities = FALSE,
  quality.type = "Phred",
  stop.on.error = TRUE,
  verbose = FALSE,
  ...
)

combineModifications(gr, ...)

incompatibleModifications(gr, x, ...)

removeIncompatibleModifications(gr, x, ...)

# S4 method for ModString
separate(x, nc.type = c("short", "nc"))

# S4 method for ModStringSet
separate(x, nc.type = c("short", "nc"))

# S4 method for GRanges
separate(x)

# S4 method for GRangesList
separate(x)

# S4 method for XString,GRanges
combineIntoModstrings(
  x,
  gr,
  with.qualities = FALSE,
  quality.type = "Phred",
  stop.on.error = TRUE,
  verbose = FALSE,
  ...
)

# S4 method for XStringSet,GRangesList
combineIntoModstrings(
  x,
  gr,
  with.qualities = FALSE,
  quality.type = "Phred",
  stop.on.error = TRUE,
  verbose = FALSE,
  ...
)

# S4 method for XStringSet,GRanges
combineIntoModstrings(
  x,
  gr,
  with.qualities = FALSE,
  quality.type = "Phred",
  stop.on.error = TRUE,
  verbose = FALSE,
  ...
)

# S4 method for GRanges
combineModifications(gr)

# S4 method for GRangesList
combineModifications(gr)

# S4 method for GRanges,XString
incompatibleModifications(gr, x)

# S4 method for GRanges,XStringSet
incompatibleModifications(gr, x)

# S4 method for GRangesList,XStringSet
incompatibleModifications(gr, x)

# S4 method for GRanges,XString
removeIncompatibleModifications(gr, x)

# S4 method for GRanges,XStringSet
removeIncompatibleModifications(gr, x)

# S4 method for GRangesList,XStringSet
removeIncompatibleModifications(gr, x)

Arguments

x

For separate: a ModString/ModStringSet or GRanges/GRangesListobject

For combineIntoModstrings: a XString and a XStringSet object.

nc.type

the type of nomenclature to be used. Either "short" or "nc". "Short" for m3C would be "m3C", "nc" for m3C would be "3C". ( default = "short")

gr

a GRanges object

with.qualities

TRUE or FALSE (default): Should the values from a score column of the GRanges object stored? If set with.qualities = TRUE, combineIntoModstrings will try to construct a QualityScaledModStringSet object.

quality.type

the type of QualityXStringSet used, if with.qualities = TRUE. Must be on of the following values: "Phred","Solexa","Illumina".

stop.on.error

For combineIntoModstrings: TRUE(default) or FALSE: Should an error be raised upon encounter of incompatible positions?

verbose

For combineIntoModstrings: TRUE or FALSE (default): Should verbose information reported on the positions filled with modifications? This settings is passed onto modifyNucleotides.

...
  • default.quality: For combineIntoModstrings: the default.quality default value for non-modified positions. (default: default.quality = 0L)

Value

for separate a GRanges object and for combineIntoModstrings a ModString* object or a QualityScaledModStringSet, if with.qualities = TRUE.

Examples

library(GenomicRanges)
# ModDNAString
seq <- ModDNAString(paste(alphabet(ModDNAString()), collapse = ""))
seq
#> 55-letter ModDNAString object
#> seq: ACGTN-+.pδO]DJeg`bU∝πI763218∉⊆⊇Rαmh×f4νX'κo()ηa⇓⇑"√/≡ζ~

gr <- separate(seq)
gr
#> GRanges object with 47 ranges and 1 metadata column:
#>        seqnames    ranges strand |         mod
#>           <Rle> <IRanges>  <Rle> | <character>
#>    [1]        1         9      + |         5pT
#>    [2]        1        10      + |         3mT
#>    [3]        1        11      + |       O4meT
#>    [4]        1        12      + |         1mT
#>    [5]        1        13      + |         dhT
#>    ...      ...       ...    ... .         ...
#>   [43]        1        51      + |        6haA
#>   [44]        1        52      + |         2mA
#>   [45]        1        53      + |      2a6haA
#>   [46]        1        54      + |       6,6mA
#>   [47]        1        55      + |       6ncmA
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

seq2 <- combineIntoModstrings(as(seq,"DNAString"),gr)
seq2
#> 55-letter ModDNAString object
#> seq: ACGTN-+.pδO]DJeg`bU∝πI763218∉⊆⊇Rαmh×f4νX'κo()ηa⇓⇑"√/≡ζ~

seq == seq2
#> [1] TRUE
# ModRNAString
seq <- ModRNAString(paste(alphabet(ModRNAString()), collapse = ""))
seq
#> 153-letter ModRNAString object
#> seq: ACGUN-+.œεξα"KO]±Γ}/£ÿ≠≈*∞[ω∏%2:B#≤Z...Ϩ≡Ϫ√ΘE=6¿(eDς9⊄«⊆I⊇8y∑WPQšÎ÷H<;ÜVυY€

gr <- separate(seq)
gr
#> GRanges object with 145 ranges and 1 metadata column:
#>         seqnames    ranges strand |         mod
#>            <Rle> <IRanges>  <Rle> | <character>
#>     [1]        1         9      + |        m1Am
#>     [2]        1        10      + |        m1Gm
#>     [3]        1        11      + |        m1Im
#>     [4]        1        12      + |     m1acp3Y
#>     [5]        1        13      + |         m1A
#>     ...      ...       ...    ... .         ...
#>   [141]        1       149      + |          xU
#>   [142]        1       150      + |       cmo5U
#>   [143]        1       151      + |      mcmo5U
#>   [144]        1       152      + |          yW
#>   [145]        1       153      + |         imG
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

# Separating RNA modifications
gr <- gr[1]
separate(gr)
#> GRanges object with 2 ranges and 1 metadata column:
#>       seqnames    ranges strand |         mod
#>          <Rle> <IRanges>  <Rle> | <character>
#>   [1]        1         9      + |          Am
#>   [2]        1         9      + |         m1A
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

# ... and combine them again (both operations work only on a subset of
# modifications)
combineModifications(separate(gr))
#> GRanges object with 1 range and 1 metadata column:
#>       seqnames    ranges strand |         mod
#>          <Rle> <IRanges>  <Rle> | <character>
#>   [1]        1         9      + |        m1Am
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

# handling incompatible modifications
seq <- RNAString("AGCU")
gr <- GRanges(c("chr1:1:+","chr1:2:+"),mod="m1A")
incompatibleModifications(gr,seq)
#> [1] FALSE  TRUE

#
removeIncompatibleModifications(gr,seq)
#> GRanges object with 1 range and 1 metadata column:
#>       seqnames    ranges strand |         mod
#>          <Rle> <IRanges>  <Rle> | <character>
#>   [1]     chr1         1      + |         m1A
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths