Correlation Between Duplicates
Estimate the correlation between duplicate spots (regularly spaced replicate spots on the same array) or between technical replicates from a series of arrays.
duplicateCorrelation(object, design=NULL, ndups=2, spacing=1, block=NULL, trim=0.15, weights=NULL)
object |
a numeric matrix of expression values, or any data object from which |
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 |
ndups |
a positive integer giving the number of times each gene is printed on an array. |
spacing |
the spacing between the rows of |
block |
vector or factor specifying a blocking variable |
trim |
the fraction of observations to be trimmed from each end of |
weights |
an optional numeric matrix of the same dimension as |
When block=NULL
, this function estimates the correlation between duplicate spots (regularly spaced within-array replicate spots).
If block
is not null, this function estimates the correlation between repeated observations on the blocking variable.
Typically the blocks are biological replicates and the repeated observations are technical replicates.
In either case, the correlation is estimated by fitting a mixed linear model by REML individually for each gene.
The function also returns a consensus correlation, which is a robust average of the individual correlations, which can be used as input for
functions lmFit
or gls.series
.
At this time it is not possible to estimate correlations between duplicate spots and between technical replicates simultaneously.
If block
is not null, then the function will set ndups=1
, which is equivalent to ignoring duplicate spots.
For this function to return statistically useful results, there must be at least two more arrays than the number of coefficients to be estimated, i.e., two more than the column rank of design
.
The function may take long time to execute as it fits a mixed linear model for each gene for an iterative algorithm. It is not uncommon for the function to return a small number of warning messages that correlation estimates cannot be computed for some individual genes. This is not a serious concern providing that there are only a few such warnings and the total number of genes is large. The consensus estimator computed by this function will not be materially affected by a small number of genes.
A list with components
consensus.correlation |
the average estimated inter-duplicate correlation. The average is the trimmed mean of the individual correlations on the atanh-transformed scale. |
cor |
same as |
atanh.correlations |
numeric vector of length |
Gordon Smyth
Smyth, G. K., Michaud, J., and Scott, H. (2005). The use of within-array replicate spots for assessing differential expression in microarray experiments. Bioinformatics 21(9), 2067-2075. [http://bioinformatics.oxfordjournals.org/content/21/9/2067] [Preprint with corrections: http://www.statsci.org/smyth/pubs/dupcor.pdf]
These functions use mixedModel2Fit
from the statmod package.
An overview of linear model functions in limma is given by 06.LinearModels.
# Simulate gene expression data for 100 probes and 6 microarrays # Microarray are in two groups # First two probes are more highly expressed in second group # Std deviations vary between genes with prior df=4 sd <- 0.3*sqrt(4/rchisq(100,df=4)) y <- matrix(rnorm(100*6,sd=sd),100,6) rownames(y) <- paste("Gene",1:100) y[1:2,4:6] <- y[1:2,4:6] + 2 design <- cbind(Grp1=1,Grp2vs1=c(0,0,0,1,1,1)) options(digits=3) # Fit with correlated arrays # Suppose each pair of arrays is a block block <- c(1,1,2,2,3,3) dupcor <- duplicateCorrelation(y,design,block=block) dupcor$consensus.correlation fit1 <- lmFit(y,design,block=block,correlation=dupcor$consensus) fit1 <- eBayes(fit1) topTable(fit1,coef=2) # Fit with duplicate probes # Suppose two side-by-side duplicates of each gene rownames(y) <- paste("Gene",rep(1:50,each=2)) dupcor <- duplicateCorrelation(y,design,ndups=2) dupcor$consensus.correlation fit2 <- lmFit(y,design,ndups=2,correlation=dupcor$consensus) dim(fit2) fit2 <- eBayes(fit2) topTable(fit2,coef=2)
Please choose more modern alternatives, such as Google Chrome or Mozilla Firefox.