plot.PHMM {aphid}R Documentation

Plot profile hidden Markov models.

Description

plot.PHMM provides a visual representation of a profile hidden Markov model.

Usage

## S3 method for class 'PHMM'
plot(
  x,
  from = "start",
  to = "end",
  just = "center",
  arrexp = 1,
  textexp = 1,
  ...
)

Arguments

x

an object of class "PHMM".

from

an integer giving the module number to start the plot sequence from. Also accepts the chracter string "start" (module 0; default).

to

an integer giving the module number to terminate the plot sequence. Also accepts the chracter string "end" (default).

just

a character string giving the justfication of the plot relative to the device. Accepted values are "left", "center" and "right".

arrexp

the expansion factor to be applied to the arrows in the plot.

textexp

the expansion factor to be applied to the text in the plot.

...

additional arguments to be passed to plot.

Details

"plot.PHMM" Plots a "PHMM" object as a directed graph with sequential modules consisting of squares, diamonds and circles representing match, insert and delete states, respectively. Modules are interconnected by directed lines with line-weights proportional to the transition probabilities between the states. Since the plotted models are generally much longer than they are high, it is usually better to output the plot to a PDF file as demonstrated in the example below.

Value

NULL (invisibly).

Author(s)

Shaun Wilkinson

References

Durbin R, Eddy SR, Krogh A, Mitchison G (1998) Biological sequence analysis: probabilistic models of proteins and nucleic acids. Cambridge University Press, Cambridge, United Kingdom.

See Also

plot.HMM

Examples

  ## Small globin alignment example from Durbin et al (1998) Figure 5.3
  data(globins)
  ## derive a profile hidden Markov model from the alignment
  globins.PHMM <- derivePHMM(globins, residues = "AMINO", seqweights = NULL)
  ## plot the PHMM
  plot(globins.PHMM, main = "Profile hidden Markov model for globins")
  ##
  ## derive a profile hidden Markov model from the woodmouse dataset in the
  ## ape package
  library(ape)
  data(woodmouse)
  woodmouse.PHMM <- derivePHMM(woodmouse)
  ## plot partial model to viewer device
  plot(woodmouse.PHMM, from = 0, to = 5)
  ## plot the entire model to a PDF in the current working directory
  
  tmpf <- tempfile(fileext = ".pdf")
  nr <- ceiling((woodmouse.PHMM$size + 2)/10)
  pdf(file = tmpf, width = 8.27, height = nr * 2)
  par(mfrow = c(nr, 1), mar = c(0, 0, 0, 0) + 0.1)
  from <- 0
  to <- 10
  for(i in 1:nr){
    plot(woodmouse.PHMM, from = from, to = to, just = "left")
    from <- from + 10
    to <- min(to + 10, woodmouse.PHMM$size + 1)
  }
  dev.off()
  

[Package aphid version 1.3.5 Index]