Creating Custom CDF's for Affymetrix Chips in Bioconductor

What?

For those who don’t know, CDF files are chip definition format files that define which probes belong to which probesets, and are necessary to use any of the standard summarization methods such as RMA, and others.

Why?

Because we can, and because custom definitions have been shown to be quite useful. See the information over at Brainarray.

Why not somewhere else?

A lot of times other people create custom CDF files based on their own criteria, and make it subsequently available for others to use (see the Brainarray for an example of what some are doing, as well as PlandbAffy)

You have a really nifty idea for a way to reorganize the probesets on an Affymetrix chip to perform a custom analysis, but you don’t want to go to the trouble of actually creating the CDF files and Bioconductor packages normally required to do the analysis, and yet you want to test and develop your analysis method.

How?

It turns out you are in luck. At least for AffyBatch objects in Bioconductor (created by calling ReadAffy), the CDF information is stored as an attached environment that can be easily hacked and modified to your hearts content. Environments in R are quite important and useful, and I wouldn’t have come up with this if I hadn’t been working in R for the past couple of years, but figured someone else might benefit from this knowledge.

The environment

In R, one can access an environment like so:

get("objName", envName) # get the value of object in the environment
ls(envName)

What is also very cool, is that one can extract the objects in an environment to a list, and also create their own environment from a list using list2env. Using this methodology, we can create our own definition of probesets that can be used by standard Bioconductor routines to summarize the probes into probesets.

A couple of disclaimers:

  • I have only tried this on 3’ expression arrays
  • There might be a better way to do this, but I couldn’t find it (let me know in the comments)

Example

require(affy)
## Loading required package: affy
## Loading required package: BiocGenerics
## Loading required package: parallel
## 
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:parallel':
## 
##     clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
##     clusterExport, clusterMap, parApply, parCapply, parLapply,
##     parLapplyLB, parRapply, parSapply, parSapplyLB
## The following objects are masked from 'package:stats':
## 
##     IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
## 
##     anyDuplicated, append, as.data.frame, basename, cbind,
##     colMeans, colnames, colSums, dirname, do.call, duplicated,
##     eval, evalq, Filter, Find, get, grep, grepl, intersect,
##     is.unsorted, lapply, lengths, Map, mapply, match, mget, order,
##     paste, pmax, pmax.int, pmin, pmin.int, Position, rank, rbind,
##     Reduce, rowMeans, rownames, rowSums, sapply, setdiff, sort,
##     table, tapply, union, unique, unsplit, which, which.max,
##     which.min
## Loading required package: Biobase
## Welcome to Bioconductor
## 
##     Vignettes contain introductory material; view with
##     'browseVignettes()'. To cite Bioconductor, see
##     'citation("Biobase")', and for packages 'citation("pkgname")'.
require(estrogen)
## Loading required package: estrogen
require(hgu95av2cdf)
## Loading required package: hgu95av2cdf
## Warning: replacing previous import 'AnnotationDbi::tail' by 'utils::tail'
## when loading 'hgu95av2cdf'
## Warning: replacing previous import 'AnnotationDbi::head' by 'utils::head'
## when loading 'hgu95av2cdf'
## 
datadir = system.file("extdata", package="estrogen")

pd = read.AnnotatedDataFrame(file.path(datadir, "estrogen.txt"), header=TRUE, sep="", row.names=1)
pData(pd)
##              estrogen time.h
## low10-1.cel    absent     10
## low10-2.cel    absent     10
## high10-1.cel  present     10
## high10-2.cel  present     10
## low48-1.cel    absent     48
## low48-2.cel    absent     48
## high48-1.cel  present     48
## high48-2.cel  present     48
celDat = ReadAffy(filenames = rownames(pData(pd)), 
                  phenoData = pd,
                  verbose=TRUE, celfile.path=datadir)
## 1 reading /software/R_libs/R351_bioc37/estrogen/extdata/low10-1.cel ...instantiating an AffyBatch (intensity a 409600x8 matrix)...done.
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/low10-1.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/low10-2.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/high10-1.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/high10-2.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/low48-1.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/low48-2.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/high48-1.cel
## Reading in : /software/R_libs/R351_bioc37/estrogen/extdata/high48-2.cel

This loads up the data, reads in the raw data, and gets it ready for us to use. Now, lets see what is in the actual CDF environment.

topProbes <- head(ls(hgu95av2cdf)) # get a list of probesets
topProbes
## [1] "100_g_at"  "1000_at"   "1001_at"   "1002_f_at" "1003_s_at" "1004_at"
exSet <- get(topProbes[1], hgu95av2cdf)
exSet
##           pm     mm
##  [1,] 175218 175858
##  [2,] 356689 357329
##  [3,] 227696 228336
##  [4,] 237919 238559
##  [5,] 275173 275813
##  [6,] 203444 204084
##  [7,] 357984 358624
##  [8,] 368524 369164
##  [9,] 285352 285992
## [10,] 304510 305150
## [11,] 159937 160577
## [12,] 223929 224569
## [13,] 282764 283404
## [14,] 270003 270643
## [15,] 303343 303983
## [16,] 389048 389688

