rma {xps}R Documentation

Robust Multi-Array Average Expression Measure

Description

This function converts a DataTreeSet into an ExprTreeSet using the robust multi-array average (RMA) expression measure.

Usage

rma(xps.data,
    filename   = character(0),
    filedir    = getwd(),
    tmpdir     = "",
    background = "pmonly",
    normalize  = TRUE,
    option     = "transcript",
    exonlevel  = "",
    xps.scheme = NULL,
    add.data   = TRUE,
    verbose    = TRUE)

xpsRMA(object, ...)

Arguments

xps.data object of class DataTreeSet.
filename file name of ROOT data file.
filedir system directory where ROOT data file should be stored.
tmpdir optional temporary directory where temporary ROOT files should be stored.
background probes used to compute background, one of ‘pmonly’, ‘mmonly’, ‘both’; for genome/exon arrays one of ‘genomic’, ‘antigenomic’
normalize logical. If TRUE normalize data using quantile normalization.
option option determining the grouping of probes for summarization, one of ‘transcript’, ‘exon’, ‘probeset’; exon arrays only.
exonlevel exon annotation level determining which probes should be used for summarization; exon/genome arrays only.
xps.scheme optional alternative SchemeTreeSet.
add.data logical. If TRUE expression data will be included as slot data.
verbose logical, if TRUE print status information.
object object of class DataTreeSet.
... the arguments described above.

Details

This function computes the RMA (Robust Multichip Average) expression measure described in Irizarry et al. for both expression arrays and exon arrays. For exon arrays it is necessary to supply the requested option and exonlevel.

Following options are valid for exon arrays:
transcript: expression levels are computed for transcript clusters, i.e. probe sets containing the same ‘transcript_cluster_id’.
exon: expression levels are computed for exon clusters, i.e. probe sets containing the same ‘exon_id’, where each exon cluster consists of one or more probesets.
probeset: expression levels are computed for individual probe sets, i.e. for each ‘probeset_id’.

Following exonlevel annotations are valid for exon arrays:
core: probesets supported by RefSeq and full-length GenBank transcripts.
metacore: core meta-probesets.
extended: probesets with other cDNA support.
metaextended: extended meta-probesets.
full: probesets supported by gene predictions only.
metafull: full meta-probesets.
affx: standard AFFX controls.
all: combination of above (including affx).

Following exonlevel annotations are valid for whole genome arrays:
core: probesets with category ‘unique’, ‘similar’ and ‘mixed’.
metacore: probesets with category ‘unique’ only.
affx: standard AFFX controls.
all: combination of above (including affx).

Exon levels can also be combined, with following combinations being most useful:
exonlevel="metacore+affx": core meta-probesets plus AFFX controls
exonlevel="core+extended": probesets with cDNA support
exonlevel="core+extended+full": supported plus predicted probesets

Exon level annotations are described in the Affymetrix whitepaper exon_probeset_trans_clust_whitepaper.pdf:
“Exon Probeset Annotations and Transcript Cluster Groupings”.

In order to use an alternative SchemeTreeSet set the corresponding SchemeSet xps.scheme.

xpsRMA is the DataSet method called by function rma, containing the same parameters.

Value

An ExprTreeSet

Note

In contrary to other implementations of RMA the expression measure is given to you in linear scale, analogously to the expression measures computed with mas5 and mas4.

It is also possible to skip background correction by setting parameter background="none".

For the analysis of many exon arrays it may be better to define a tmpdir, since this will store only the results in the main file and not e.g. background and normalized intensities, and thus will reduce the file size of the main file. For quantile normalization memory should not be an issue, however medianpolish depends on RAM unless you are using a temporary file.

Parameter exonlevel determines not only which probes are used for medianpolish, but also the probes used for background calculation and for quantile normalization. If you want to use seperate probes for background calculation, quantile normalization and medianpolish summarization, you can pass a numeric vector containing three integer values corresponding to the respective exonlevel, e.g. you can use exonlevel=c(16316,8252,8252), see function exonLevel for more details.

