1 Introduction

The signatureSearchData package provides access to the reference data used by the associated signatureSearch software package. The latter allows to search with a query gene expression signature (GES) a database of treatment GESs to identify cellular states sharing similar expression responses (connections). This way one can identify drugs or gene knockouts that induce expression phenotypes similar to a sample of interest. The resulting associations may lead to novel functional insights how perturbagens of interest interact with biological systems.

Currently, signatureSearchData includes GES data from the CMap (Connectivity Map) and LINCS (Library of Network-Based Cellular Signatures) projects that are largely based on drug and genetic perturbation experiments performed on variable numbers of human cell lines (Lamb et al. 2006; Subramanian et al. 2017). In signatureSearchData these data sets have been preprocessed to be compatible with the different gene expression signature search (GESS) algorithms implemented in signatureSearch. The preprocessed data types include but are not limited to normalized gene expression values (e.g. intensity values), log fold changes (LFC) and Z-scores, p-values or FDRs of differentially expressed genes (DEGs), rankings based on selected preprocessing routines or sets of top up/down-regulated DEGs.

The CMap data were downloaded from the CMap project site (Version build02). The latter is a collection of over 7,000 gene expression profiles (signatures) obtained from perturbation experiments with 1,309 drug-like small molecules on five human cancer cell lines. The Affymetrix Gene Chip technology was used to generate the CMAP2 data set.

In 2017, the LINCS Consortium generated a similar but much larger data set where the total number of gene expression signatures was scaled up to over one million. This was achieved by switching to a much more cost effective gene expression profiling technology called L1000 assay (Peck et al. 2006; Edgar, Domrachev, and Lash 2002). The current set of perturbations covered by the LINCS data set includes 19,811 drug-like small molecules applied at variable concentrations and treatment times to ~70 human non-cancer (normal) and cancer cell lines. Additionally, it includes several thousand genetic perturbagens composed of gene knockdown and over-expression experiments.

The data structures and search algorithms used by signatureSearch and signatureSearchData are designed to work with most genome-wide expression data including hybridization-based methods, such as Affymetrix or L1000, as well as sequencing-based methods, such as RNA-Seq. Currently, signatureSearchData does not include preconfigured RNA-Seq reference data mainly due to the lack of large-scale perturbation studies (e.g. drug-based) available in the public domain that are based on RNA-Seq. This situation may change in the near future once the technology has become more affordable for this purpose.

2 Install and Load Package

signatureSearchData is a R/Bioconductor package and can be installed using BiocManager::install().

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("signatureSearchData")

After the package is installed, it can be loaded in an R session as follows.

library(signatureSearchData)

3 Explore Data Sets

A summary of the data sets provided by the signatureSearchData package can be obtained with the query function of the ExperimentHub package. The information is stored in an object of class ExperimentHub, here assigned to ssd.

library(ExperimentHub)
eh <- ExperimentHub()
ssd <- query(eh, c("signatureSearchData"))
ssd
## ExperimentHub with 12 records
## # snapshotDate(): 2019-10-22 
## # $dataprovider: Broad Institute, GO, DrugBank, Broad Institute, STITCH
## # $species: Homo sapiens
## # $rdataclass: character, data.frame, environment, list
## # additional mcols(): taxonomyid, genome, description,
## #   coordinate_1_based, maintainer, rdatadateadded, preparerclass,
## #   tags, rdatapath, sourceurl, sourcetype 
## # retrieve records with, e.g., 'object[["EH3223"]]' 
## 
##            title       
##   EH3223 | cmap        
##   EH3224 | cmap_expr   
##   EH3225 | cmap_rank   
##   EH3226 | lincs       
##   EH3227 | lincs_expr  
##   ...      ...         
##   EH3230 | goAnno_drug 
##   EH3231 | GO_DATA     
##   EH3232 | GO_DATA_drug
##   EH3233 | taurefList  
##   EH3234 | ES_NULL

