Title: | Secondary Structure Plotting for RNA |
---|---|
Description: | Functions for creating and manipulating RNA secondary structure plots. |
Authors: | Jim Maher [ths], JP Bida [aut, cre], Jonathan Price [ctb] |
Maintainer: | JP Bida <[email protected]> |
License: | GPL-3 |
Version: | 1.3 |
Built: | 2024-11-12 02:54:00 UTC |
Source: | https://github.com/jlp-bioinf/rrna |
Set of functions for creating and manipulating RNA secondary structure plots from CT files or bracket notations.
Package: | RRNA |
Type: | Package |
Version: | 1.0 |
Date: | 2015-07-27 |
License: | GPL-3 |
JP Bida Maintainer: JP Bida <[email protected]
### Create a CT file from bracket notation ct=makeCt("(((...(((...)))...(((...)))...)))","AAAUUUCCCAAAGGGUUUAAAGGGUUUCCCUUU") coord=ct2coord(ct) RNAPlot(coord,hl=c("GGGUUU","AAAUUU"),seqcols=c(2,4),labTF=TRUE)
### Create a CT file from bracket notation ct=makeCt("(((...(((...)))...(((...)))...)))","AAAUUUCCCAAAGGGUUUAAAGGGUUUCCCUUU") coord=ct2coord(ct) RNAPlot(coord,hl=c("GGGUUU","AAAUUU"),seqcols=c(2,4),labTF=TRUE)
Internal Function used in coord2comp to generate components from secondary structures
addBreaks(d,structure)
addBreaks(d,structure)
d |
coordinates generated from ct2coord |
structure |
Any character sequence |
Returns the character sequence provided with breaks indicated by "-"
JP Bida
Given a coordinate file with multiple RNA secondary structures, it aligns all folds such that n1 is at position (x,y) and n2 has its y coordinate equal to y
alignCoord(data, n1, n2, x, y)
alignCoord(data, n1, n2, x, y)
data |
R data frame containing the coordinates for plotting a given secondary structure |
n1 |
Nucleotide position that will be translated to (x,y) |
n2 |
Nucleotide position that will have its y coordinate equal to y |
x |
x coordinate that n1 will be translated to |
y |
y coordinate that n1 will be translated to |
Returns a data frame containing fold coordinates.
JP Bida
### Create two RNA secondary structures #### ct1=makeCt( "((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ct2=makeCt( "((((...(((((((....))))))).((..((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ### Create a coordinate file #### dat1=ct2coord(ct1) ### Each RNA fold needs its own id ### dat1$id=1 #### Create a coordinate file #### dat2=ct2coord(ct2) ### Each RNA fold needs its own id ### dat2$id=2 dat=rbind(dat1,dat2) adat=alignCoord(dat,1,46,0,0) ### Plot the aligned RNA folds #### RNAPlot(adat[adat$id==1,]) l=length(adat$seq[adat$id==2]) RNAPlot(adat[adat$id==2,],modspec=TRUE,modp=c(1:l),modcol=rep(4,l),mod=rep(16,l),add=TRUE)
### Create two RNA secondary structures #### ct1=makeCt( "((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ct2=makeCt( "((((...(((((((....))))))).((..((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ### Create a coordinate file #### dat1=ct2coord(ct1) ### Each RNA fold needs its own id ### dat1$id=1 #### Create a coordinate file #### dat2=ct2coord(ct2) ### Each RNA fold needs its own id ### dat2$id=2 dat=rbind(dat1,dat2) adat=alignCoord(dat,1,46,0,0) ### Plot the aligned RNA folds #### RNAPlot(adat[adat$id==1,]) l=length(adat$seq[adat$id==2]) RNAPlot(adat[adat$id==2,],modspec=TRUE,modp=c(1:l),modcol=rep(4,l),mod=rep(16,l),add=TRUE)
Generates and RNA secondary structure plot from a CT file. Removes pseudoKnots automatically and allows them to be drawn back in with pseudoTF=TRUE.
aptPlotCT(file, ranges = 0, add = FALSE, hl = NULL, seqcols = NULL, seqTF = FALSE, labTF = FALSE, nt = FALSE, dp = 0.5, modspec = FALSE, modp = NULL, mod = NULL, modcol = NULL, tsize = 0.5, main = "", pseudoTF = FALSE, pseudo_nums = NULL, ticks = NULL, ticksTF = FALSE )
aptPlotCT(file, ranges = 0, add = FALSE, hl = NULL, seqcols = NULL, seqTF = FALSE, labTF = FALSE, nt = FALSE, dp = 0.5, modspec = FALSE, modp = NULL, mod = NULL, modcol = NULL, tsize = 0.5, main = "", pseudoTF = FALSE, pseudo_nums = NULL, ticks = NULL, ticksTF = FALSE )
file |
CT file name |
ranges |
A data frame containing the ranges of sequence positions that should be highlighted with given colors.
|
add |
Should the new plot be added to an existing plot TRUE/FALSE |
hl |
Takes an array of sequences and highlights them with seqcol
|
seqcols |
Colors that should be used to highlight the sequences given in hl |
seqTF |
Unused argument included for backward compatibility |
labTF |
TRUE/FALSE plot the legend |
nt |
TRUE/FALSE plot the nucleotide sequence on the secondary structure |
dp |
Floating point value to determine how far from the coordinates the nucleotide sequence should be plotted. Values between 0 and 5 usually work best. |
modspec |
TRUE/FALSE modify specific positions in the secondary structure. Used in combination with modp,mod,and modcol. This allows you to change the shape and color of nucleotide in the secondary structure. |
modp |
Array defining the specific positions to be modified in the plot
|
mod |
Array defining the pch values to be plotted at the positions given by modp.
|
modcol |
Array of color values to be used for plotting at the positions defined by modp in the secondary structure.
|
tsize |
Text size used for plotting the nucleotide sequence in the secondary structure. Only applicable when nt=TRUE. Values between 0.1 and 4 work well. |
main |
Title used for the plot when labTF is set to TRUE. |
pseudoTF |
Plot pseudo knot sequences |
pseudo_nums |
indices of the nucleotides included in pseudoknots |
ticksTF |
TRUE/FALSE include ticks |
ticks |
Positions where the ticks should be drawn. These are sequence positions in the RNA molecule |
Returns and R plot object
JP Bida
### PseudoKnots ### pk= makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext = ".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile,row.names=FALSE,col.names=TRUE) aptPlotCT(ctfile,ticksTF=TRUE,ticks=seq(1,60,by=5),pseudoTF=TRUE,pseudo_nums=c(19,20,43,42))
### PseudoKnots ### pk= makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext = ".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile,row.names=FALSE,col.names=TRUE) aptPlotCT(ctfile,ticksTF=TRUE,ticks=seq(1,60,by=5),pseudoTF=TRUE,pseudo_nums=c(19,20,43,42))
Given a bracket notation for RNA secondary structure and an index of a ")" bracket type the backward function will find the "(" bracket that closes the ")" at the given index.
backward(stc, i)
backward(stc, i)
stc |
Array of brackets and dots.
|
i |
index giving the position of a bracket |
returns the index of the bracket that closes the bracket at the given index
JP Bida
a=unlist(strsplit("(((...)))...((..))","")) ind=backward(a,7)
a=unlist(strsplit("(((...)))...((..))","")) ind=backward(a,7)
A bpl file can be created from a given coordinate file for inputing into other RNA visulatization programs
bplfile(dat, name)
bplfile(dat, name)
dat |
Coordinate file created by ct2coord or loadCoords functions
|
name |
Name of the file outputed |
Creates the file with the given "name"
JP Bida
ct=makeCt("((((...))))","AAAACCCUUUU") ### Create the coordinate file ### dat=ct2coord(ct) bplfile(dat,tempfile(fileext = ".bpl"))
ct=makeCt("((((...))))","AAAACCCUUUU") ### Create the coordinate file ### dat=ct2coord(ct) bplfile(dat,tempfile(fileext = ".bpl"))
Given an integer N the function returns N (x,y) coordinates for a polygon with N sides each of length 1. This is used to plot the loops in an RNA structure
circleCoord(n)
circleCoord(n)
n |
Integer determining the number of sides |
Data frame with columns x,y defining coordinates of the polygons
JP Bida
pts=circleCoord(10) plot(pts$x,pts$y)
pts=circleCoord(10) plot(pts$x,pts$y)
Split RNA secondary structures into closed components that can be used for 3D assembly into larger molecules.
coord2comp(coord)
coord2comp(coord)
coord |
coordinates generated from ct2coord |
Data frame containing substructures
JP Bida
Creates a coordinate file from a CT file that has been loaded into a data frame
ct2coord(input)
ct2coord(input)
input |
Data frame representing a ct file. Created from makeCt or loadCt |
Returns a coordinate file for the secondary structure represented in the CT file
Pseudoknots sometimes cause trouble
JP Bida
ct=makeCt("((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) coord=ct2coord(ct) RNAPlot(coord)
ct=makeCt("((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) coord=ct2coord(ct) RNAPlot(coord)
Knet files are used as inputs for KnetFold secondary structure prediction program
ct2knet(file, ind = 0)
ct2knet(file, ind = 0)
file |
Name of the CT file being converted to KnetFold file |
ind |
Index used to relabel sequence indexes |
Retuns a string containing the contains of the knet file
JP Bida
pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext=".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile,row.names=FALSE,col.names=TRUE) ### Convert CT file to Knet ### out=ct2knet(ctfile,0)
pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext=".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile,row.names=FALSE,col.names=TRUE) ### Convert CT file to Knet ### out=ct2knet(ctfile,0)
Split RNA secondary structures into closed components that can be used for 3D assembly into larger molecules.
ct2ss(dat, cleanup)
ct2ss(dat, cleanup)
dat |
Dataframe from ctfile generated by loadCT |
cleanup |
TRUE/FALSE indicating if CT file should be cleaned up removing pseudoknots and multiple pairings preserving the most base-pairs in a structure compatible with bracket notation. |
Secondary structure in bracket notation
JP Bida
Given a bracket notation for RNA secondary structure and an index of a "(" bracket type the forward function will find the ")" bracket that closes the "(" at the given index.
forward(stc, i)
forward(stc, i)
stc |
a=unlist(strsplit("(((...)))...((..))","")) |
i |
Interger index |
Integer index
JP Bida
a=unlist(strsplit("(((...)))...((..))","")) ind=forward(a,1)
a=unlist(strsplit("(((...)))...((..))","")) ind=forward(a,1)
Given hydrogen bond data from a PDB file generate closed structures that can be used to assemble larger 3D structures
fuzzy2comp(dat)
fuzzy2comp(dat)
dat |
Dataframe containing hydrogen bond data generated from the RSIM tertiary structure program |
RNA secondary structure components
JP Bida
Given hydrogen bond data generated from a pdb file using the RSIM 3D structure prediction program, generate a bracket notation that does not contain pseudoknots that can be used as input for the RNAPlot function.
fuzzy2ct(dat)
fuzzy2ct(dat)
dat |
Dataframe with columns containing the data pdb - name of pdb file model - model number in the file pos1 - nucleotide 1 position in hydrogen bonded pair pos2 - nucleotide 2 position in hydrogen bonded pair face1 - face of NT 1 participating in the bond face2 - face of NT 2 participating in the bond h-bonds - probility assigned to number of hydrogen bonds |
returns secondary structure in bracket notation for the 3D structure
JP Bida
Generates coordinates for a loop in a secondary structure. Internal function used by RNAPlot.
genCords(loop, p1, p2, input, vn)
genCords(loop, p1, p2, input, vn)
loop |
List contianing a data frame that has the subset of nucleotides in a given loop |
p1 |
The position of the first nucleotide in the loop |
p2 |
The position of the second nucleotide in the loop |
input |
The data frame contianing the coordinate file for the entire RNA secondary structure |
vn |
A flag that flips over y axis if vn = 1. |
Returns a set of points
JP Bida
### This is an internal function ###
### This is an internal function ###
Adds BREAK to bracket notation in secondary structure indicating 5' and 3' discontinuities at a base-pair
genDB(map)
genDB(map)
map |
Dataframe generated by a custom script |
Returns inputed map with breaks indicated
JP Bida
Coordinate files can be created from the viennaRNA library.
loadCoords(filename)
loadCoords(filename)
filename |
Name of the coordinate file being loaded |
Data frame containing the coordinate file
JP Bida
The RRNAFold program generates the coordinate files used by RRNA
https://github.com/jpbida/ViennaScripts
### Create a test coordinate file using ct2coord ### ct=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) coord=ct2coord(ct) ### add an id ### coord$id=1 ### write out test file ### crfile = tempfile(fileext=".cr") write.table( coord[,c('id','x','y','seq','num','bound')], col.names=FALSE,row.names=FALSE,sep=",",file=crfile ) ### Read in the coordinate file ## input=loadCoords(crfile) ### Plot the file using RNAPlot ## RNAPlot(input)
### Create a test coordinate file using ct2coord ### ct=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) coord=ct2coord(ct) ### add an id ### coord$id=1 ### write out test file ### crfile = tempfile(fileext=".cr") write.table( coord[,c('id','x','y','seq','num','bound')], col.names=FALSE,row.names=FALSE,sep=",",file=crfile ) ### Read in the coordinate file ## input=loadCoords(crfile) ### Plot the file using RNAPlot ## RNAPlot(input)
A variety of RNA secondary structure prediction programs produce CT files. You can load these CT files into R using the loadCT function.
loadCt(file)
loadCt(file)
file |
The name of the CT file being loaded |
Returns at data frame containing the CT file data
JP Bida
### Create a CT file with PseudoKnots ### pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext = ".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile, row.names=FALSE,col.names=TRUE) ctfile=loadCt(ctfile) ### Before using ct2coord you need to remove the pseudo knots ### l=pseudoKnot(ctfile) dat=l[[2]] cd=ct2coord(dat) RNAPlot(cd)
### Create a CT file with PseudoKnots ### pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 ### Create a CT file for testing ### ctfile = tempfile(fileext = ".ct") write.table(pk[,c(1,4,2,3,6,5)],file=ctfile, row.names=FALSE,col.names=TRUE) ctfile=loadCt(ctfile) ### Before using ct2coord you need to remove the pseudo knots ### l=pseudoKnot(ctfile) dat=l[[2]] cd=ct2coord(dat) RNAPlot(cd)
Used by RNAPlot to get the length of a loop
loopLength(input, start)
loopLength(input, start)
input |
CT file
|
start |
Position of the first nucleotide in the the loop |
Retuns a list contianing the output and stems
JP Bida
ct=makeCt("((((...((((..))))..((((...)))).))))","AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA") out=loopLength(ct,4)
ct=makeCt("((((...((((..))))..((((...)))).))))","AAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAAA") out=loopLength(ct,4)
Given an RNA secondary structure in bracket notation containing no pseudoKnots this function creates an R data frame that represents the secondary structures CT file.
makeCt(struct, seq)
makeCt(struct, seq)
struct |
Bracket notation.
|
seq |
String containing the RNA sequence
|
Returns a data frame representing the bracket notaiton secondary structure in a CT file like format.
JP Bida
st="(((((....)))))..((..))" seq="AUAAUUAAAAAAAACCCCCAAA" ct=makeCt(st,seq)
st="(((((....)))))..((..))" seq="AUAAUUAAAAAAAACCCCCAAA" ct=makeCt(st,seq)
Given all hydrogen binding data from a PDB file identify the secondary structure that maximizes pairings using a greedy approach.
pairScores(PairDF)
pairScores(PairDF)
PairDF |
Dataframe generated from PDB file using RSIM. |
Returns PairDF file with pairings removed that do not contribute to the bracket notation.
JP Bida
internal function used to remove pseudoKnots before calling ct2coord
pseudoKnot(ctDat)
pseudoKnot(ctDat)
ctDat |
R data frame representing a CT file for RNA secondary structure |
Returns a list with the first item being a list of pseudoKnots and the second item being a CT file data frame with all pseudoKnots removed from the structure
JP Bida
pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC") pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 l=pseudoKnot(pk) ## Positions of removed pseudo knots ## removed=l[[1]] ### clean ct file that can be used by ct2coord ### ct=l[[2]]
pk=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC") pk$bound[pk$pos==20]=42 pk$bound[pk$pos==19]=43 pk$bound[pk$pos==43]=19 pk$bound[pk$pos==42]=20 l=pseudoKnot(pk) ## Positions of removed pseudo knots ## removed=l[[1]] ### clean ct file that can be used by ct2coord ### ct=l[[2]]
Given fold data from loadFolds or ct2coords RNAPlot plots the secondary structure
RNAPlot(data, ranges = 0, add = FALSE, hl = NULL, seqcols = NULL, seqTF = FALSE, labTF = FALSE, nt = FALSE, dp = 0.5, modspec = FALSE, modp = NULL, mod = NULL, modcol = NULL, tsize = 0.5, main = "", pointSize = 2, lineWd = 2)
RNAPlot(data, ranges = 0, add = FALSE, hl = NULL, seqcols = NULL, seqTF = FALSE, labTF = FALSE, nt = FALSE, dp = 0.5, modspec = FALSE, modp = NULL, mod = NULL, modcol = NULL, tsize = 0.5, main = "", pointSize = 2, lineWd = 2)
data |
R data frame containing the coordinates for plotting a given secondary structure
|
ranges |
A data frame containing the ranges of sequence positions that should be highlighted with given colors.
|
add |
Should the new plot be added to an existing plot TRUE/FALSE |
hl |
Takes an array of sequences and highlights them with seqcol
|
seqcols |
Colors that should be used to highlight the sequences given in hl |
seqTF |
Unused argument included for backward compatibility |
labTF |
TRUE/FALSE plot the legend |
nt |
TRUE/FALSE plot the nucleotide sequence on the secondary structure |
dp |
Floating point value to determine how far from the coordinates the nucleotide sequence should be plotted. Values between 0 and 5 usually work best. |
modspec |
TRUE/FALSE modify specific positions in the secondary structure. Used in combination with modp,mod,and modcol. This allows you to change the shape and color of nucleotide in the secondary structure. |
modp |
Array defining the specific positions to be modified in the plot
|
mod |
Array defining the pch values to be plotted at the positions given by modp.
|
modcol |
Array of color values to be used for plotting at the positions defined by modp in the secondary structure.
|
tsize |
Text size used for plotting the nucleotide sequence in the secondary structure. Only applicable when nt=TRUE. Values between 0.1 and 4 work well. |
main |
Title used for the plot when labTF is set to TRUE. |
pointSize |
The size of points plotted in the secondary structure. Values betwen 0.1-5 work well. |
lineWd |
Line width for base pairings and backbone of secondary structures. |
Returns a generic R plot that can be used with the jpeg, postscript, etc. functions.
JP Bida
## Create a CT file from bracket notation and sequence ### ct=makeCt( "((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ## Create a coordinate file based on the CT file ### dat=ct2coord(ct) ### Create a plot of the secondary structure ### RNAPlot(dat) ### Plot positions 1:4 as green and 43:46 circles ## ### and show the legend ranges=data.frame(min=c(1,43),max=c(4,46),col=c(2,3), desc=c("Region 1","Region 2") ) RNAPlot(dat,ranges,labTF=TRUE) ### Highlight the sequences CUCCU and CCCCAAA ### RNAPlot(dat,hl=c("CUCCU","CCCCAAA"),seqcol=c(2,4),labTF=TRUE,main="RNA Molecule") ### Modify specific positions #### RNAPlot( dat, modspec=TRUE, modp=c(1:4,43:46),mod=c(17,17,15,15,16,16,16,16), modcol=c(rep(2,2),rep(3,2),rep(4,4)) ) ### RNA Plot with nucleotides ### RNAPlot(dat,nt=TRUE) ### RNA plot with nucleotides RNAPlot( dat,nt=TRUE,modspec=TRUE,modp=c(1:4,43:46), mod=c(17,17,15,15,16,16,16,16), modcol=c(rep(2,2),rep(3,2),rep(4,4)) ) ### RNA Plot wiht nucleotides and dots ### RNAPlot(dat) RNAPlot(dat,nt=TRUE,add=TRUE,dp=0.75)
## Create a CT file from bracket notation and sequence ### ct=makeCt( "((((...(((((((....)))))))...((((...))))...))))", "CCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC" ) ## Create a coordinate file based on the CT file ### dat=ct2coord(ct) ### Create a plot of the secondary structure ### RNAPlot(dat) ### Plot positions 1:4 as green and 43:46 circles ## ### and show the legend ranges=data.frame(min=c(1,43),max=c(4,46),col=c(2,3), desc=c("Region 1","Region 2") ) RNAPlot(dat,ranges,labTF=TRUE) ### Highlight the sequences CUCCU and CCCCAAA ### RNAPlot(dat,hl=c("CUCCU","CCCCAAA"),seqcol=c(2,4),labTF=TRUE,main="RNA Molecule") ### Modify specific positions #### RNAPlot( dat, modspec=TRUE, modp=c(1:4,43:46),mod=c(17,17,15,15,16,16,16,16), modcol=c(rep(2,2),rep(3,2),rep(4,4)) ) ### RNA Plot with nucleotides ### RNAPlot(dat,nt=TRUE) ### RNA plot with nucleotides RNAPlot( dat,nt=TRUE,modspec=TRUE,modp=c(1:4,43:46), mod=c(17,17,15,15,16,16,16,16), modcol=c(rep(2,2),rep(3,2),rep(4,4)) ) ### RNA Plot wiht nucleotides and dots ### RNAPlot(dat) RNAPlot(dat,nt=TRUE,add=TRUE,dp=0.75)
Rotates a point a given angle around a given center point.
rotateS(x2, y2, x0, y0, ang)
rotateS(x2, y2, x0, y0, ang)
x2 |
x coordinate of the position being rotated |
y2 |
y coordinate of the position being rotated |
x0 |
x coordinate of the center of rotation |
y0 |
y coordinate of the center of rotation |
ang |
rotatation angle in radians |
Returns a rotated point
JP Bida
### Rotate a point 90 degress ### rotateS(0,1,0,0,pi/2)
### Rotate a point 90 degress ### rotateS(0,1,0,0,pi/2)
Rotates a set of points around a center point a given number of radians
rotateV(x2, y2, x0, y0, ang)
rotateV(x2, y2, x0, y0, ang)
x2 |
Vector containing x coordinates being rotated |
y2 |
Vector containing y coordinates being rotated |
x0 |
x coordinate of center of rotation |
y0 |
y coordinate of center of rotation |
ang |
Angle of rotation given in radians |
set of rotated points
JP Bida
x=c(1,0,-1,0) y=c(0,1,0,-1) pts=rotateV(x,y,0,0,pi/4)
x=c(1,0,-1,0) y=c(0,1,0,-1) pts=rotateV(x,y,0,0,pi/4)
internal function that generates coordinates for an RNA secondary structure stem
stemCords(input, p1, p2, x1, y1, x2, y2, x3, y3)
stemCords(input, p1, p2, x1, y1, x2, y2, x3, y3)
input |
ct file as data frame |
p1 |
index of nucleotide in first base pair of the stem |
p2 |
index of nucleotide in first base pair of the stem |
x1 |
x coordinate of p1 |
y1 |
y coordinate of p1 |
x2 |
x coordinate of p2 |
y2 |
y coordinate of p2 |
x3 |
direction vector x component |
y3 |
direction vector y component |
set of points
This is an internal function not recommend for use out side of the ct2coord function
JP Bida
### Internal Function ###
### Internal Function ###
Given a coordinate file, a point, and an angle in radians transformFold rotates the fold around the given point the given number of radians.
transformFold(dat, x0, y0, ang)
transformFold(dat, x0, y0, ang)
dat |
Coordinate file containing multiple RNA folds |
x0 |
x coordinate of center of rotation |
y0 |
y coordinate of center of rotation |
ang |
angle of rotation in radians |
dat frame containing the rotated coordinates
JP Bida
ct=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC") c1=ct2coord(ct) RNAPlot(c1) c2=transformFold(c1,0,0,pi/2) c3=transformFold(c2,0,0,pi/2) c4=transformFold(c3,0,0,pi/2) RNAPlot(c2,add=TRUE) RNAPlot(c3,add=TRUE) RNAPlot(c4,add=TRUE)
ct=makeCt("((((...(((((((.........)))))))...((((.........))))...))))", "AAAAAAAACCCCCCCCAAAGGGGGGGAUUACCCCUCCUUUAAAAGGGUUUUCCCCCCC") c1=ct2coord(ct) RNAPlot(c1) c2=transformFold(c1,0,0,pi/2) c3=transformFold(c2,0,0,pi/2) c4=transformFold(c3,0,0,pi/2) RNAPlot(c2,add=TRUE) RNAPlot(c3,add=TRUE) RNAPlot(c4,add=TRUE)
internal function to translate points
translate(x1, y1, x2, y2)
translate(x1, y1, x2, y2)
x1 |
x coordinates being translated |
y1 |
y coordinates being translated |
x2 |
dx for translation |
y2 |
dy for translation |
set of points
JP Bida
## Internal Function ##
## Internal Function ##
Given a CT file checks for correct column names, all nucleotides have positions greater than 0, unique rows per nucleotide, multiple pairings, duplicate pairings so that the CT file can be used to generate a secondary structure using bracket notation.
validCT(ct)
validCT(ct)
ct |
Dataframe generated from loadCT |
Returns messages indicating issues with the CT file that would prevent it from generating a secondary structure in bracket notation
JP Bida