printtipWeights            package:limma            R Documentation

_S_u_b-_a_r_r_a_y _Q_u_a_l_i_t_y _W_e_i_g_h_t_s

_D_e_s_c_r_i_p_t_i_o_n:

     Estimates relative quality weights for each sub-array in a
     multi-array experiment.

_U_s_a_g_e:

     printtipWeights(object, design = NULL, weights = NULL, method = "genebygene", layout, maxiter = 50, tol = 1e-10, trace=FALSE)

_A_r_g_u_m_e_n_t_s:

  object: object of class 'numeric', 'matrix', 'MAList', 'marrayNorm',
          or 'exprSet' containing log-ratios or log-values of
          expression for a series of spotted microarrays.

  design: the design matrix of the microarray experiment, with rows
          corresponding to arrays and columns to coefficients to be
          estimated.  Defaults to the unit vector meaning that the
          arrays are treated as replicates.

 weights: optional numeric matrix containing prior weights for each
          spot.

  method: character string specifying the estimating algorithm to be
          used. Choices are '"genebygene"' and '"reml"'.

  layout: list specifying the dimensions of the spot matrix and the
          grid matrix. For details see 'PrintLayout-class'.

 maxiter: maximum number of iterations allowed.

     tol: convergence tolerance.

   trace: logical variable. If true then output diagnostic information
          at each iteration of '"reml"' algorithm.

_D_e_t_a_i_l_s:

     The relative reliability of each sub-array (print-tip group) is
     estimated by measuring how well the expression values for that
     sub-array follow the linear model.

     The method described in Ritchie et al (2006) and implemented in 
     the 'arrayWeights' function is adapted for this purpose.A
     heteroscedastic model is fitted to the expression values for  each
     gene by calling the function 'lm.wfit'.  The dispersion model  is
     fitted to the squared residuals from the mean fit, and is set up
     to  have sub-array specific coefficients, which are updated in
     either full REML  scoring iterations, or using an efficient
     gene-by-gene update algorithm.   The final estimates of the
     sub-array variances are converted to weights.

     The data object 'object' is interpreted as for 'lmFit'. In
     particular, the arguments 'design', 'weights' and 'layout' will 
     be extracted from the data 'object' if available and do not
     normally need to  be set explicitly in the call; if any of these
     are set in the call then they  will over-ride the slots or
     components in the data 'object'.

_V_a_l_u_e:

     A matrix of sub-array weights which can be passed to 'lmFit'.

_A_u_t_h_o_r(_s):

     Matthew Ritchie and Gordon Smyth

_R_e_f_e_r_e_n_c_e_s:

     Ritchie, M. E., Diyagama, D., Neilson, van Laar, R., J., Dobrovic,
     A., Holloway, A., and Smyth, G. K. (2006). Empirical array quality
     weights in the analysis of microarray data. BMC Bioinformatics 7,
     261. <URL: http://www.biomedcentral.com/1471-2105/7/261/abstract>

_S_e_e _A_l_s_o:

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

_E_x_a_m_p_l_e_s:

     library(sma)
     # Subset of data from ApoAI case study in Limma User's Guide
     data(MouseArray)
     # Avoid non-positive intensities
     RG <- backgroundCorrect(mouse.data, method="half")
     MA <- normalizeWithinArrays(RG, mouse.setup)
     MA <- normalizeBetweenArrays(MA, method="Aq")
     targets <- data.frame(Cy3=I(rep("Pool",6)),Cy5=I(c("WT","WT","WT","KO","KO","KO")))
     design <- modelMatrix(targets, ref="Pool")
     subarrayw <- printtipWeights(MA, design, layout=mouse.setup)
     fit <- lmFit(MA, design, weights=subarrayw)
     fit2 <- contrasts.fit(fit, contrasts=c(-1,1))
     fit2 <- eBayes(fit2)
     # Use of sub-array weights increases the significance of the top genes
     topTable(fit2)
     # Create an image plot of sub-array weights from each array
     zlim <- c(min(subarrayw), max(subarrayw))
     par(mfrow=c(3,2), mai=c(0.1,0.1,0.3,0.1))
     for(i in 1:6) 
             imageplot(subarrayw[,i], layout=mouse.setup, zlim=zlim, main=paste("Array", i))