The titles of the data sets can be returned with ssd$title.

ssd$title
##  [1] "cmap"               "cmap_expr"          "cmap_rank"         
##  [4] "lincs"              "lincs_expr"         "dtlink_db_clue_sti"
##  [7] "goAnno"             "goAnno_drug"        "GO_DATA"           
## [10] "GO_DATA_drug"       "taurefList"         "ES_NULL"

More detailed information about each data set can be returned as a list, below subsetted to 10th entry with [10].

as.list(ssd)[10]
## $EH3232
## ExperimentHub with 1 record
## # snapshotDate(): 2019-10-22 
## # names(): EH3232
## # package(): signatureSearchData
## # $dataprovider: GO
## # $species: Homo sapiens
## # $rdataclass: environment
## # $rdatadateadded: 2019-10-22
## # $title: GO_DATA_drug
## # $description: GO annotation environment after drug mappings
## # $taxonomyid: 9606
## # $genome: GRCh38
## # $sourcetype: MySQL
## # $sourceurl: https://bioconductor.org/packages/release/data/annotation...
## # $sourcesize: NA
## # $tags: c("ExperimentHub", "ExperimentData") 
## # retrieve record with 'object[["EH3232"]]'

Details about the usage of ExperimentHub can be found in its vignettes here.

4 LINCS Signature Database

The L1000 assay, used for generating the LINCS data, measures the expression of 978 landmark genes and 80 control genes by loading amplified mRNA populations onto beads and then detecting their abundance with a fluorescent-based method (Peck et al. 2006). The expression of 11,350 additional genes is imputed from the landmark genes by using as training data a large collection of Affymetrix gene chips (Edgar, Domrachev, and Lash 2002).

The LINCS data have been pre-processed by the Broad Institute to 5 different levels and are available for download from GEO. Level 1 data are the raw mean fluorescent intensity values that come directly from the Luminex scanner. Level 2 data are the expression intensities of the 978 landmark genes. They have been normalized and used to impute the expression of the additional 11,350 genes, forming Level 3 data. A robust z-scoring procedure was used to generate differential expression values from the normalized profiles (Level 4). Finally, a moderated z-scoring procedure was applied to the replicated samples of each experiment (mostly 3 replicates) to compute a weighted average signature (Level 5). For a more detailed description of the preprocessing methods used by the LINCS project, readers want to refer to the LINCS user guide.

Disregarding replicates, the LINCS data set contains 473,647 signatures with unique cell type and treatment combinations. This includes 19,811 drug-like small molecules tested on different cell lines at multiple concentrations and treatment times. In addition to compounds, several thousand genetic perturbations (gene knock-downs and over expressions) have been tested. Currently, the data described in this vignette are restricted to signatures of small molecule treatments across different cells lines. However, users have the option to assemble any custom collection of the LINCS data. For consistency, only signatures at one specific concentration (10\(\mu\)M) and one time point (24h) have been selected for each small molecule in the default collection. These choices are similar to the conditions used in primary high-throughput compound screens of cell lines. Since the selected compound concentrations and treatment duration have not been tested by LINCS across all cell types yet, a subset of compounds had to be selected that best met the chosen treatment requirements. This left us with 8,104 compounds that were uniformly tested at the chosen concentration and treatment time, but across variable numbers of cell lines. The total number of expression signatures meeting this requirement is 45,956, while the total number of cell lines included in this data set is 30.

4.1 Z-scores from ExperimentHub

The LINCS sub-dataset, filtered and assembled according to the above criteria, can be downloaded from Bioconductor’s ExperimentHub as HDF5 file. In the example below, the path to this file is assigned to a character vector called lincs_path. A summary of the content of the HDF5 file can be returned with the h5ls function. Note, due to the large size of the LINCS data set, its download takes too much time to evaluate the following code section during the build time of this vignette.

