rlm.series {limma}R Documentation

Robust Linear Model for Series of Microarrays

Description

Fit linear models for each gene to a series of microarrays. Fit is by robust M-estimation.

Usage

rlm.series(M,design=NULL,ndups=1,spacing=1,weights=NULL,...)

Arguments

M a numeric matrix. Usually the log-ratios of expression for a series of cDNA microarrrays with rows corresponding to genes and columns to arrays.
design the design matrix of the microarray experiment, with rows corresponding to arrays and columns to comparisons to be estimated. The number of rows must match the number of columns of M. Defaults to the unit vector meaning that the arrays are treated as replicates.
ndups a positive integer giving the number of times each gene is printed on an array. nrow(M) must be divisible by ndups.
spacing the spacing between the rows of M corresponding to duplicate spots, spacing=1 for consecutive spots
weights an optional numeric matrix of the same dimension as M containing weights for each spot. If it is of different dimension to M, it will be filled out to the same size.
... any other arguments are passed to rlm.default.

Details

The linear model is fit for each gene by calling the function rlm from the MASS library.

Warning: don't use weights with this function unless you understand how rlm treats weights. The treatment of weights is different from that of lm.series and gls.series.

Value

A list with components

coefficients numeric matrix containing the estimated coefficients for each linear model. Same number of rows as M, same number of columns as design.
stdev.unscaled numeric matrix conformal with coef containing the unscaled standard deviations for the coefficient estimators. The standard errors are given by stdev.unscaled * sigma.
sigma numeric vector containing the residual standard deviation for each gene.
df.residual numeric vector giving the degrees of freedom corresponding to sigma.

Author(s)

Gordon Smyth

See Also

rlm.

An overview of linear model functions in limma is given by 5.LinearModels.

Examples

#  Simulate gene expression data,
#  6 microarrays and 100 genes with one gene differentially expressed in first 3 arrays
M <- matrix(rnorm(100*6,sd=0.3),100,6)
M[1,1:3] <- M[1,1:3] + 2
#  Design matrix includes two treatments, one for first 3 and one for last 3 arrays
design <- cbind(First3Arrays=c(1,1,1,0,0,0),Last3Arrays=c(0,0,0,1,1,1))
fit <- rlm.series(M,design=design)
eb <- ebayes(fit)
#  Large values of eb$t indicate differential expression
qqt(eb$t[,1],df=fit$df+eb$df.prior)
abline(0,1)

[Package limma version 1.6.7 Index]