We can see here that the first probe set 100_g_at has 16 perfect-match and mis-match probes in associated with it.

Lets summarize the original data using RMA.

rma1 <- exprs(rma(celDat))
## Background correcting
## Normalizing
## Calculating Expression
head(rma1)
##           low10-1.cel low10-2.cel high10-1.cel high10-2.cel low48-1.cel
## 100_g_at     9.642896    9.741496     9.537036     9.353625    9.591697
## 1000_at     10.398169   10.254362    10.003971     9.903528   10.374866
## 1001_at      5.717613    5.881008     5.859563     5.954028    5.960540
## 1002_f_at    5.512596    5.801807     5.571065     5.608132    5.390064
## 1003_s_at    7.783927    8.007975     8.037999     7.835120    7.926487
## 1004_at      7.289162    7.603670     7.488539     7.771506    7.521789
##           low48-2.cel high48-1.cel high48-2.cel
## 100_g_at     9.570590     9.475796     9.530655
## 1000_at     10.033520    10.345066     9.863321
## 1001_at      6.020889     5.981080     6.285192
## 1002_f_at    5.494511     5.508104     5.630107
## 1003_s_at    8.138870     7.994937     8.233338
## 1004_at      7.599544     7.456149     7.675171

Now lets get the data as a list, and then create a new environment to be used for summarization.

allSets <- ls(hgu95av2cdf)
allSetDat <- mget(allSets, hgu95av2cdf)

allSetDat[1]
## $`100_g_at`
##           pm     mm
##  [1,] 175218 175858
##  [2,] 356689 357329
##  [3,] 227696 228336
##  [4,] 237919 238559
##  [5,] 275173 275813
##  [6,] 203444 204084
##  [7,] 357984 358624
##  [8,] 368524 369164
##  [9,] 285352 285992
## [10,] 304510 305150
## [11,] 159937 160577
## [12,] 223929 224569
## [13,] 282764 283404
## [14,] 270003 270643
## [15,] 303343 303983
## [16,] 389048 389688
hgu2 <- list2env(allSetDat)
celDat@cdfName <- "hgu2"

rma2 <- exprs(rma(celDat))
## Background correcting
## Normalizing
## Calculating Expression
head(rma2)
##           low10-1.cel low10-2.cel high10-1.cel high10-2.cel low48-1.cel
## 100_g_at     9.642896    9.741496     9.537036     9.353625    9.591697
## 1000_at     10.398169   10.254362    10.003971     9.903528   10.374866
## 1001_at      5.717613    5.881008     5.859563     5.954028    5.960540
## 1002_f_at    5.512596    5.801807     5.571065     5.608132    5.390064
## 1003_s_at    7.783927    8.007975     8.037999     7.835120    7.926487
## 1004_at      7.289162    7.603670     7.488539     7.771506    7.521789
##           low48-2.cel high48-1.cel high48-2.cel
## 100_g_at     9.570590     9.475796     9.530655
## 1000_at     10.033520    10.345066     9.863321
## 1001_at      6.020889     5.981080     6.285192
## 1002_f_at    5.494511     5.508104     5.630107
## 1003_s_at    8.138870     7.994937     8.233338
## 1004_at      7.599544     7.456149     7.675171

What about removing the MM columns? RMA only uses the PM, so it should still work.

allSetDat <- lapply(allSetDat, function(x){
  x[,1, drop=F]
})

allSetDat[1]
## $`100_g_at`
##           pm
##  [1,] 175218
##  [2,] 356689
##  [3,] 227696
##  [4,] 237919
##  [5,] 275173
##  [6,] 203444
##  [7,] 357984
##  [8,] 368524
##  [9,] 285352
## [10,] 304510
## [11,] 159937
## [12,] 223929
## [13,] 282764
## [14,] 270003
## [15,] 303343
## [16,] 389048
hgu3 <- list2env(allSetDat)
celDat@cdfName <- "hgu3"
rma3 <-exprs(rma(celDat))
## Background correcting
## Normalizing
## Calculating Expression
head(rma3)
##           low10-1.cel low10-2.cel high10-1.cel high10-2.cel low48-1.cel
## 100_g_at     9.642896    9.741496     9.537036     9.353625    9.591697
## 1000_at     10.398169   10.254362    10.003971     9.903528   10.374866
## 1001_at      5.717613    5.881008     5.859563     5.954028    5.960540
## 1002_f_at    5.512596    5.801807     5.571065     5.608132    5.390064
## 1003_s_at    7.783927    8.007975     8.037999     7.835120    7.926487
## 1004_at      7.289162    7.603670     7.488539     7.771506    7.521789
##           low48-2.cel high48-1.cel high48-2.cel
## 100_g_at     9.570590     9.475796     9.530655
## 1000_at     10.033520    10.345066     9.863321
## 1001_at      6.020889     5.981080     6.285192
## 1002_f_at    5.494511     5.508104     5.630107
## 1003_s_at    8.138870     7.994937     8.233338
## 1004_at      7.599544     7.456149     7.675171