library(ExperimentHub); library(rhdf5)
eh <- ExperimentHub()
query(eh, c("signatureSearchData", "lincs"))
lincs_path <- eh[['EH3226']]
rhdf5::h5ls(lincs_path) 

In this case the loaded data instance includes moderated Z-scores from DE analyses of 12,328 genes for 8,140 compound treatments across a total of 30 cell lines corresponding to 45,956 expression signatures. This data set can be used by all set-based and correlation-based GESS methods provided by the signatureSearch package.

4.2 Z-scores from GEO

The following explains how to generate the above LINCS data object from scratch. This also illustrates how to filter the LINCS level 5 data in other ways.

4.2.1 Download Level 5 Data

Download and unzip the following files from GEO entry GSE92742:

  • GSE92742_Broad_LINCS_Level5_COMPZ.MODZ_n473647x12328.gctx.gz
  • GSE92742_Broad_LINCS_gene_info.txt.gz
  • GSE92742_Broad_LINCS_sig_info.txt.gz

The following code examples assume that the downloaded datasets are stored in a sub-directory called data. All paths in this vignette are given relative to the present working directory of a user’s R session.

4.2.2 Filter Signatures

The following selects LINCS Level 5 signatures of compound treatments at a concentration of 10\(\mu\)M and a treatment time of 24 hours. Note, the import command below may issue a warning message that can be ignored.

meta42 <- readr::read_tsv("./data/GSE92742_Broad_LINCS_sig_info.txt") 
dose <- meta42$pert_idose[7]
## filter rows by 'pert_type' as compound, 10uM concentration, and 24h treatment time
meta42_filter <- sig_filter(meta42, pert_type="trt_cp", dose=dose, time="24 h") # 45956 X 14

4.2.3 Z-score Data in HDF5

Next, the large Z-score matrix of expression signatures is imported step-wise in subsets of manageable sizes and then appended to an HDF5 file (here lincs.h5). In this vignette, the latter is referred to as the LINCS Z-score database. Since the size of the full matrix is several GBs in size, it would consume too much memory to be read into R at once. Reading the matrix in smaller batches and appending them to an HDF5 file is much more memory efficient. Subsequently, the readHDF5chunk function, defined by the signatureSearchData package, imports user-definable subsets of the data from the HDF5 file into a SummarizedExperiment object, here assigned to se.

library(signatureSearch)
gctx <- "./data/GSE92742_Broad_LINCS_Level5_COMPZ.MODZ_n473647x12328.gctx"
gctx2h5(gctx, cid=meta42_filter$sig_id, new_cid=meta42_filter$pert_cell_factor,
        h5file="./data/lincs.h5", by_ncol=5000, overwrite=TRUE)
se <- readHDF5chunk(h5file="./data/lincs.h5", colindex=1:5000)

4.3 Intensities from ExperimentHub

The LINCS Level 3 data can be downloaded from GEO the same way as described above for the Level 5 data. The Level 3 data contain normalized gene expression values across all treatments and cell lines used by LINCS. The Level 3 signatures were filtered using the same dosage and duration criteria as the Level 5 data. The biological replicate information included in the Level 3 data were collapsed to mean values. Subsequently, the resulting matrix of mean expression values was written to an HDF5 file. The latter is referred to as lincs_expr database containing 38,824 signatures for a total of 5,925 small molecule treatments and 30 cell lines. Although the LINCS Level 3 and 5 data are filtered here the same way, the number of small molecules represented in the Level 3 data (5,925) is smaller than in the Level 5 data (8,140). The reason for this inconsistency is most likely that the Level 3 dataset, downloadable from GEO, is incomplete.

The filtered and processed LINCS Level3 data (lincs_expr) can be loaded from Bioconductor’s ExperimentHub interface as follows.

library(ExperimentHub)
eh <- ExperimentHub()
query(eh, c("signatureSearchData", "lincs_expr"))
lincs_expr_path <- eh[['EH3227']]

