kooperberg               package:limma               R Documentation

_K_o_o_p_e_r_b_e_r_g _M_o_d_e_l-_B_a_s_e_d _B_a_c_k_g_r_o_u_n_d _C_o_r_r_e_c_t_i_o_n

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

     This function uses a Bayesian model to  background correct GenePix
     microarray data.

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

     kooperberg(RG, a=TRUE, layout=RG$printer, verbose=TRUE)

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

      RG: an RGList of GenePix data, read in using 'read.maimages',
          with 'other.columns=c("F635 SD","B635 SD","F532 SD","B532
          SD","B532 Mean","B635 Mean","F Pixels","B Pixels")'.

       a: logical.  If 'TRUE', the 'a' parameters in the model
          (equation 3 and 4) are estimated for each slide.  If 'FALSE'
          the 'a' parameters are set to unity.

  layout: list containing print layout with components 'ngrid.r',
          'ngrid.c', 'nspot.r' and 'nspot.c'.  Defaults to
          'RG$printer'.

 verbose: logical.  If 'TRUE', progress is reported to standard output.

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

     This function is for use with GenePix data and is designed to cope
     with the problem of large numbers of negative intensities and
     hence missing values on the log-intensity scale. It avoids missing
     values in most cases and at the same time dampens down the
     variability of log-ratios for low intensity spots. See Kooperberg
     et al (2002) for more details.

     'kooperberg' uses the foreground and background intensities,
     standard deviations and number of pixels to compute empirical
     estimates of the model  parameters as described in equation 2 of
     Kooperberg et al (2002).

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

     An 'RGList' containing the components 

       R: matrix containing the background adjusted intensities for the
          red channel for each spot for each array

       G: matrix containing the background adjusted intensities for the
          green channel for each spot for each array

 printer: list containing print layout

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

     Matthew Ritchie

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

     Kooperberg, C., Fazzio, T. G., Delrow, J. J., and Tsukiyama, T.
     (2002) Improved background correction for spotted DNA microarrays.
     _Journal of Computational Biology_ *9*, 55-66.

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

     04.Background gives an overview of background correction functions
     defined in the LIMMA package.

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

     #  This is example code for reading and background correcting GenePix data
     #  given GenePix Results (gpr) files in the working directory (data not
     #  provided).
     ## Not run: 
     genepixFiles <- dir(pattern="*\\.gpr$") # get the names of the GenePix image analysis output files in the current directory
     RG <- read.maimages(genepixFiles, source="genepix", other.columns=c("F635 SD","B635 SD","F532 SD","B532 SD","B532 Mean","B635 Mean","F Pixels","B Pixels"))
     RGmodel <- kooperberg(RG)
     MA <- normalizeWithinArrays(RGmodel)
     ## End(Not run)

