R: Methods implementing Lattice displays for flow data
levelplot,formula,flowSet-method
R Documentation
Methods implementing Lattice displays for flow data
Description
Various methods implementing multipanel visualizations for flow data using
infrastructure provided in the lattice package. The original generics for
these methods are defined in lattice, and these S4 methods (mostly) dispatch
on a formula and the data argument which must be of class
flowSet or flowFrame. The formula has to be fairly basic:
conditioning can be done using phenodata variables and channel names (the
colnames slot) can be used as panel variables. See examples below for
sample usage.
a formula describing the structure of the plot and the variables to
be used in the display.
data
a flowSet object that serves as a source of data.
xlab,ylab
Labels for data axes, with suitable defaults taken from the
formula
as.table,contour,labels
These arguments are passed unchanged to the
corresponding methods in lattice, and are listed here only because they
provide different defaults. See documentation for the original methods for
details.
n
the number of bins on each axis to be used when evaluating the
density
f.value,distribution
number of points used in Q-Q plot, and the
reference distribution used. See qqmath for
details.
reorder.by
a function, which is applied to each column. The columns
are ordered by the results. Reordering can be suppressed by setting this to
NULL.
time
A character string giving the name of the column recording time.
exclude.time
logical, specifying whether to exclude the time variable
from a scatter plot matrix or parallel coordinates plot. It is rarely
meaningful not to do so.
filter
flowCore filter
...
more arguments, usually passed on to the underlying lattice
methods.
Details
Not all standard lattice arguments will have the intended effect, but many
should. For a fuller description of possible arguments and their effects,
consult documentation on lattice (Trellis docs would also work for the
fundamentals).
Methods
qqmath
signature(x = "formula", data = "flowSet"): creates
theoretical quantile plots of a given channel, with one or more samples per
panel
levelplot
signature(x = "formula", data = "flowSet"): similar
to the xyplot method, but plots estimated density (using
kde2d) with a common z-scale and an optional color
key.
parallel
signature(x = "flowFrame", data = "missing"): draws a
parallel coordinates plot of all channels (excluding time, by default) of a
flowFrame object. This is rarely useful without transparency, but
that is currently only possible with the pdf device (and
perhaps the aqua device as well).
Examples
data(GvHD)
qqmath( ~ `FSC-H` | factor(Patient), GvHD,
grid = TRUE, type = "l",
f.value = ppoints(100))
## contourplot of bivariate density:
require(colorspace)
YlOrBr <- c("#FFFFD4", "#FED98E", "#FE9929", "#D95F0E", "#993404")
colori <- colorRampPalette(YlOrBr)
levelplot(asinh(`SSC-H`) ~ asinh(`FSC-H`) | Visit + Patient, GvHD, n = 20,
col.regions = colori(50), main = "Contour Plot")
## parallel coordinate plots
parallel(GvHD[["s6a01"]])
## Not run:
## try with PDF device
parallel(GvHD[["s7a01"]], alpha = 0.01)
## End(Not run)
Results
R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)
R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.
R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.
Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.
> library(flowViz)
Loading required package: flowCore
Loading required package: lattice
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/flowViz/lattice-methods.Rd_%03d_medium.png", width=480, height=480)
> ### Name: levelplot,formula,flowSet-method
> ### Title: Methods implementing Lattice displays for flow data
> ### Aliases: lattice-methods levelplot,formula,flowSet-method
> ### parallel,flowFrame,missing-method parallel,formula,flowSet-method
> ### qqmath,formula,flowSet-method
> ### Keywords: dplot methods
>
> ### ** Examples
>
> data(GvHD)
>
> qqmath( ~ `FSC-H` | factor(Patient), GvHD,
+ grid = TRUE, type = "l",
+ f.value = ppoints(100))
>
>
> ## contourplot of bivariate density:
>
> require(colorspace)
Loading required package: colorspace
> YlOrBr <- c("#FFFFD4", "#FED98E", "#FE9929", "#D95F0E", "#993404")
> colori <- colorRampPalette(YlOrBr)
> levelplot(asinh(`SSC-H`) ~ asinh(`FSC-H`) | Visit + Patient, GvHD, n = 20,
+ col.regions = colori(50), main = "Contour Plot")
>
>
>
> ## parallel coordinate plots
>
> parallel(GvHD[["s6a01"]])
>
> ## Not run:
> ##D
> ##D ## try with PDF device
> ##D parallel(GvHD[["s7a01"]], alpha = 0.01)
> ##D
> ## End(Not run)
>
>
>
>
>
> dev.off()
null device
1
>