In this case the loaded lincs_expr instance includes mean expression values of 12,328 genes for 5,925 compound treatments across a total of 30 cell lines. This data set can be used by all correlation-based GESS methods provided by the signatureSearch package.

4.4 Intensities from GEO

The following steps explain how to generate the above data set from scratch. This also illustrates how to filter the LINCS Level 3 data in other ways.

4.4.1 Download Level 3 Data

Download and unzip the following files from GEO entry GSE92742:

  • GSE92742_Broad_LINCS_Level3_INF_mlr12k_n1319138x12328.gctx.gz
  • GSE92742_Broad_LINCS_gene_info.txt.gz
  • GSE92742_Broad_LINCS_inst_info.txt.gz

As above, the following code examples assume that the downloaded datasets are stored in a sub-directory called data. All paths in this vignette are given relative to the present working directory of a user’s R session.

4.4.2 Filter Signatures

The following selects LINCS Level 3 signatures of compound treatments at a concentration of 10\(\mu\)M and a treatment time of 24 hours.

inst42 <- readr::read_tsv("./data/GSE92742_Broad_LINCS_inst_info.txt") 
inst42_filter <- inst_filter(inst42, pert_type="trt_cp", dose=10, dose_unit="um", 
                             time=24, time_unit="h") # 169,795 X 13

4.4.3 Mean Intensities in HDF5

Next, mean expression values are calculated among biological replicates and then appended in batches to the corresponding HDF5 file.

# It takes some time
meanExpr2h5(gctx="./data/GSE92742_Broad_LINCS_Level3_INF_mlr12k_n1319138x12328.gctx",
            inst=inst42_filter, h5file="./data/lincs_expr.h5") # 12328 X 38824

5 CMap2 Signature Database

CMap2 (Version build02) contains GESs for 1,309 drugs and eight cell lines that were generated with Affymetrix Gene Chips as expression platform. In some cases this includes drug treatments at different concentrations and time points. For consistency, the CMap2 data was reduced to drug treatments with concentrations and time points that are comparable to those used for the above LINCS data. CMap2 data can be downloaded from GEO or its project site either in raw format or as rank transformed matrix. The ranks are based on DEG analyses of drug treatments (drug vs. no-drug) where the resulting Z-scores were used to generate the rank matrix. The latter was used here and is referred to as rankMatrix. The Affymetrix probe set identifiers stored in the row name slot of this matrix were translated into gene identifies. To obtain a matrix with unique gene identifiers, the ranks for genes represented by more than one probe set were averaged and then re-ranked accordingly. This final gene level rank matrix, referred to as cmap_rank, contains rank profiles for 12,403 genes from 1,309 compound treatments in up to 5 cells corresponding to a total of 3,587 treatment signatures. This matrix can be used for all GESS methods in the signatureSearch package that are compatible with rank data, such as the gess_cmap method.

5.1 Rank Matrix from ExperimentHub

The cmap_rank data can be downloaded from Bioconductor’s ExperimentHub as HDF5 file. Since CMap2 is much smaller than LINCS, it can be imported in its entirety into a SummarizedExperiment object (here assigned to se) without excessive memory requirements.

library(ExperimentHub)
eh <- ExperimentHub()
query(eh, c("signatureSearchData", "cmap_rank"))
cmap_rank_path <- eh[["EH3225"]]
se <- readHDF5chunk(h5file=cmap_rank_path, colindex=1:3587)

5.2 Rank Matrix from Sources

The following steps explain how to generate the above CMap2 rank data set from scratch.

5.2.1 Download Rank Data

The rankMatrix can be downloaded from the CMap project site here. The specific file to download from this site is rankMatrix.txt.zip. As before, it should be saved and unzipped in the data directory of a user’s R session.

5.2.2 Filter instances

The following selects from rankMatrix for each compound the chosen treatment concentration and time point. This is achieved with help of the experiment annotation file cmap_instances_02.txt, also available from the CMap project site. Since this file is relatively small it has been included in the signatureSearchData package from where it can be loaded into R as shown below.