Author(s)

Christian Stratowa

References

Rafael. A. Irizarry, Benjamin M. Bolstad, Francois Collin, Leslie M. Cope, Bridget Hobbs and Terence P. Speed (2003), Summaries of Affymetrix GeneChip probe level data Nucleic Acids Research 31(4):e15

Bolstad, B.M., Irizarry R. A., Astrand M., and Speed, T.P. (2003), A Comparison of Normalization Methods for High Density Oligonucleotide Array Data Based on Bias and Variance. Bioinformatics 19(2):185-193

Irizarry, RA, Hobbs, B, Collin, F, Beazer-Barclay, YD, Antonellis, KJ, Scherf, U, Speed, TP (2003) Exploration, Normalization, and Summaries of High Density Oligonucleotide Array Probe Level Data. Biostatistics .Vol. 4, Number 2: 249-264

See Also

express

Examples

## first, load ROOT scheme file and ROOT data file
scheme.test3 <- root.scheme(paste(.path.package("xps"),"schemes/SchemeTest3.root",sep="/"))
data.test3 <- root.data(scheme.test3, paste(.path.package("xps"),"rootdata/DataTest3_cel.root",sep="/"))

data.rma <- rma(data.test3,"tmp_Test3RMA",tmpdir="",background="pmonly",normalize=TRUE,verbose=FALSE)

## get data.frame
expr.rma <- validData(data.rma)
head(expr.rma)

## plot results
if (interactive()) {
boxplot(data.rma)
boxplot(log2(expr.rma))
}

rm(scheme.test3, data.test3)
gc()

## Not run: 
## examples using Affymetrix human tissue dataset (see also xps/examples/script4exon.R)
## first, load ROOT scheme file and ROOT data file from e.g.:
scmdir <- "/Volumes/GigaDrive/CRAN/Workspaces/Schemes"
datdir <- "/Volumes/GigaDrive/CRAN/Workspaces/ROOTData"

## 1. example - expression array, e.g. HG-U133_Plus_2:
scheme.u133p2 <- root.scheme(paste(scmdir,"Scheme_HGU133p2_na25.root",sep="/"))
data.u133p2   <- root.data(scheme.u133p2, paste(datdir,"HuTissuesU133P2_cel.root",sep="/"))

workdir <- "/Volumes/GigaDrive/CRAN/Workspaces/Exon/hutissues/u133p2"
data.rma <- rma(data.u133p2,"MixU133P2RMA",filedir=workdir,tmpdir="",
                background="pmonly",normalize=TRUE)

## 2. example - whole genome array, e.g. HuGene-1_0-st-v1:
scheme.genome <- root.scheme(paste(scmdir,"Scheme_HuGene10stv1r3_na25.root",sep="/"))
data.genome   <- root.data(scheme.genome, paste(datdir,"HuTissuesGenome_cel.root",sep="/"))

workdir <- "/Volumes/GigaDrive/CRAN/Workspaces/Exon/hutissues/hugene"
data.g.rma <- rma(data.genome,"HuGeneMixRMAMetacore",filedir=workdir,tmpdir="",
                  background="antigenomic",normalize=T,exonlevel="metacore+affx")

## 3. example - exon array, e.g. HuEx-1_0-st-v2:
scheme.exon <- root.scheme(paste(scmdir,"Scheme_HuEx10stv2r2_na25.root",sep="/"))
data.exon   <- root.data(scheme.exon, paste(datdir,"HuTissuesExon_cel.root",sep="/"))

workdir <- "/Volumes/GigaDrive/CRAN/Workspaces/Exon/hutissues/exon"
data.x.rma <- rma(data.exon,"MixRMAMetacore",filedir=workdir,tmpdir="",background="antigenomic",
                  normalize=T,option="transcript",exonlevel="metacore")
## End(Not run)

[Package xps version 1.2.10 Index]