Last data update: 2014.03.03

R: PLOT CONTROL WELLS IN RTCA DATA
controlViewR Documentation

PLOT CONTROL WELLS IN RTCA DATA

Description

A convenience function to plot sample wells with control wells on an E-plate in RTCA system. To use the function the phenoData field of the RTCA object must contain a field named “GeneSymbol”.

Usage

controlView(rtca, genesymbol = c("Allstar", "COPB2", "GFP", "mock", "PLK1", "WEE1"), cols, ylim, smooth = FALSE, group = TRUE, ylab = "Normalized cell index", xlab = "Time interval (hour)", drawsd = TRUE, normline = TRUE, ncol = 1, legendpos = "topleft", pData.column="GeneSymbol",...)

Arguments

rtca

An object of RTCA.To use the function, the phenoData must contain a column which name is specified by the pData.column parameter.

genesymbol

character, gene symbols to be plotted.

cols

character, colors used by the provided gene symbols

ylim

y-axis lim

smooth

logical, whether the RTCA object should be smoothed before plotting

group

logical. If ‘group’ is set to TRUE, wells with the same GeneSymbol will be summarized and plotted. For instance, these could be biological replicates. Otherwise each well is plotted separatedly

ylab

y axis label

xlab

x axis label

drawsd

logical, should the error bar be drawn to represent standard deviation?

normline

logical, should the base-time indicated by a line? See ratioTransform for the concept of the base-time

ncol

integer, legend column number

legendpos

character, legend position

pData.column

The column which the genesymbol parameter will be matched with

...

other parameters passed to the plot function

Details

The function is often called to draw sample and control in one plot.

Value

NULL, the function is called for its side effect

Author(s)

Jitao David Zhang jitao_david.zhang@roche.com

See Also

RTCA

Examples

require(RTCA)
  
ofile <- system.file("extdata/testOutput.csv", package="RTCA")
pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA")

pData <- read.csv(pfile, sep="\t", row.names="Well")
metaData <- data.frame(labelDescription=c(
"Rack number",
"siRNA catalogue number",
"siRNA gene symbol",
"siRNA EntrezGene ID",
"siRNA targeting accession"
))

phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData)
x <- parseRTCA(ofile, phenoData=phData)

controlView(x, genesymbol=c("mock","COPB2","PLK1"),ylim=c(0,2))

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(RTCA)
Loading required package: Biobase
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: RColorBrewer
Loading required package: gtools
> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/RTCA/controlView.Rd_%03d_medium.png", width=480, height=480)
> ### Name: controlView
> ### Title: PLOT CONTROL WELLS IN RTCA DATA
> ### Aliases: controlView
> ### Keywords: hplot
> 
> ### ** Examples
> 
> require(RTCA)
>   
> ofile <- system.file("extdata/testOutput.csv", package="RTCA")
> pfile <- system.file("extdata/testOutputPhenoData.csv", package="RTCA")
> 
> pData <- read.csv(pfile, sep="\t", row.names="Well")
> metaData <- data.frame(labelDescription=c(
+ "Rack number",
+ "siRNA catalogue number",
+ "siRNA gene symbol",
+ "siRNA EntrezGene ID",
+ "siRNA targeting accession"
+ ))
> 
> phData <- new("AnnotatedDataFrame", data=pData, varMetadata=metaData)
> x <- parseRTCA(ofile, phenoData=phData)
Read 245 items
> 
> controlView(x, genesymbol=c("mock","COPB2","PLK1"),ylim=c(0,2))
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>