path <- system.file("extdata", "cmap_instances_02.txt", package="signatureSearchData")
cmap_inst <- read.delim(path, check.names=FALSE) 
inst_id <- cmap_inst$instance_id[!duplicated(paste(cmap_inst$cmap_name, cmap_inst$cell2, sep="_"))]
rankM <- read.delim("./data/rankMatrix.txt", check.names=FALSE, row.names=1) # 22283 X 6100
rankM_sub <- rankM[, as.character(inst_id)]
colnames(rankM_sub) <- unique(paste(cmap_inst$cmap_name, cmap_inst$cell2, "trt_cp", sep="__"))

5.2.3 Annotated Gene Level Data

5.2.3.1 Annotation Information

The following generates annotation information for Affymetirx probe set identifiers. Note, the three different Affymetrix chip types (HG-U133A, HT_HG-U133A, U133AAofAv2) used by CMap2 share nearly all probe set identifiers, meaning it is possible to use the same annotation package (here hgu133a.db) for all three.

library(hgu133a.db)
myAnnot <- data.frame(ACCNUM=sapply(contents(hgu133aACCNUM), paste, collapse=", "), 
                      SYMBOL=sapply(contents(hgu133aSYMBOL), paste, collapse=", "), 
                      UNIGENE=sapply(contents(hgu133aUNIGENE), paste, collapse=", "), 
                      ENTREZID=sapply(contents(hgu133aENTREZID), paste, collapse=", "), 
                      ENSEMBL=sapply(contents(hgu133aENSEMBL), paste, collapse=", "), 
                      DESC=sapply(contents(hgu133aGENENAME), paste, collapse=", "))
saveRDS(myAnnot, "./data/myAnnot.rds")

5.2.3.2 Gene Level Data

The probe2gene function transforms probe set to gene level data. If genes are represented by several probe sets then their mean intensities are used.

rankM_sub_gene <- probe2gene(rankM_sub, myAnnot) 

5.2.4 Store Rank Matrix in HDF5 file

The sub-setted rankMatrix is written to an HDF5 file, referred to as cmap_rank database.

matrix2h5(rankM_sub_gene, "./data/cmap_rank.h5", overwrite=TRUE) # 12403 X 3587
rhdf5::h5ls("./data/cmap_rank.h5")
## Read in cmap_rank.h5 as SummarizedExperiment object by chunks
cmap_rank_se <- readHDF5chunk("./data/cmap_rank.h5", colindex=1:5)

5.3 Intensities from ExperimentHub

To search CMap2 with signatureSearch's correlation based GESS methods (gess_cor), normalized gene expression values (here intensities) are required where the biological replicate information has been collapsed to mean values. For this, the cmap_expr database has been created from CEL files, which are the raw data of the Affymetrix technology. To obtain normalized expression data, the CEL files were downloaded from the CMap project site, and then processed with the MAS5 algorithm. Gene level expression data was generated the same way as described above. Next, the gene expression values for different concentrations and treatment times of each compound and cell were averaged. Subsequently, the expression matrix was saved to an HDF5 file, referred to as the cmap_expr database. It represents mean expression values of 12,403 genes for 1,309 compound treatments in up to 5 cells (3,587 signatures in total).

The cmap_expr database can be downloaded as HDF5 file from Bioconductor’s ExperimentHub as follows.

library(ExperimentHub)
eh <- ExperimentHub()
query(eh, c("signatureSearchData", "cmap_expr"))
cmap_expr_path <- eh[["EH3224"]]

This data set can be used by all correlation-based GESS methods provided by the signatureSearch package.

5.4 Intensities from Sources

How to generate the above cmap_expr database from scratch is explained in the Supplementary Material section of this vignette (see Section 8).

6 Custom Signature Databases

