motifCircos: plot sequence logo stacks with a radial phylogenic tree and...

Description Usage Arguments Value Author(s) See Also Examples

Description

plot sequence logo stacks with a radial phylogenic tree and multiple color rings. The difference from plotMotifStackWithRadialPhylog is that it has more color setting and one more group of pfms.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
motifCircos(
  phylog,
  pfms = NULL,
  pfms2 = NULL,
  R = 2.5,
  r.tree = 1,
  col.tree.bg = NULL,
  col.tree.bg.alpha = 1,
  cnodes = 0,
  labels.nodes = names(phylog$nodes),
  clabel.nodes = 0,
  r.leaves = NA,
  cleaves = 1,
  labels.leaves = names(phylog$leaves),
  clabel.leaves = 1,
  col.leaves = rep("black", length(labels.leaves)),
  col.leaves.bg = NULL,
  col.leaves.bg.alpha = 1,
  r.pfms = NA,
  r.pfms2 = NA,
  r.rings = 0,
  col.rings = list(),
  col.inner.label.circle = NULL,
  inner.label.circle.width = 0.02,
  col.outer.label.circle = NULL,
  outer.label.circle.width = 0.02,
  draw.box = FALSE,
  clockwise = FALSE,
  init.angle = if (clockwise) 90 else 0,
  angle = 360,
  pfmNameSpliter = ";",
  rcpostfix = "(RC)",
  motifScale = c("linear", "logarithmic", "none"),
  ic.scale = TRUE,
  plotIndex = FALSE,
  IndexCol = "black",
  IndexCex = 0.8,
  groupDistance = NA,
  groupDistanceLineCol = "red",
  plotAxis = FALSE
)

Arguments

phylog

an object of class phylog

pfms

a list of objects of class pfm

pfms2

a list of objects of class pfm

R

radius of canvas

r.tree

half width of the tree

col.tree.bg

a vector of colors for tree background

col.tree.bg.alpha

a alpha value [0, 1] of colors for tree background

cnodes

a character size for plotting the points that represent the nodes, used with par("cex")*cnodes. If zero, no points are drawn

labels.nodes

a vector of strings of characters for the nodes labels

clabel.nodes

a character size for the nodes labels, used with par("cex")*clabel.nodes. If zero, no nodes labels are drawn

r.leaves

width of the leaves

cleaves

a character size for plotting the points that represent the leaves, used with par("cex")*cleaves. If zero, no points are drawn

labels.leaves

a vector of strings of characters for the leaves labels

clabel.leaves

a character size for the leaves labels, used with par("cex")*clavel.leaves

col.leaves

a vector of colors for leaves labels

col.leaves.bg

a vector of colors for background of leaves labels

col.leaves.bg.alpha

alpha value [0, 1] for the colors of backgroud of leaves labels

r.pfms

width of the pfms

r.pfms2

width of the pfms2

r.rings

a vector of width of color rings

col.rings

a list of color rings

col.inner.label.circle

a vector of colors for inner cirlce of pfms

inner.label.circle.width

width for inner circle of pfms

col.outer.label.circle

a vector of colors for outer circle of pfms

outer.label.circle.width

width for outer circle of pfms

draw.box

if TRUE draws a box around the current plot with the function box()

clockwise

a logical value indicating if slices are drawn clockwise or counter clockwise

init.angle

number specifying the starting angle (in degrees) for the slices. Defaults to 0 (i.e., ‘3 o’clock') unless clockwise is true where init.angle defaults to 90 (degrees), (i.e., ‘12 o’clock')

angle

number specifying the angle (in degrees) for phylogenic tree. Defaults 360

pfmNameSpliter

spliter when name of pfms/pfms2 contain multiple node of labels.leaves

rcpostfix

the postfix for reverse complements

motifScale

the scale of logo size

ic.scale

logical. If TRUE, the height of each column is proportional to its information content. Otherwise, all columns have the same height.

plotIndex

logical. If TRUE, will plot index number in the motifLogo which can help user to describe the motifLogo

IndexCol

The color of the index number when plotIndex is TRUE.

IndexCex

The cex of the index number when plotIndex is TRUE.

groupDistance

show groupDistance on the draw

groupDistanceLineCol

groupDistance line color, default: red

plotAxis

logical. If TRUE, will plot distance axis.

Value

none

Author(s)

Jianhong Ou

See Also

plotMotifStackWithRadialPhylog

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
if(interactive() || Sys.getenv("USER")=="jianhongou"){
    library("MotifDb")
    matrix.fly <- query(MotifDb, "Dmelanogaster")
    motifs <- as.list(matrix.fly)
    motifs <- motifs[grepl("Dmelanogaster-FlyFactorSurvey-", 
                            names(motifs), fixed=TRUE)]
    names(motifs) <- gsub("Dmelanogaster_FlyFactorSurvey_", "", 
                gsub("_FBgn[0-9]+$", "", 
                  gsub("[^a-zA-Z0-9]","_", 
                     gsub("(_[0-9]+)+$", "", names(motifs)))))
    motifs <- motifs[unique(names(motifs))]
    pfms <- sample(motifs, 50)
    hc <- clusterMotifs(pfms)
    library(ade4)
    phylog <- ade4::hclust2phylog(hc)
    leaves <- names(phylog$leaves)
    pfms <- pfms[leaves]
    pfms <- mapply(pfms, names(pfms), FUN=function(.ele, .name){
                 new("pfm",mat=.ele, name=.name)})
    pfms <- DNAmotifAlignment(pfms, minimalConsensus=3)
    library(RColorBrewer)
    color <- brewer.pal(12, "Set3")
    motifCircos(phylog, pfms, cleaves = 0.5, clabel.leaves = 0.7, 
                     col.tree.bg=rep(color, each=5), 
                     col.leaves=rep(color, each=5),
                      r.rings=c(0.02, 0.03, 0.04), 
                      col.rings=list(sample(colors(), 50), 
                                     sample(colors(), 50), 
                                     sample(colors(), 50)))
  }

motifStack documentation built on Nov. 8, 2020, 7:43 p.m.