crossmeta streamlines the cross-platform effect size and pathway meta-analysis of microarray data. For the analysis, you will need a list of Affymetrix, Illumina, and/or Agilent GSE numbers from GEO. All 21 species in the current homologene build are supported.


Obtaining Raw Data

Search GEO to find relevant microarray data for the meta-analysis. For this example, I searched for the PI3K inhibitor LY-294002. The search was filtered as follows:

This search produced 35 hits from which five were chosen. In practice, it is best to select all GSEs matching the necessary criteria (supplementary raw data files available, multiple samples for each treatment, and from either Affymetrix, Illumina, or Agilent arrays).

After identifying GSEs for the meta-analysis, download and decompress the raw data as follows. For raw Illumina data only, you may need to check that it has the correct format (and edit it if needed).

library(crossmeta)

# specify where data will be downloaded
data_dir <- file.path(getwd(), "data", "LY")

# gather all GSEs
gse_names  <- c("GSE9601", "GSE15069", "GSE50841", "GSE34817", "GSE29689")

# gather Illumina GSEs (see 'Checking Raw Illumina Data')
illum_names <- c("GSE50841", "GSE34817", "GSE29689")

# download raw data
# get_raw(gse_names, data_dir)

Checking Raw Illumina Data

It is difficult to automate loading raw Illumina data files because they lack a standardized format. crossmeta will attempt to fix the headers of raw Illumina data files so that they can be loaded. If crossmeta fails, you will have to edit the headers of the raw Illumina data files yourself or omit the offending studies.

To edit raw Illumina data headers, I recommend that you download and set Sublime Text 2 as your default text editor. It has very nice regular expression capabilities. Here is a good regular expression cheat-sheat.

Raw illumina files will be in data_dir in a seperate folder for each GSE. They are usually .txt files and include non-normalized in their name. Ensure the following:

Also ensure that column names have the following format:

To open these files one at a time with your default text editor:

# this is why we gathered Illumina GSEs
open_raw_illum(illum_names, data_dir)

To illustrate (crossmeta fixes these particular headers automatically):

For GSE50841:

For GSE34817:

For GSE29689:

A bioconductor data package (lydata) is available where all the above was performed. This data package will be used for subsequent demonstrations.

library(lydata)

# location of raw data
data_dir <- system.file("extdata", package = "lydata")

Loading and Annotating Data

After downloading the raw data, it must be loaded and annotated. The necessary bioconductor annotation data packages will be downloaded as needed.

# reloads if previously called
esets <- load_raw(gse_names, data_dir)

If crossmeta does not know the bioconductor annotation data package for a given platform, entry will be requested. Use the given platform (GPL) to search online for the correct package name. For example, entry was requested for GSE29689 (platform GPL6883). A search for GPL6883 on GEO identified the title for this platform as Illumina HumanRef-8 v3.0 expression beadchip. A subsequent search on Bioconductor for illuminahuman identified the appropriate package as illuminaHumanv3.

If you can’t find the bioconductor annotation data package (or none exists), type enter with a blank entry and crossmeta will attempt annotation using entrez gene ids in the feature data from the study in question. If entrez ids are absent, annotation will fail. To proceed, I reccommend that you add entrez ids and then use crossmeta to map these to hgnc symbols. As an example, if annotation had failed for GSE15069:

library(Biobase)
library(AnnotationDbi)

# check feature data to see what columns are available
head(fData(esets$GSE15069))

# if using RStudio
# View(fData(esets$GSE15069))

# annotation package for appropriate species
library(org.Mm.eg.db)

# map from accession number to entrez gene ids
acnums  <- as.character(fData(esets$GSE15069)$GB_ACC)
enids   <- mapIds(org.Mm.eg.db, acnums, "ENTREZID", "ACCNUM")

# add 'GENE_ID' column with entrez ids
fData(esets$GSE15069)$GENE_ID <- enids 

# use crossmeta to map from entrez gene ids to homologous hgnc symbol
esets$GSE15069 <- symbol_annot(esets$GSE15069)

# to overwrite saved eset (to avoid repeating above)
saveRDS(esets$GSE15069, file.path(data_dir, "GSE15069", "GSE15069_eset.rds"))

Differential Expression

After loading and annotating the data, you may proceed with differential expression analysis:

anals <- diff_expr(esets, data_dir)

You will be prompted to type group names and select samples for each control and test group that you wish to compare. Once done for a study, click Done and you will be prompted to do the same for the next study. You can re-select a previous control group by selecting the group name in the drop down box. You can also view and delete current contrasts in the Contrasts tab.

If you need to look up experimental details for the study (e.g. group membership) Click the GEO button in the top left corner.

For Illumina samples, If only titles are shown (no accessions), then samples from raw Illumina data could not be matched confidently to samples from the processed Illumina data. In this case, the sample titles are from the raw data and you may need to look in the individual sample records on GEO to identify which sample belongs to which group (not always possible).

If you need to re-run the above analysis, you can avoid having to reselect samples and rename groups:

# load auto-saved results of previous call to diff_expr
prev <- load_diff(gse_names, data_dir)