Custom databases of GESs can be built with the build_custom_db function provided by the signatureSearch package. For this the user provides custom genome-wide gene expression data (e.g. for drug, disease or genetic perturbations) in a data.frame or matrix. The gene expression data can be most types of the pre-processed gene expression values described under section 1.4 of the signatureSearch vignette.

7 Additional Datasets

The signatureSearchData package also contains several annotation datasets, such as drug-target information of small molecules. They are required for signatureSearch's functional enrichment analysis (FEA) routines. Currently, most of these annotation data were downloaded from the following databases:

8 Supplemental Material

8.1 CMap2 Intensities from Sources

The following steps explain how to generate the cmap_expr database in subsection 5.3 from scratch. They are intended for expert users and have been included here for reproduciblity reasons.

8.1.1 Directory Structure

The large number of files processed in the next steps are organized in two sub-directories of a user’s R session. Input files will be stored in a data directory, while all output files will be written to a results directory.

dir.create("data"); dir.create("data/CEL"); dir.create("results") 

8.1.2 Download of CEL Files

The getCmapCEL function will download the 7,056 CEL files from the CMap project site, and save each of them to a subdirectory named CEL under data. Since this download step will take some time, the argument rerun has been assigned FALSE in the below function call to avoid running it accidentally. To execute the download, the argument rerun needs to be assigned TRUE. If the raw data are not needed, users can skip this time consuming step and work with the preprocessed cmap_expr database downloaded from the ExperimentHub instead.

getCmapCEL(rerun=FALSE) 

8.1.3 Determine Chip Type

The CMAP data set is based on three different Affymetrix chip types (HG-U133A, HT_HG-U133A and U133AAofAv2). The following extracts the chip type information from the downloaded CEL files and stores the information in an rds file with the path ./data/chiptype.rds.

library(affxparser)
celfiles <- list.files("./data/CEL", pattern=".CEL$")
chiptype <- sapply(celfiles, function(x) affxparser::readCelHeader(paste0("data/CEL/", x))$chiptype)
saveRDS(chiptype, "./data/chiptype.rds")

8.1.4 Normalization of CEL Files

The following processes the CEL files from each chip type separately using the MAS5 normalization algorithm. The results will be written to 3 subdirectores under data that are named after the chip type names. To reduce the memory consumption of this step, the CEL files are normalized in batches of 200. The normalization takes about 10 hours without parallelization. To save time, this process can be easily accelerated on a computer cluster.

chiptype <- readRDS("./data/chiptype.rds")
chiptype_list <- split(names(chiptype), as.character(chiptype))
normalizeCel(chiptype_list, batchsize=200, rerun=FALSE) 

Next the results from each chip type are assembled in a data frame. After this all three of these data frames are combined to a single one, here named mas5df.

chiptype_dir <- unique(chiptype)
combineResults(chiptype_dir, rerun=FALSE)
mas5df <- combineNormRes(chiptype_dir, norm_method="MAS5")

8.1.5 Gene Level Data

After moving the myAnnot.rds file from above into the data directory, the probe2gene function is used to transforms probe set to gene level data. If genes are represented by several probe sets then their mean intensities are used.

myAnnot <- readRDS("./data/myAnnot.rds") 
mas5df <- probe2gene(mas5df, myAnnot) 
saveRDS(mas5df,"./data/mas5df.rds")

8.1.6 Average Intensities Across Replicates

The following averages the normalized gene expression values for different concentrations, treatment times and replicates of compounds and cell types.

mas5df <- readRDS("./data/mas5df.rds") # dim: 12403 x 7056
path <- system.file("extdata", "cmap_instances_02.txt", package="signatureSearchData")
cmap_inst <- read.delim(path, check.names=FALSE) 
cmap_drug_cell_expr <- meanExpr(mas5df, cmap_inst) # dim: 12403 X 3587
saveRDS(cmap_drug_cell_expr, "./data/cmap_drug_cell_expr.rds")

8.1.7 Mean Intensities in HDF5

The normalized and averaged expression values are saved to an HDF5 file, referred to as cmap_expr database.

