## ----style, echo = FALSE, results = 'asis'---------------------------------
BiocStyle::markdown(css.files = c('custom.css'))

## ---- echo=FALSE-----------------------------------------------------------
suppressPackageStartupMessages({
  library(Modstrings)
  library(GenomicRanges)
})

## ---- eval=FALSE-----------------------------------------------------------
#  library(Modstrings)
#  library(GenomicRanges)

## ----object_creation2------------------------------------------------------
r <- RNAString("ACGUG")
mr2 <- modifyNucleotides(r,5,"m7G")
mr2
mr3 <- modifyNucleotides(r,5,"7G",nc.type = "nc")
mr3

## ----object_creation3------------------------------------------------------
mr4 <- ModRNAString(paste0("ACGU",alphabet(ModRNAString())[33]))
mr4

## ----object_creation4------------------------------------------------------
gr <- GRanges("1:5", mod = "m7G")
mr5 <- combineIntoModstrings(r, gr)
mr5

## ----object_creation5------------------------------------------------------
rs <- RNAStringSet(list(r,r,r,r,r))
names(rs) <- paste0("Sequence", seq_along(rs))
gr2 <- GRanges(seqnames = names(rs)[c(1,1,2,3,3,4,5,5)],
               ranges = IRanges(start = c(4,5,5,4,5,5,4,5),width = 1),
               mod = c("D","m7G","m7G","D","m7G","m7G","D","m7G"))
gr2
mrs <- combineIntoModstrings(rs, gr2)
mrs

## ----object_separation-----------------------------------------------------
gr3 <- separate(mrs)
rs2 <- RNAStringSet(mrs)
gr3
rs2

## ----object_comparison_teaser----------------------------------------------
r <- RNAString("ACGUG")
mr2 <- modifyNucleotides(r,5,"m7G")
r == RNAString(mr2)

## ----object_comparison-----------------------------------------------------
r == ModRNAString(r)
r == mr
rs == ModRNAStringSet(rs)
rs == c(mrs[1:3],rs[4:5])

## ----object_conversion-----------------------------------------------------
RNAString(mr)

## ----object_qual-----------------------------------------------------------
qmrs <- QualityScaledModRNAStringSet(mrs,
                                     PhredQuality(c("!!!!h","!!!!h","!!!!h",
                                                    "!!!!h","!!!!h")))
qmrs

## ----object_qual_sep_combine-----------------------------------------------
qgr <- separate(qmrs)
qgr
combineIntoModstrings(mrs,qgr, with.qualities = TRUE)

## ----object_io-------------------------------------------------------------
writeModStringSet(mrs, file = "test.fasta")
# note the different function name. Otherwise empty qualities will be written
writeQualityScaledModStringSet(qmrs, file = "test.fastq")
mrs2 <- readModRNAStringSet("test.fasta", format = "fasta")
mrs2
qmrs2 <- readQualityScaledModRNAStringSet("test.fastq")
qmrs2

## ----object_io_unlink, include=FALSE---------------------------------------
unlink("test.fasta")
unlink("test.fastq")

## ----object_pattern--------------------------------------------------------
matchPattern("U7",mr)
vmatchPattern("D7",mrs)
mrl <- unlist(mrs)
matchLRPatterns("7ACGU","U7ACG",100,mrl)

## ----cars------------------------------------------------------------------
# read the lines
test <- readLines(system.file("extdata","test.fasta",package = "Modstrings"),
                      encoding = "UTF-8")
head(test,2)
# keep every second line as sequence, the other one as name
names <- test[seq.int(from = 1, to = 104, by = 2)]
seq <- test[seq.int(from = 2, to = 104, by = 2)]
# sanitize input. This needs to be adapt to the individual case
names <- gsub(" ","_",
              gsub("> ","",
                   gsub(" \\| ","-",
                        names)))
seq <- gsub("-","",gsub("_","",seq))
names(seq) <- names

# sanitize special characters to Modstrings equivalent
seq <- sanitizeFromModomics(seq)
seq <- ModRNAStringSet(seq)
seq

# convert the contained modifications into a tabular format
separate(seq)

## ----sessioninfo-----------------------------------------------------------
sessioninfo::session_info()