# supply prev to diff_expr
# anals <- diff_expr(esets, data_dir, prev_anals=prev)

Multi-dimensional scaling plots are generated for each GSE. These plots are generated from expression data with the effects of surrogate variables removed.

Non-GUI Selection

Using the graphical user interface for sample selection can be error prone and time consuming for GSEs with a large number of samples. For these cases, you may prefer to specify group membership of samples using sample information from the study in question. As an example:

library(Biobase)

# load eset
gse_name  <- c("GSE34817")
eset <- load_raw(gse_name, data_dir)

# inspect pData of eset
# View(pData(eset$GSE34817))  # if using RStudio
head(pData(eset$GSE34817))    # otherwise

# get group info from pData (differs based on eset)
group <- pData(eset$GSE34817)$characteristics_ch1.1

# make group names concise and valid
group <- gsub("treatment: ", "", group)
group <- make.names(group)

# add group to eset pData
pData(eset$GSE34817)$group <- group

# setup selections
sel <- setup_prev(eset, contrasts = "LY-DMSO")

# run differential expression analysis
# anal <- diff_expr(eset, data_dir, prev_anal = sel)

Tissue Sources

Pathway and effect-size meta-analyses can be performed for all studies and seperately for each tissue source. To do so, first add tissue sources and specify any sources that should be paired (treated as the same source for subsequent meta-analyses).

# run GUI to add tissue sources
# anals <- add_sources(prev, data_dir)

# for usage details
?add_sources

Effect Size Meta-Analysis

After tissue sources have been specified for differential expression analyses, overall effect sizes can be determined through meta-analysis. The meta-analysis method determines an overall standardized effect size and false discovery rate for each gene that is present in a specified fraction of studies (30% by default).

# re-load previous analyses if need to
prev <- load_diff(gse_names, data_dir)
anals <- diff_expr(esets, data_dir, prev_anals = prev)
## Warning: Not all samples grouped: see https://support.bioconductor.org/p/73107/
## #73109

## Warning: Not all samples grouped: see https://support.bioconductor.org/p/73107/
## #73109

## Warning: Not all samples grouped: see https://support.bioconductor.org/p/73107/
## #73109

## Warning: Not all samples grouped: see https://support.bioconductor.org/p/73107/
## #73109

## Warning: Not all samples grouped: see https://support.bioconductor.org/p/73107/
## #73109
# perform meta analyses by tissue source
es_res <- es_meta(anals, by_source = TRUE)

# for explanation of values
?es_meta

The meta-analysis method was adapted from GeneMeta. Differences include the use of moderated unbiased effect sizes calculated by metaMA and the allowance for genes measured in only a fraction of studies.

sessionInfo()
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] limma_3.50.0        Biobase_2.54.0      BiocGenerics_0.40.0
## [4] lydata_1.19.0       crossmeta_1.20.0   
## 
## loaded via a namespace (and not attached):
##  [1] locfit_1.5-9.4         Rcpp_1.0.7             lattice_0.20-45       
##  [4] png_0.1-7              Biostrings_2.62.0      digest_0.6.28         
##  [7] mime_0.12              R6_2.5.1               GenomeInfoDb_1.30.0   
## [10] stats4_4.1.1           RSQLite_2.2.8          evaluate_0.14         
## [13] sva_3.42.0             httr_1.4.2             zlibbioc_1.40.0       
## [16] rlang_0.4.12           data.table_1.14.2      annotate_1.72.0       
## [19] miniUI_0.1.1.1         jquerylib_0.1.4        blob_1.2.2            
## [22] S4Vectors_0.32.0       Matrix_1.3-4           rmarkdown_2.11        
## [25] splines_4.1.1          shinyjs_2.0.0          BiocParallel_1.28.0   
## [28] stringr_1.4.0          RCurl_1.98-1.5         bit_4.0.4             
## [31] shiny_1.7.1            compiler_4.1.1         httpuv_1.6.3          
## [34] xfun_0.27              metaMA_3.1.2           mgcv_1.8-38           
## [37] htmltools_0.5.2        KEGGREST_1.34.0        GenomeInfoDbData_1.2.7
## [40] edgeR_3.36.0           matrixStats_0.61.0     IRanges_2.28.0        
## [43] XML_3.99-0.8           crayon_1.4.1           later_1.3.0           
## [46] bitops_1.0-7           grid_4.1.1             nlme_3.1-153          
## [49] jsonlite_1.7.2         xtable_1.8-4           lifecycle_1.0.1       
## [52] DBI_1.1.1              magrittr_2.0.1         stringi_1.7.5         
## [55] cachem_1.0.6           XVector_0.34.0         genefilter_1.76.0     
## [58] promises_1.2.0.1       fdrtool_1.2.16         bslib_0.3.1           
## [61] ellipsis_0.3.2         vctrs_0.3.8            tools_4.1.1           
## [64] bit64_4.0.5            SMVar_1.3.3            parallel_4.1.1        
## [67] fastmap_1.1.0          survival_3.2-13        yaml_2.2.1            
## [70] AnnotationDbi_1.56.0   memoise_2.0.0          knitr_1.36            
## [73] sass_0.4.0