cmap_drug_cell_expr <- readRDS("./data/cmap_drug_cell_expr.rds")
## match colnames to '(drug)__(cell)__(factor)' format
colnames(cmap_drug_cell_expr) <- gsub("(^.*)_(.*$)", "\\1__\\2__trt_cp", 
                                      colnames(cmap_drug_cell_expr)) 
matrix2h5(cmap_drug_cell_expr, "./data/cmap_expr.h5", overwrite=TRUE)
## Read in cmap_expr.h5 by chunks
cmap_expr_se <- readHDF5chunk("./data/cmap_expr.h5", colindex=1:5)

9 Session Info

sessionInfo()
## R version 3.6.1 (2019-07-05)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.10-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        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] parallel  stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] ExperimentHub_1.12.0      AnnotationHub_2.18.0     
## [3] BiocFileCache_1.10.2      dbplyr_1.4.2             
## [5] BiocGenerics_0.32.0       signatureSearchData_1.0.0
## 
## loaded via a namespace (and not attached):
##   [1] backports_1.1.5               fastmatch_1.1-0              
##   [3] plyr_1.8.4                    igraph_1.2.4.1               
##   [5] lazyeval_0.2.2                GSEABase_1.48.0              
##   [7] splines_3.6.1                 BiocParallel_1.20.0          
##   [9] GenomeInfoDb_1.22.0           ggplot2_3.2.1                
##  [11] urltools_1.7.3                digest_0.6.22                
##  [13] htmltools_0.4.0               GOSemSim_2.12.0              
##  [15] viridis_0.5.1                 GO.db_3.10.0                 
##  [17] magrittr_1.5                  memoise_1.1.0                
##  [19] limma_3.42.0                  readr_1.3.1                  
##  [21] annotate_1.64.0               graphlayouts_0.5.0           
##  [23] matrixStats_0.55.0            R.utils_2.9.0                
##  [25] enrichplot_1.6.0              prettyunits_1.0.2            
##  [27] colorspace_1.4-1              blob_1.2.0                   
##  [29] rappdirs_0.3.1                ggrepel_0.8.1                
##  [31] xfun_0.10                     dplyr_0.8.3                  
##  [33] crayon_1.3.4                  RCurl_1.95-4.12              
##  [35] jsonlite_1.6                  graph_1.64.0                 
##  [37] genefilter_1.68.0             zeallot_0.1.0                
##  [39] survival_3.1-7                glue_1.3.1                   
##  [41] polyclip_1.10-0               gtable_0.3.0                 
##  [43] zlibbioc_1.32.0               XVector_0.26.0               
##  [45] DelayedArray_0.12.0           Rhdf5lib_1.8.0               
##  [47] scales_1.0.0                  DOSE_3.12.0                  
##  [49] DESeq_1.38.0                  DBI_1.0.0                    
##  [51] Rcpp_1.0.3                    viridisLite_0.3.0            
##  [53] xtable_1.8-4                  progress_1.2.2               
##  [55] gridGraphics_0.4-1            bit_1.1-14                   
##  [57] europepmc_0.3                 preprocessCore_1.48.0        
##  [59] stats4_3.6.1                  htmlwidgets_1.5.1            
##  [61] httr_1.4.1                    fgsea_1.12.0                 
##  [63] RColorBrewer_1.1-2            pkgconfig_2.0.3              
##  [65] XML_3.98-1.20                 R.methodsS3_1.7.1            
##  [67] farver_1.1.0                  ggplotify_0.0.4              
##  [69] tidyselect_0.2.5              rlang_0.4.1                  
##  [71] reshape2_1.4.3                later_1.0.0                  
##  [73] AnnotationDbi_1.48.0          munsell_0.5.0                
##  [75] BiocVersion_3.10.1            tools_3.6.1                  
##  [77] visNetwork_2.0.8              RSQLite_2.1.2                
##  [79] ggridges_0.5.1                evaluate_0.14                
##  [81] stringr_1.4.0                 fastmap_1.0.1                
##  [83] yaml_2.2.0                    knitr_1.25                   
##  [85] bit64_0.9-7                   tidygraph_1.1.2              
##  [87] purrr_0.3.3                   ggraph_2.0.0                 
##  [89] RBGL_1.62.1                   mime_0.7                     
##  [91] R.oo_1.23.0                   DO.db_2.9                    
##  [93] xml2_1.2.2                    compiler_3.6.1               
##  [95] curl_4.2                      interactiveDisplayBase_1.24.0
##  [97] affyio_1.56.0                 tibble_2.1.3                 
##  [99] tweenr_1.0.1                  geneplotter_1.64.0           
## [101] stringi_1.4.3                 lattice_0.20-38              
## [103] Matrix_1.2-17                 vctrs_0.2.0                  
## [105] GSEAlm_1.46.0                 signatureSearch_1.0.0        
## [107] pillar_1.4.2                  lifecycle_0.1.0              
## [109] BiocManager_1.30.9            triebeard_0.3.0              
## [111] data.table_1.12.6             cowplot_1.0.0                
## [113] bitops_1.0-6                  httpuv_1.5.2                 
## [115] GenomicRanges_1.38.0          qvalue_2.18.0                
## [117] R6_2.4.0                      affy_1.64.0                  
## [119] promises_1.1.0                gridExtra_2.3                
## [121] gCMAP_1.30.0                  IRanges_2.20.0               
## [123] MASS_7.3-51.4                 assertthat_0.2.1             
## [125] rhdf5_2.30.0                  SummarizedExperiment_1.16.0  
## [127] Category_2.52.1               S4Vectors_0.24.0             
## [129] GenomeInfoDbData_1.2.2        hms_0.5.2                    
## [131] clusterProfiler_3.14.0        grid_3.6.1                   
## [133] tidyr_1.0.0                   rvcheck_0.1.6                
## [135] rmarkdown_1.16                ggforce_0.3.1                
## [137] Biobase_2.46.0                shiny_1.4.0

