## ----include=TRUE, echo=FALSE, results = 'asis', warning=FALSE---------------- library(knitr) library(pander) fun_descr_xnastring <- data.frame("function" = c('XNAString', 'XNAReverseComplement', 'predictMfeStructure', 'predictDuplexStructure', 'XNAStringToHelm'), description = c('Create XNAString object by passing at least base (sugar, name, backbone, target, conjugates and dictionary are optional).', 'Take XNAString object and return reverse complement of base slot.', 'Take XNAString object and apply RNAfold_MFE function from ViennaRNA package on base slot (if single stranded molecule).', 'Take XNAString object and apply RNAcofold_MFE function from ViennaRNA package on base slot (if double stranded molecule) or duplicate base and apply RNAcofold_MFE function from ViennaRNA (if single stranded molecule)', 'Change XNAString/XNAStringSet object to helm notation.')) #m <- mtcars[1:3, 1:4] pandoc.table(fun_descr_xnastring, keep.line.breaks = TRUE, style='rmarkdown', justify = 'left', caption = "XNAString class methods", split.table = Inf) # style = 'rmarkdown' ## ----include=TRUE, echo=FALSE, results = 'asis', warning=FALSE---------------- library(knitr) library(pander) fun_descr_xnastringset <- data.frame("function" = c('XNAStringSet', 'set2List', 'set2Dt', 'dt2Set', '[', '[[' ), description = c('Create XNAStringSet object by passing XNAString objects list.', 'Change XNAStringSet object to list of XNAString objects.', 'Change XNAStringSet object to data.table.', 'Change data.table (or data.frame) to XNAStringSet object', 'Extract single/multiple XNAString objects (XNAStringSet object returned) by passing index/indexes number or name/names.', 'Extract single XNAString object (XNAString object returned) by passing index number or name.' )) #m <- mtcars[1:3, 1:4] pandoc.table(fun_descr_xnastringset, keep.line.breaks = TRUE, style='rmarkdown', justify = 'left', caption = "XNAStringSet class methods", split.table = Inf) # style = 'rmarkdown' ## ----include=TRUE, echo=FALSE, results = 'asis', warning=FALSE---------------- library(knitr) library(pander) fun_descr <- data.frame("function" = c( 'name / name<-', 'base / base<-', 'sugar / sugar<-', 'backbone / backbone<-', 'target / target<-', 'conjugate5 / conjugate5<-', 'conjugate3 / conjugate3<-', 'secondary_structure / secondary_structure<-', 'duplex_structure / duplex_structure<-', 'dictionary / dictionary<-', 'compl_dictionary / compl_dictionary<-', 'XNAStringFromHelm', 'XNAPairwiseAlignment', 'XNAMatchPattern', 'XNAVmatchPattern', 'XNAMatchPDict', 'XNAAlphabetFrequency', 'XNADinucleotideFrequency', 'mimir2XnaDict'), description = c( 'Extract / overwrite name slot.', 'Extract / overwrite base slot.', 'Extract / overwrite sugar slot.', 'Extract / overwrite backbone slot.', 'Extract / overwrite target slot.', 'Extract / overwrite conjugate5 slot.', 'Extract / overwrite conjugate3 slot.', 'Extract / overwrite secondary_structure slot.', 'Extract / overwrite duplex_structure slot.', 'Extract / overwrite dictionary slot.', 'Extract / overwrite compl_dictionary slot.', 'Change helm notation to XNAString/XNAStringSet object.', 'Inherited from Biostrings package. Solve global/local/ends-free alignment problems.', 'Inherited from Biostrings package. Find/count all the occurrences of a given pattern (typically short) in a reference sequence (typically long). Support mismatches and indels.', 'Inherited from Biostrings package. Find/count all the occurrences of a given pattern (typically short) in a set of reference sequences. Support mismatches and indels.', 'Inherited from Biostrings package. Find/count all the occurrences of a set of patterns in a reference sequence. Support a small number of mismatches.', 'Tabulate the letters and count frequency for nucleotides.', 'Tabulate the letters and count frequency for dinucleotides.', 'Reformat mimir table to XNA dictionary standards')) #m <- mtcars[1:3, 1:4] pandoc.table(fun_descr, keep.line.breaks = TRUE, style='rmarkdown', justify = 'left', caption = "XNAString and XNAStringSet class methods", split.table = Inf) # style = 'rmarkdown' ## ----include=TRUE, echo=FALSE, results = 'asis', warning=FALSE---------------- library(knitr) library(pander) fun_descr_other <- data.frame("function" = c('concatDict', 'mimir2XnaDict'), description = c('Concatanate custom HELM-symbol dictionary with built-in HELM-symbol xna_dictionary.', 'Rewrite dictionary table to standard format.' )) #m <- mtcars[1:3, 1:4] pandoc.table(fun_descr_other, keep.line.breaks = TRUE, style='rmarkdown', justify = 'left', caption = "XNAString class methods", split.table = Inf) # style = 'rmarkdown' ## ----include=FALSE, echo=FALSE------------------------------------------------ library(XNAString) library(Biostrings) XNAString::xna_dictionary ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'ATCG', sugar = 'OOOO') obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'ATCG', sugar = 'OOOO', backbone = 'SOS', target = Biostrings::DNAStringSet('TAGC'), conjugate3 = "[5gn2c6]", name = "oligo1") obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'ATCG', sugar = 'OOOO', backbone = 'SOS', target = Biostrings::DNAStringSet(c('TAGC', 'TATC'))) obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = c('ATCG', 'TAGC'), sugar = c('OOOO', 'OOOO'), backbone = c('SOS','SOS'), target = Biostrings::DNAStringSet(c('TAGC', 'TATC'))) obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- d1 <- Biostrings::DNAString(x = "ACGATCG") obj <- XNAString(base = d1) obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString(base = 'ATCG') XNAString(base = 'ATCG', default_sugar = 'F', default_backbone = 'X') XNAString(base = Biostrings::DNAString('ATCG')) XNAString(base = Biostrings::DNAString('ATCG'), default_backbone = 'O', default_sugar = 'R') ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'ATCG', sugar = 'FODL', backbone = 'SOO') base(obj) base(obj) <- 'CTGA' base(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString::complementary_bases obj1 <- XNAString(base = "ACEGTTGGT", sugar = 'FODDDDDDD', conjugate3 = 'TAG') XNAString::XNAReverseComplement(obj1) ## ----include=TRUE, echo=TRUE-------------------------------------------------- custom_compl_dict <- rbind( XNAString::complementary_bases[seq(1, 5),], data.table::data.table( base = 'U', target = 'R', compl_target = 'T' ) ) obj2 <- XNAString(base = "ACGCUUA", sugar = 'DDDDDDD', compl_dictionary = custom_compl_dict) obj2 ## ----include=TRUE, echo=TRUE-------------------------------------------------- #create custom complementary dictonary with complementary bases coded as IUPAC compl_dict <- XNAString::complementary_bases compl_dict[base == "G"]$target <- "Y" compl_dict[base == "T"]$target <- "R" # if you have T in your base sequence compl_dict[base == "U"]$target <- "R" # if you have U in your base sequence compl_dict xna <- XNAString::XNAString(base = "ACGTACGT", sugar = "DDDDDDDD",compl_dictionary = compl_dict) xna ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'GAGAGGGAACCAGGCAGGGACCGCAGACAACA', sugar = 'FODLMFODLMFODLMFODLMFODLMFFFFFFF') XNAString::predictMfeStructure(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = Biostrings::DNAStringSet(c('GAGAGGGAACCAGGCAGGGACCGCAGACAACA', 'GAGAGGGAACCAGGCAGGGACCGCAGACAACA'))) XNAString::predictDuplexStructure(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = Biostrings::DNAString('GAGAGGGAACCAGGCAGGGACCGCAGACAACA'), sugar = 'FODLMFODLMFODLMFODLMFODLMFFFFFFF') XNAString::predictDuplexStructure(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = 'ATCG', sugar = 'FODD') XNAString_obj2 <- XNAString(base = 'TTCT', sugar = 'FOLL') XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAStringSet_obj ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = c('ATCG', 'TAGC'), sugar = c('OOOO', 'OOOO'), backbone = c('SOS','SOS'), target = Biostrings::DNAStringSet(c('TAGC', 'TACC'))) XNAString_obj2 <- XNAString(base = c('GGCG', 'TATC'), sugar = c('OOOO', 'OOOO'), target = Biostrings::DNAStringSet(c('CCGC', 'TATG'))) XNAStringSet_siRNA <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAStringSet_siRNA ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAStringSet(base = c("ATGCT","TGCAT","ATATG")) XNAStringSet(base = c("ATGCT","TGCAT","ATATG"), default_sugar = 'R', default_backbone = 'X') XNAStringSet(base = c("ATGCT","TGCAT","ATATG"),sugar = c("DDDDD","DDDDD","DDDDD")) XNAStringSet(base= list(c('TT', 'GG'), c('TG', 'GT'), c('TG')), sugar = list(c('FF', 'FO'), c('OO', 'OF'), c('OO')), backbone =list(c('X', 'X'), c('X', 'X'), c('X'))) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAStringSet(base= c('TT', 'GG'), sugar = c('FF', 'FF')) compl_dict <- XNAString::complementary_bases compl_dict[base == "G"]$target <- "Y" compl_dict[base == "T"]$target <- "R" # if you have T in your base sequence compl_dict[base == "U"]$target <- "R" # if you have U in your base sequence XNAStringSet(base= c('TT', 'GG'), sugar = c('FF', 'FF'), compl_dict = compl_dict) # compl_dict in use only if target empty XNAStringSet(base= c('TT', 'GG'), sugar = c('FF', 'FF'), target = c('AA', 'CC'), compl_dict = compl_dict) ## ----include=TRUE, echo=TRUE-------------------------------------------------- dt <- data.table::data.table(base = c('TT', 'GG')) out1 <- dt2Set(dt) out1 dt <- data.table::data.table(base = c('TT', 'GG'), default_sugar = 'R', default_backbone = 'X') out2 <- dt2Set(dt) out2 dt <- data.table::data.table(base= list(c('TT', 'GG'), c('TG', 'GT'), c('TG')), sugar = list(c('FF', 'FO'), c('OO', 'OF'), c('OO')), backbone =list(c('X', 'X'), c('X', 'X'), c('X'))) out3 <- dt2Set(dt) out3 ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(name = 'oligo1', base = 'ATCG', sugar = 'FODD') XNAString_obj2 <- XNAString(name = 'oligo2',base = 'TTCT', sugar = 'FOLL') XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) set2List(XNAStringSet_obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAStringSet_obj[2] XNAStringSet_obj['oligo2'] XNAStringSet_obj[c(1,2)] XNAStringSet_obj[c('oligo1', 'oligo2')] XNAStringSet_obj[[2]] XNAStringSet_obj[['oligo2']] ## ----include=TRUE, echo=TRUE-------------------------------------------------- set2Dt(XNAStringSet_obj, slots =c('name', 'base', 'sugar', 'backbone', 'target', 'conjugate5', 'conjugate3')) ## ----include=TRUE, echo=TRUE-------------------------------------------------- base(XNAStringSet_siRNA, 1) base(XNAStringSet_siRNA, 2) base(XNAStringSet_siRNA, 1) <- c('CTGA', 'CCCC') base(XNAStringSet_siRNA, 2) <- c('TTTT', 'TTTT') XNAStringSet_siRNA ## ----include=TRUE, echo=TRUE-------------------------------------------------- dt <- data.table::data.table(base= c('AEACACACACACAEAE', 'AAETGCTETGTATTTTTE'), sugar = c('LLLDDDDDDDDDDLLL', 'LLLDDDDDDDDDDDLLLL'), backbone =c('SSSSSSSSSSSSSSS', 'SSSSSSSSSSSSSSSSS')) dt dt2Set(dt) ## ----include=TRUE, echo=TRUE-------------------------------------------------- helm <- "CHEM1{[5gn2c6]}|RNA1{P.[dR](C)P.[dR](A)P.[LR](G)[sP].[LR](A)[sP].[LR](G)[sP].[LR](A)[sP].[dR](A)[sP].[dR](G)[sP].[dR](G)[sP].[dR](C)[sP].[dR](A)[sP].[dR](C)[sP].[dR](A)[sP].[dR](G)[sP].[dR](A)[sP].[LR]([5meC])[sP].[LR](G)[sP].[LR](G)}$CHEM1,RNA1,1:R2-1:R1$$$V2.0" XNAStringFromHelm(helm, name ='oligo1') ## ----include=TRUE, echo=TRUE-------------------------------------------------- helm <- "RNA1{[mR](C)P.[mR](C)P.[fR](C)P.[mR](C)P.[fR](C)P.[mR](G)P.[fR](C)P.[mR](C)P.[fR](G)P.[mR](T)P.[fR](G)P.[mR](G)P.[fR](T)P.[mR](T)P.[fR](C)P.[mR](A)P.[fR](T)P.[mR](A)P.[fR](A)}|RNA2{[fR](T)P.[fR](T)P.[mR](A)P.[fR](T)P.[mR](G)P.[fR](A)P.[mR](A)P.[fR](C)P.[mR](C)P.[fR](A)P.[mR](C)P.[fR](G)P.[mR](G)P.[fR](C)P.[mR](A)P.[fR](G)P.[mR](G)P.[fR](G)P.[mR](G)P.[fR](C)P.[mR](G)}$RNA1,RNA2,2:pair-56:pair|RNA1,RNA2,5:pair-53:pair|RNA1,RNA2,8:pair-50:pair|RNA1,RNA2,11:pair-47:pair|RNA1,RNA2,14:pair-44:pair|RNA1,RNA2,17:pair-41:pair|RNA1,RNA2,20:pair-38:pair|RNA1,RNA2,23:pair-35:pair|RNA1,RNA2,26:pair-32:pair|RNA1,RNA2,29:pair-29:pair|RNA1,RNA2,32:pair-26:pair|RNA1,RNA2,35:pair-23:pair|RNA1,RNA2,38:pair-20:pair|RNA1,RNA2,41:pair-17:pair|RNA1,RNA2,44:pair-14:pair|RNA1,RNA2,47:pair-11:pair|RNA1,RNA2,50:pair-8:pair|RNA1,RNA2,53:pair-5:pair|RNA1,RNA2,56:pair-2:pair$$$V2.0" XNAStringFromHelm(helm) ## ----include=TRUE, echo=TRUE-------------------------------------------------- helm <- "CHEM1{[5gn2c6]}|RNA1{P.[dR](C)P.[dR](A)P.[LR](T)[sP].[LR](T)[sP].[LR](G)[sP].[dR](A)[sP].[dR](A)[sP].[dR](T)[sP].[dR](A)[sP].[dR](A)[sP].[dR](G)[sP].[dR](T)[sP].[dR](G)[sP].[dR](G)[sP].[dR](A)[sP].[LR](T)[sP].[LR](G)[sP].[LR](T)}$CHEM1,RNA1,1:R2-1:R1$$$V2.0" XNAStringFromHelm(helm) ## ----include=TRUE, echo=TRUE-------------------------------------------------- helm <- c("CHEM1{[5gn2c6]}|RNA1{P.[dR](C)P.[dR](A)P.[LR](G)[sP].[LR](A)[sP].[LR](G)[sP].[LR](A)[sP].[dR](A)[sP].[dR](G)[sP].[dR](G)[sP].[dR](C)[sP].[dR](A)[sP].[dR](C)[sP].[dR](A)[sP].[dR](G)[sP].[dR](A)[sP].[LR]([5meC])[sP].[LR](G)[sP].[LR](G)}$CHEM1,RNA1,1:R2-1:R1$$$V2.0", "RNA1{[LR](T)[sP].[LR](G)[sP].[dR](T)[sP].[dR](G)[sP].[LR](T)[sP].[LR](G)[sP].[dR](T)[sP].[dR](G)[sP].[LR](T)[sP].[LR](G)[sP].[dR](T)[sP].[dR](G)[sP].[LR](T)[sP].[LR](G)[sP].[LR](T)}$$$$V2.0") XNAStringFromHelm(helm, name =c('oligo1', 'oligo2')) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'GAGTTACTTGCCAAET', sugar = 'LLLDMDDDDDDDDLLL', backbone = 'XXXXXXXXXXXXXX2') XNAStringToHelm(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString( base = Biostrings::DNAStringSet(c("CCCC", "GGGG")), sugar = c("OOFO", "FFOF"), backbone = c("OOO", "OOO"), target = '', conjugate3 = "", conjugate5 = "") XNAStringToHelm(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString( base = c("CCCCEGC", "UUAUGAT"), sugar = c("OOFOFOF", "FFOFOFO"), backbone = c("OOOOOO", "OOOOOO"), target = '', conjugate3 = "[5gn2c6]", conjugate5 = "") XNAStringToHelm(obj) ## ----include=TRUE, echo=TRUE-------------------------------------------------- obj <- XNAString(base = 'ATCGATATATATACACATGTATGATG', sugar = 'OOOODDDDDDDDDDDDDDDDDDDDDD', target = DNAStringSet(c('TAGCTATATATATGTGTACATACTAC', 'TAGCTAGATATATGTGTACATACTAC'))) subject <- 'ATCGATATATATACACATGTATGATGTAGCTATATATATGTGTACATACTACATCGATATATATACACATGTATGATG' ## ----include=TRUE, echo=TRUE-------------------------------------------------- substitutionMatrix <- Biostrings::nucleotideSubstitutionMatrix() XNAString::XNAPairwiseAlignment(pattern = obj, subject = subject, substitutionMatrix = substitutionMatrix) ## ----include=TRUE, echo=TRUE-------------------------------------------------- substitutionMatrix <- Biostrings::nucleotideSubstitutionMatrix() XNAString::XNAPairwiseAlignment(pattern = obj, subject = subject, type = "local", substitutionMatrix = substitutionMatrix) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString::XNAMatchPattern(pattern = obj, subject = subject) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString::XNAMatchPattern(pattern = obj, subject = subject, target.number = 2) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString::XNAMatchPattern(pattern = obj, subject = subject, target.number = 2, max.mismatch = 1) ## ----include=TRUE, echo=TRUE-------------------------------------------------- subject <- c('ATCGATATATATACACATGTATGATGTAGCTATATATATGTGTACATACTACATCGATATATATACACATGTATGATG', 'ATCGATATATATACACATGTATGATGTAGCTATATATATGTGTGCGACTACATCGATATATATACACATGTATGATG') XNAString::XNAVmatchPattern(pattern = obj, subject) ## ----include=TRUE, echo=TRUE-------------------------------------------------- subject <- 'ATCGATATATATACACATGTATGATGTAGCTATATATATGTGTACATACTACATCGATATATATACACATGTATGATG' XNAString::XNAMatchPDict(pdict = obj, subject = subject) ## ----include=TRUE, echo=TRUE-------------------------------------------------- xnastring_obj <- XNAString(name = 'oligo1', base = c('AEEE'), sugar = c('FFOO'), target = DNAStringSet('TTT')) XNAString::XNAAlphabetFrequency(obj = xnastring_obj, slot = 'base') XNAString::XNAAlphabetFrequency(obj = xnastring_obj, slot = 'base', as.prob = TRUE) XNAString::XNAAlphabetFrequency(obj = xnastring_obj, slot = 'base', base_only = TRUE) XNAString::XNAAlphabetFrequency(obj = xnastring_obj, slot = 'base', letters = c('A', 'C')) ## ----include=TRUE, echo=TRUE-------------------------------------------------- xnastring_obj <- XNAString(name = 'oligo1', base = c('AAEC', 'ECTA'), sugar = c('FFOO', 'DDLM')) XNAString::XNAAlphabetFrequency(obj = xnastring_obj, slot = 'sugar', matrix_nbr = 2) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = 'ATCG', sugar = 'FODD') XNAString_obj2 <- XNAString(base = 'TTCT', sugar = 'FOLL') XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAString::XNAAlphabetFrequency(obj = XNAStringSet_obj, slot = 'base') ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = c('ATCG', 'TAGC'), sugar = c('OFOO', 'ODDF'), backbone = c('SOS','SOS')) XNAString_obj2 <- XNAString(base = c('GGCG', 'TATC'), sugar = c('OOOO', 'OOFO')) XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAString::XNAAlphabetFrequency(obj = XNAStringSet_obj, slot = 'sugar', matrix_nbr = 2) ## ----include=TRUE, echo=TRUE-------------------------------------------------- xnastring_obj <- XNAString(name = 'oligo1', base = c('GCGC'), sugar = c('FODL'), target = DNAStringSet('TTTT')) XNAString::XNADinucleotideFrequency(obj = xnastring_obj, slot = 'base') XNAString::XNADinucleotideFrequency(obj = xnastring_obj, slot = 'base', as.prob = TRUE) XNAString::XNADinucleotideFrequency(obj = xnastring_obj, slot = 'base', base_only = TRUE) XNAString::XNADinucleotideFrequency(obj = xnastring_obj, slot = 'base', double_letters = c('GC', 'CU')) ## ----include=TRUE, echo=TRUE-------------------------------------------------- my_dict <- data.table::data.table(type = c(rep('base', 3), rep('sugar', 2), rep('backbone', 2)), symbol = c('G', 'E', 'A', 'F', 'O', 'S', 'B')) xnastring_obj <- XNAString_obj1 <- XNAString( base = c('AGGE', 'EEEA'), sugar = c('FFFO', 'OOOO'), backbone = c('SBS', 'SBS'), dictionary = my_dict ) XNAString::XNADinucleotideFrequency(obj = xnastring_obj, slot = 'sugar', matrix_nbr = 2) ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = 'ATCG', sugar = 'FODD') XNAString_obj2 <- XNAString(base = 'TTCT', sugar = 'FOLL') XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAString::XNADinucleotideFrequency(obj = XNAStringSet_obj, slot = 'base') ## ----include=TRUE, echo=TRUE-------------------------------------------------- XNAString_obj1 <- XNAString(base = c('ATCG', 'TAGC'), sugar = c('OFOO', 'ODDF'), backbone = c('SOS','SOS')) XNAString_obj2 <- XNAString(base = c('GGCG', 'TATC'), sugar = c('OOOO', 'OOFO')) XNAStringSet_obj <- XNAStringSet(objects = list(XNAString_obj1, XNAString_obj2)) XNAString::XNADinucleotideFrequency(obj = XNAStringSet_obj, slot = 'sugar', matrix_nbr = 2) ## ----include=TRUE, echo=TRUE-------------------------------------------------- dt <- data.table::data.table(HELM = c("([PPG])", "[fR]", "[srP]"), TS_BASE_SEQ = c("F", NA, NA), TS_SUGAR_SEQ = c(NA, NA, 'F'), TS_BACKBONE_SEQ = c(NA, 'S', NA)) dt mimir2XnaDict(dt, 'TS_BASE_SEQ', 'TS_SUGAR_SEQ', 'TS_BACKBONE_SEQ') ## ----include=TRUE, echo=TRUE-------------------------------------------------- my_dict <- data.table::data.table(HELM = c('[[B]]'), type = c('base'), symbol = c('B')) concatDict(my_dict) ## ----include=TRUE, echo=TRUE, error=TRUE-------------------------------------- my_dict <- data.table::data.table(HELM = c('[[U]]'), type = c('base'), symbol = c('U')) concatDict(my_dict) ## ----include = FALSE---------------------------------------------------------- knitr::opts_chunk$set( collapse = TRUE, comment = "#>" ) ## ----include=TRUE, echo=TRUE, error=TRUE-------------------------------------- helm <- "CHEM1{[5gn2c6]}|RNA1{P.[dR](C)P.[dR](A)P.[LR](G)[sP].[LR](A)[sP].[LR](G)[sP].[LR](A)[sP].[dR](A)[sP].[dR](G)[sP].[dR](G)[sP].[dR](C)[sP].[dR](A)[sP].[dR](C)[sP].[dR](A)[sP].[dR](G)[sP].[dR](A)[sP].[LR]([5meC])[sP].[LR](G)[sP].[LR](G)}$CHEM1,RNA1,1:R2-1:R1$$$V2.0" xna_obj <- XNAStringFromHelm(helm, name ='oligo1') xna_obj subject <- 'ATCGATATATATACACCGTCTGTGCCTTCTCACTACATCGAG' substitutionMatrix <- Biostrings::nucleotideSubstitutionMatrix() XNAString::XNAPairwiseAlignment(pattern = obj, subject = subject, substitutionMatrix = substitutionMatrix) XNAString::XNAMatchPattern(pattern = xna_obj, subject = subject) ## ----include=TRUE, echo=TRUE, error=TRUE-------------------------------------- base <- list(c('ATCGATAT', 'ATCGATAT'), c('TGGGGGTGC', 'ATCGGGAT'), c('CCCTAGTA')) set_obj <- XNAStringSet(base = base) set_obj XNAAlphabetFrequency(obj = set_obj, slot = 'base') XNAAlphabetFrequency(obj = set_obj, slot = 'base', matrix_nbr = 2) XNAAlphabetFrequency(obj = set_obj, slot = 'base', as.prob = TRUE) XNADinucleotideFrequency(obj = set_obj, slot = 'base', double_letters = c('AT', 'GA', 'GT')) XNADinucleotideFrequency(obj = set_obj, slot = 'base', base_only = TRUE) ## ----include=TRUE, echo=TRUE, error=TRUE-------------------------------------- #create custom complementary dictonary with complementary bases coded as IUPAC compl_dict <- XNAString::complementary_bases compl_dict[base == "G"]$target <- "Y" compl_dict[base == "T"]$target <- "R" # if you have T in your base sequence compl_dict[base == "U"]$target <- "R" # if you have U in your base sequence compl_dict xna <- XNAString::XNAStringSet(base = c("ACGTACG", "CCCGTAC", "AATACTT"), compl_dict = compl_dict) xna ## ----sessionInfo, echo=FALSE-------------------------------------------------- sessionInfo()