What if we only want to use the first 5 probesets?

allSetDat <- allSetDat[1:5]
hgu4 <- list2env(allSetDat)
celDat@cdfName <- "hgu4"
celDat
## AffyBatch object
## size of arrays=640x640 features (21 kb)
## cdf=hgu4 (5 affyids)
## number of samples=8
## number of genes=5
## annotation=hgu95av2
## notes=
rma4 <- exprs(rma(celDat))
## Background correcting
## Normalizing
## Calculating Expression
rma4
##           low10-1.cel low10-2.cel high10-1.cel high10-2.cel low48-1.cel
## 100_g_at     9.463007    9.554665     9.449050     9.401976    9.447562
## 1000_at     10.182753   10.009785    10.009785     9.970396   10.102424
## 1001_at      5.943840    6.005177     5.944015     6.089531    6.237329
## 1002_f_at    5.787166    5.846225     5.816964     5.814798    5.763175
## 1003_s_at    7.750877    7.769401     7.913021     7.864052    7.860778
##           low48-2.cel high48-1.cel high48-2.cel
## 100_g_at     9.457986     9.401366     9.431078
## 1000_at     10.009785    10.197065     9.889555
## 1001_at      6.147957     6.189200     6.206669
## 1002_f_at    5.763175     5.740757     5.755085
## 1003_s_at    7.917565     7.862614     7.928691
dim(rma4)
## [1] 5 8

Custom CDF

To generate our custom CDF, we are going to set our own names, and take random probes from all of the probes on the chip. The actual criteria of which probes should be together can be made using any method the author chooses.

maxIndx <- 640*640

customCDF <- lapply(seq(1,100), function(x){
  tmp <- matrix(sample(maxIndx, 20), nrow=20, ncol=1)
  colnames(tmp) <- "pm"
  return(tmp)
})

names(customCDF) <- seq(1, 100)

hgu5 <- list2env(customCDF)
celDat@cdfName <- "hgu5"
rma5 <- exprs(rma(celDat))
## Background correcting
## Normalizing
## Calculating Expression
head(rma5)
##     low10-1.cel low10-2.cel high10-1.cel high10-2.cel low48-1.cel
## 1      7.129135    7.098127     7.068482     7.233563    7.374926
## 10     6.032985    5.979999     5.986230     5.932740    5.826887
## 100    6.682842    6.456283     6.495291     6.519710    6.659279
## 11     7.851354    7.915515     7.902346     8.095335    7.902186
## 12     6.467172    6.447098     6.597328     6.417236    6.591901
## 13     7.093439    7.165885     7.160012     7.123641    7.146212
##     low48-2.cel high48-1.cel high48-2.cel
## 1      7.058663     7.007227     7.056640
## 10     6.046502     6.332338     6.120624
## 100    6.611351     6.525982     6.488238
## 11     7.948811     8.141210     7.987187
## 12     6.474024     6.469687     6.303738
## 13     7.175856     7.408415     7.137689

I hope this information is useful to someone else. I know it made my life a lot easier.

sessionInfo()
## R version 3.5.1 (2018-07-02)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
## 
## Matrix products: default
## BLAS: /software/R-351/lib/R/lib/libRblas.so
## LAPACK: /software/R-351/lib/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] parallel  stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] hgu95av2cdf_2.18.0  estrogen_1.26.0     affy_1.58.0        
## [4] Biobase_2.42.0      BiocGenerics_0.28.0
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.2            AnnotationDbi_1.42.1  knitr_1.24           
##  [4] magrittr_1.5          IRanges_2.14.12       zlibbioc_1.28.0      
##  [7] bit_1.1-14            blob_1.1.1            stringr_1.4.0        
## [10] tools_3.5.1           xfun_0.8              DBI_1.0.0            
## [13] htmltools_0.3.6       bit64_0.9-7           yaml_2.2.0           
## [16] digest_0.6.20         preprocessCore_1.42.0 bookdown_0.9         
## [19] affyio_1.50.0         S4Vectors_0.18.3      memoise_1.1.0        
## [22] evaluate_0.14         RSQLite_2.1.1         rmarkdown_1.15       
## [25] blogdown_0.10         stringi_1.4.3         compiler_3.5.1       
## [28] BiocInstaller_1.30.0  stats4_3.5.1

Originally published 2013/07/13, moved to http://rmflight.github.io on 2013/12/04.