10 Funding

This project is funded by NIH grants U19AG02312 and U24AG051129 awarded by the National Institute on Aging (NIA). Subcomponents of the environment are based on methods developed by projects funded by NSF awards ABI-1661152 and PGRP-1810468. The High-Performance Computing (HPC) resources used for optimizing and applying the code of this project were funded by NIH and NSF grants 1S10OD016290-01A1 and MRI-1429826, respectively.

References

Edgar, Ron, Michael Domrachev, and Alex E Lash. 2002. “Gene Expression Omnibus: NCBI Gene Expression and Hybridization Array Data Repository.” Nucleic Acids Res. 30 (1):207–10. https://www.ncbi.nlm.nih.gov/pubmed/11752295.

Lamb, Justin, Emily D Crawford, David Peck, Joshua W Modell, Irene C Blat, Matthew J Wrobel, Jim Lerner, et al. 2006. “The Connectivity Map: Using Gene-Expression Signatures to Connect Small Molecules, Genes, and Disease.” Science 313 (5795):1929–35. http://dx.doi.org/10.1126/science.1132939.

Peck, David, Emily D Crawford, Kenneth N Ross, Kimberly Stegmaier, Todd R Golub, and Justin Lamb. 2006. “A Method for High-Throughput Gene Expression Signature Analysis.” Genome Biol. 7 (7):R61. http://dx.doi.org/10.1186/gb-2006-7-7-r61.

Subramanian, Aravind, Rajiv Narayan, Steven M Corsello, David D Peck, Ted E Natoli, Xiaodong Lu, Joshua Gould, et al. 2017. “A Next Generation Connectivity Map: L1000 Platform and the First 1,000,000 Profiles.” Cell 171 (6):1437–1452.e17. http://dx.doi.org/10.1016/j.cell.2017.10.049.