R version: R version 4.4.0 beta (2024-04-15 r86425)
Bioconductor version: 3.19
Package version: 2.24.0
In this article we introduce a complete workflow for a typical (Affymetrix) microarray analysis. Data import, preprocessing, differential expression and enrichment analysis are discussed.
The data set used (1) is from a paper studying the differences in gene expression in inflamed and non-inflamed tissue. 14 patients suffering from Ulcerative colitis (UC) and 15 patients with Crohn’s disease (CD) were tested, and from each patient inflamed and non-inflamed colonic mucosa tissue was obtained via a biopsy. This is a typical clinical data set consisting of 58 arrays in total. Our aim is to analyze differential expression (DE) between the tissues. Our results show a substantial overlap with the results of the original paper.
The workflow is wrapped in a package called maEndToEnd.
You can install the package via the BiocManager.
if (!require("BiocManager"))
    install.packages("BiocManager")
BiocManager::install("maEndToEnd")Currently, the workflow has been available since Bioconductor version 3.8, which was first released in October 2018.
Once the workflow package has been successfully installed, we can use
a call to library() in order to load it. This will also load
all the other packages neccessary to run the workflow.
suppressPackageStartupMessages({library("maEndToEnd")})   
   groupGOTerms:    GOBPTerm, GOMFTerm, GOCCTerm environments built.Below, you find a list of packages that are required by the workflow. Some Helper/Styling packages have been commented here, as they are not strictly neccesary to execute the workflow.
#General Bioconductor packages
    library(Biobase)
    library(oligoClasses)
     
#Annotation and data import packages
    library(ArrayExpress)
    library(pd.hugene.1.0.st.v1)
    library(hugene10sttranscriptcluster.db)
     
#Quality control and pre-processing packages
    library(oligo)
    library(arrayQualityMetrics)
     
#Analysis and statistics packages
    library(limma)
    library(topGO)
    library(ReactomePA)
    library(clusterProfiler)
     
#Plotting and color options packages
    library(gplots)
    library(ggplot2)
    library(geneplotter)
    library(RColorBrewer)
    library(pheatmap)
    library(enrichplot)
     
#Formatting/documentation packages
   #library(rmarkdown)
   #library(BiocStyle)
    library(dplyr)
    library(tidyr)
#Helpers:
    library(stringr)
    library(matrixStats)
    library(genefilter)
    library(openxlsx)
   #library(devtools)The first step of the analysis is to download the raw data CEL files. These files are produced by the array scanner software and contain the measured probe intensities. The data we use have been deposited at ArrayExpress and have the accession code E-MTAB-2967.
We will store these files in the directory raw_data_dir which defaults to a cache location determined by R. Alternatively, you can provide any location you like here.
raw_data_dir <- tools::R_user_dir(package = "maEndToEnd", which = "cache")
if (!dir.exists(raw_data_dir)) {
    dir.create(raw_data_dir)
}Each ArrayExpress data set has a landing page summarizing the data set,
and we use the getAEfunction from the ArrayExpress Bioconductor package to obtain
the ftp links to the raw data files (Data from Palmieri et. al. on ArrayEpress).
With the code below, we download the raw data (also including annotation data) from
ArrayExpress (2)
by using the getAE-function. The data are saved in the raw_data_dir
created above. The names of the downloaded files are returned as a list.
anno_AE <- getAE("E-MTAB-2967", path = raw_data_dir, type = "raw")
download.file(url = "https://www.ebi.ac.uk/biostudies/files/E-MTAB-2967/E-MTAB-2967.sdrf.txt",
              destfile = file.path(raw_data_dir, "E-MTAB-2967.sdrf.txt"), 
              mode = "wb")We will now have a closer look at the data we downloaded from ArrayExpress
Each dataset at ArrayExpress is stored according to the MAGE-TAB (MicroArray Gene Expression Tabular) specifications as a collection of tables bundled with the raw data. The MAGE-TAB format specifies up to five different types of files:
Other than the raw data files, the IDF and the SDRF file are important for us. The IDF file contains top level information about the experiment including title, description, submitter contact details and protocols. The SDRF file contains essential information on the experimental samples, e.g. the experimental group(s) they belong to.
Before we move on to the actual raw data import, we will briefly introduce the
ExpressionSet class contained in the Biobase package.
It is commonly used to store microarray data in Bioconductor.
Genomic data can be very complex, usually consisting of a number of different components, e.g. information on the experimental samples, annotation of genomic features measured as well as the experimental data itself. In Bioconductor, the approach is taken that these components should be stored in a single structure to easily manage the data.
The package Biobase contains standardized data structures
to represent genomic data. The ExpressionSet class is designed
to combine several different sources of information (i.e. as contained in the
various MAGE-TAB files) into a single convenient
structure. An ExpressionSet can be manipulated (e.g., subsetted, copied),
and is the input to or output of many Bioconductor functions.
The data in an ExpressionSet consist of:
The ExpressionSet class coordinates all of these data, so that one does not
have to worry about the details. However, one should keep in mind that the rownames of the phenoData have to match the column names of the assay data, while the row names of the assay data have to match the
row names of the featureData.
This is illustrated in Figure 1.
Figure 1: Structure of Bioconductor’s ExpressionSet class
You can use the functions pData and fData to extract
the sample and feature annotation, respectively, from an ExpressionSet.
The function exprs will return the expression data itself as a matrix.
We import the SDRF file with the read.delim function from the raw data
folder in order to obtain the sample annotation.
The sample names are given in the column Array.Data.File of the SDRF data table and will be used as rownames for the SDRF file.
We turn the SDRF table into an AnnotatedDataFrame
from the Biobase package that we
will need later to create an ExpressionSet for our data (3).
sdrf_location <- file.path(raw_data_dir, "E-MTAB-2967.sdrf.txt")
SDRF <- read.delim(sdrf_location)
rownames(SDRF) <- SDRF$Array.Data.File
SDRF <- AnnotatedDataFrame(SDRF)We now create the Expression Set object
raw_data, which contains array data, pheno data
(from the SDRF file) as well as the information
of the chip annotation package used.
The analysis of Affymetrix arrays starts with CEL files. These are the result of the processing of the raw image files using the Affymetrix software and contain estimated probe intensity values. Each CEL file additionally contains some metadata, such as a chip identifier.
We use the function read.celfiles from the oligo
package (4) to import the files:
raw_data <- oligo::read.celfiles(filenames = file.path(raw_data_dir, 
                                                SDRF$Array.Data.File),
                                    verbose = FALSE, phenoData = SDRF)
stopifnot(validObject(raw_data))This automatically creates an ExpressionSet, fills the sections “array data” with the data from the CEL files and uses the correct chip annotation package, in this case pd.hugene.1.0.st.v1 (the chip-type is also stored in the .CEL files).
Furthermore, we specified our AnnotatedDataFrame “SDRF”
created earlier from the SDRF file as phenoData. Thus,
we had to make sure to import the CEL files
in the order that corresponds to the SDRF table — to enforce this,
we used the column Array.Data.File of the SDRF table as the filenames
argument.
Finally, we checked whether the object created is valid (e.g. sample names match between the different tables).
We now have a first look on the raw data.
The pData function of the Biobase package directly accesses the phenoData in the ExpressionSet raw_data. With the head() function, we can view the
first six lines of the table. We have a look at the columns included
and retain only those columns that are related to the experimental
factors of interest.
head(Biobase::pData(raw_data))               Source.Name Characteristics.individual. Characteristics.organism.
   164_I_.CEL        164_I                         164              Homo sapiens
   164_II.CEL       164_II                         164              Homo sapiens
   183_I.CEL         183_I                         183              Homo sapiens
   183_II.CEL       183_II                         183              Homo sapiens
   2114_I.CEL       2114_I                        2114              Homo sapiens
   2114_II.CEL     2114_II                        2114              Homo sapiens
               Characteristics.disease. Characteristics.organism.part.
   164_I_.CEL           Crohn's disease                          colon
   164_II.CEL           Crohn's disease                          colon
   183_I.CEL            Crohn's disease                          colon
   183_II.CEL           Crohn's disease                          colon
   2114_I.CEL           Crohn's disease                          colon
   2114_II.CEL          Crohn's disease                          colon
                Characteristics.phenotype. Material.Type Protocol.REF
   164_I_.CEL  non-inflamed colonic mucosa organism part P-MTAB-41361
   164_II.CEL      inflamed colonic mucosa organism part P-MTAB-41361
   183_I.CEL   non-inflamed colonic mucosa organism part P-MTAB-41361
   183_II.CEL      inflamed colonic mucosa organism part P-MTAB-41361
   2114_I.CEL  non-inflamed colonic mucosa organism part P-MTAB-41361
   2114_II.CEL     inflamed colonic mucosa organism part P-MTAB-41361
               Protocol.REF.1 Extract.Name Protocol.REF.2 Labeled.Extract.Name
   164_I_.CEL    P-MTAB-41363        164_I   P-MTAB-41364         164_I:Biotin
   164_II.CEL    P-MTAB-41363       164_II   P-MTAB-41364        164_II:Biotin
   183_I.CEL     P-MTAB-41363        183_I   P-MTAB-41364         183_I:Biotin
   183_II.CEL    P-MTAB-41363       183_II   P-MTAB-41364        183_II:Biotin
   2114_I.CEL    P-MTAB-41363       2114_I   P-MTAB-41364        2114_I:Biotin
   2114_II.CEL   P-MTAB-41363      2114_II   P-MTAB-41364       2114_II:Biotin
                 Label Protocol.REF.3 Assay.Name Technology.Type Array.Design.REF
   164_I_.CEL  biotin    P-MTAB-41366     164_I_     array assay       A-AFFY-141
   164_II.CEL  biotin    P-MTAB-41366     164_II     array assay       A-AFFY-141
   183_I.CEL   biotin    P-MTAB-41366      183_I     array assay       A-AFFY-141
   183_II.CEL  biotin    P-MTAB-41366     183_II     array assay       A-AFFY-141
   2114_I.CEL  biotin    P-MTAB-41366     2114_I     array assay       A-AFFY-141
   2114_II.CEL biotin    P-MTAB-41366    2114_II     array assay       A-AFFY-141
               Term.Source.REF Protocol.REF.4 Array.Data.File
   164_I_.CEL     ArrayExpress   P-MTAB-41367      164_I_.CEL
   164_II.CEL     ArrayExpress   P-MTAB-41367      164_II.CEL
   183_I.CEL      ArrayExpress   P-MTAB-41367       183_I.CEL
   183_II.CEL     ArrayExpress   P-MTAB-41367      183_II.CEL
   2114_I.CEL     ArrayExpress   P-MTAB-41367      2114_I.CEL
   2114_II.CEL    ArrayExpress   P-MTAB-41367     2114_II.CEL
                                                                                   Comment..ArrayExpress.FTP.file.
   164_I_.CEL  ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
   164_II.CEL  ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
   183_I.CEL   ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
   183_II.CEL  ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
   2114_I.CEL  ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
   2114_II.CEL ftp://ftp.ebi.ac.uk/pub/databases/microarray/data/experiment/MTAB/E-MTAB-2967/E-MTAB-2967.raw.1.zip
               Factor.Value.disease.     Factor.Value.phenotype.
   164_I_.CEL        Crohn's disease non-inflamed colonic mucosa
   164_II.CEL        Crohn's disease     inflamed colonic mucosa
   183_I.CEL         Crohn's disease non-inflamed colonic mucosa
   183_II.CEL        Crohn's disease     inflamed colonic mucosa
   2114_I.CEL        Crohn's disease non-inflamed colonic mucosa
   2114_II.CEL       Crohn's disease     inflamed colonic mucosaThe columns of interest for us are the following:
We now subselect the corresponding columns:
Biobase::pData(raw_data) <- Biobase::pData(raw_data)[, c("Source.Name",
                                     "Characteristics.individual.",
                                     "Factor.Value.disease.",
                                     "Factor.Value.phenotype.")]The first step after the initial data import is the quality control of the data.
Here we check for outliers and try to see whether the data clusters as expected,
e.g. by the experimental conditions. The expression intensity values are in the
assayData sub-object “exprs” and can be accessed by the exprs(raw_data)
function. The rows represent the microarray probes, i.e. 
the single DNA locations on the chip, while the columns
represent one microarray, i.e. a sample of inflamed
and non-inflamed tissue of every patient, respectively.
Biobase::exprs(raw_data)[1:5, 1:5]     164_I_.CEL 164_II.CEL 183_I.CEL 183_II.CEL 2114_I.CEL
   1       4496       5310      4492       4511       2872
   2        181        280       137        101         91
   3       4556       5104      4379       4608       2972
   4        167        217        99         79         82
   5         89        110        69         58         47For quality control, we take the log2 of Biobase::exprs(raw_data),
as expression data is commonly analyzed on a logarithmic scale.
We then perform a principal component analysis (PCA) and plot it (Figure 2). Every point in the plot represents one sample, with the colour indicating the mucosa type (inflamed vs non-inflamed) and the shape indicating the disease (UC or CD).
exp_raw <- log2(Biobase::exprs(raw_data))
PCA_raw <- prcomp(t(exp_raw), scale. = FALSE)
percentVar <- round(100*PCA_raw$sdev^2/sum(PCA_raw$sdev^2),1)
sd_ratio <- sqrt(percentVar[2] / percentVar[1])
dataGG <- data.frame(PC1 = PCA_raw$x[,1], PC2 = PCA_raw$x[,2],
                    Disease = pData(raw_data)$Factor.Value.disease.,
                    Phenotype = pData(raw_data)$Factor.Value.phenotype.,
                    Individual = pData(raw_data)$Characteristics.individual.)
ggplot(dataGG, aes(PC1, PC2)) +
      geom_point(aes(shape = Disease, colour = Phenotype)) +
  ggtitle("PCA plot of the log-transformed raw expression data") +
  xlab(paste0("PC1, VarExp: ", percentVar[1], "%")) +
  ylab(paste0("PC2, VarExp: ", percentVar[2], "%")) +
  theme(plot.title = element_text(hjust = 0.5))+
  coord_fixed(ratio = sd_ratio) +
  scale_shape_manual(values = c(4,15)) + 
  scale_color_manual(values = c("darkorange2", "dodgerblue4"))
Figure 2: PCA plot of the log–transformed raw expression data
The PCA plot (Figure 2, performed on the log-intensity scale) of the raw data shows that the first principal component differentiates between the diseases. This means that the disease type is a major driver of gene expression differences. This might hinder our analysis, as we want to analyze the differential expression between inflamed and non-inflamed tissues, independently of the disease a person suffers from.
We also represent the probe intensities via a boxplot graph with one box per
individual microarray. (Figure 3). Note that the
oligo::boxplot function, i.e. the boxplot function of the oligo package, can take
expression sets as argument. It accesses the expression data and performs a
log2-transformation by default. We therefore can use raw_data as argument here.
oligo::boxplot(raw_data, target = "core", 
               main = "Boxplot of log2-intensitites for the raw data")
Figure 3: Intensity boxplots of the log2–transformed raw data
When looking at the boxplot (Figure 3), we see that the intensity distributions of the individual arrays are quite different, indicating the need for an appropriate normalization, which we will discuss next.
Until now, we have only performed a very basic quality control; more elaborate quality control plots are available in the package arrayQualityMetrics (5). The package produces an html report, containing the quality control plots together with a description of their aims and an identification of possible outliers. We do not discuss this tool in detail here, but simply provide the code below, which creates a report for our raw data.
arrayQualityMetrics(expressionset = raw_data,
    outdir = tempdir(),
    force = TRUE, do.logtransform = TRUE,
    intgroup = c("Factor.Value.disease.", "Factor.Value.phenotype."))After the initial import and quality assessment, the next step in processing of microarray data is background adjustment. This is essential because a proportion of the measured probe intensities are due to non-specific hybridization and the noise in the optical detection system. Therefore, observed intensities need to be adjusted to give accurate measurements of specific hybridization.
Normalization across arrays is needed in order to be able to compare measurements from different array hybridizations due to many obscuring sources of variation. These include different efficiencies of reverse transcription, labeling or hybridization reactions, physical problems with the arrays, reagent batch effects, and laboratory conditions.
After normalization, summarization is necessary to be done because on the Affymetrix platform, transcripts are represented by multiple probes, that is multiple locations on the array. For each gene, the background-adjusted and normalized intensities of all probes need to be summarized into one quantity that estimates an amount proportional to the amount of RNA transcript.
After the summarization step, the summarized data can be annotated with various information, e.g. gene symbols and ENSEMBL gene identifiers. There is an annotation database available from Bioconductor for our platform, namely the package hugene10sttranscriptcluster.db.
You can view its content like this:
head(ls("package:hugene10sttranscriptcluster.db"))   [1] "hugene10sttranscriptcluster"           
   [2] "hugene10sttranscriptcluster.db"        
   [3] "hugene10sttranscriptclusterACCNUM"     
   [4] "hugene10sttranscriptclusterALIAS2PROBE"
   [5] "hugene10sttranscriptclusterCHR"        
   [6] "hugene10sttranscriptclusterCHRLENGTHS"Additional information is available from the reference manual of the package. Essentially, the package provides a mapping from the transcript cluster identifiers to the various annotation data.
Traditionally, Affymetrix arrays (the so-called 3’ IVT arrays) were probeset based: a certain fixed group of probes were part of a probeset which represented a certain gene or transcript (note however, that a gene can be represented by multiple probesets).
The more recent “Gene” and “Exon” Affymetrix arrays are exon based and hence there are two levels of summarization to get to the gene level. The “probeset” summarization leads to the exon level. The gene / transcript level is given by “transcript clusters”. Hence, the appropriate annotation package for our chip type is called hugene10sttranscriptcluster.db.
“Gene” arrays were created as affordable versions of the “Exon” arrays, by only taking the “good” probes from the Exon array. Initially on the Exon array, at least four probes were part of one “Exon”. With the thinned out “Gene” array, many probesets were made up of three or fewer probes. This is visualized in Figure 4: Single probesets are indicated by single colours; probesets representing one gene are indicated by a colour shade: e.g., all yellow probes belong to one Exon, and all yellow, orange and red probesets belong to one gene:
Figure 4: Visualization of the difference between “Exon” type array (left) and “Gene” type array (right)
On the left side, we see plenty of probes for each Exon / probeset (i.e. each colour): therefore, a summarization on the probeset / exon level makes sense. In the gene type array, however, only a small proportion of the original probes per probeset is included. Thus, a summarization on the probeset / exon level is not recommended for “Gene” arrays but nonetheless possible by using the hugene10stprobeset.db annotation package.
Note that furthermore, there are also no longer designated match/mismatch probes present on “Gene” and “Exon” type chips. The mismatch probe was initially intended as base-level for background correction, but hasn’t prevailed due to more elaborate background correction techniques that do not require a mismatch probe.
The package oligo allows us to perform background correction, normalization and summarization in one single step using a deconvolution method for background correction, quantile normalization and the RMA (robust multichip average) algorithm for summarization.
This series of steps as a whole is commonly referred to as RMA algorithm, although strictly speaking RMA is merely a summarization method (6–8).
Before calibrating and evaluating the data, we want to perform another quality control procedure, namely Relative Log Expression (RLE), as described in the article by Gandolfo et al (9). To this end, we first perform an RMA without prior normalization:
palmieri_eset <- oligo::rma(raw_data, target = "core", normalize = FALSE)   Background correcting
   Calculating ExpressionFurther details on the RMA algorithm will be provided after RLE analysis, when the “full” RMA is carried out, including normalization.
The RLE is performed by calculating the median log2 intensity of every transcript across all arrays.
We do this by calculating the row medians of exprs(palmieri_eset),
as the transcripts are represented by the rows and the single microarrays
by the columns.
Note that we do not have to apply the log2 manually, as the output data of the RMA function is in log2 scale by default.
We then substract this transcript median intensity from every
transcript intensity via the sweep function.
We finally reshape the data into a format in which we can use to create a boxplot for each array, as before:
row_medians_assayData <- 
  Biobase::rowMedians(as.matrix(Biobase::exprs(palmieri_eset)))
RLE_data <- sweep(Biobase::exprs(palmieri_eset), 1, row_medians_assayData)
RLE_data <- as.data.frame(RLE_data)
RLE_data_gathered <- 
  tidyr::gather(RLE_data, patient_array, log2_expression_deviation)
ggplot2::ggplot(RLE_data_gathered, aes(patient_array,
                                       log2_expression_deviation)) + 
  geom_boxplot(outlier.shape = NA) + 
  ylim(c(-2, 2)) + 
  theme(axis.text.x = element_text(colour = "aquamarine4", 
                                  angle = 60, size = 6.5, hjust = 1 ,
                                  face = "bold"))
Figure 5: Boxplot for the RLE values
Note that the y-axis now displays for each microarray the deviation of expression intensity from the median expression of the respective single transcripts across arrays.
Boxes with a larger extension therefore indicate an unusually high deviation from the median in a lot of transcripts, suggesting that these arrays are different from most of the others in some way.
Boxes that are shifted in y-direction indicate a systematically higher or lower expression of the majority of transcripts in comparison to most of the other arrays. This could be caused by quality issues or batch effects.
Therefore, if shape and median of a given box varies too much from the bulk, they should be inspected and potentially removed.
By inspecting the boxplot in Figure 5, five arrays could be considered as outliers: 2826_I, 2826_II, 3262_II, 3302_II and 3332_II are negatively y-shifted.
We will keep these samples in mind for heatmap cluster analysis later on in the workflow. Arrays that are confirmed to be outliers by heatmap analysis could be removed for subsequent analysis.
Now, we can apply the full RMA algorithm to our data in order to background-correct, normalize and summarize:
palmieri_eset_norm <- oligo::rma(raw_data, target = "core")   Background correcting
   Normalizing
   Calculating ExpressionThe parameter target defines the degree of summarization, the
default option of which is “core”, using transcript clusters containing
“safely” annotated genes. For summarization on the exon level (not recommended for Gene
arrays), one can use “probeset” as the target option.
Although other methods for background correction and normalization exist,
RMA is usually a good default choice.
RMA shares information across arrays and
uses the versatile quantile normalization method that
will make the array intensity distributions match. However, it is preferable
to apply it only after outliers have been removed.
The quantile normalization algorithm used by RMA
works by replacing values by the average of identically
ranked (within a single chip) values across arrays. A more detailed
description can be found on the Wikipedia page.
An alternative to quantile normalization would be the vsn algorithm,that performs background correction and normalization by robustly shifting and scaling intensity values within arrays before log-transforming them. This is less “severe” than quantile normalization (10).
A generic model for the value of the intensity \(Y\) of a single probe on a microarray is given by
\[ Y = B + \alpha \cdot S \]
where B is a random quantity due to background noise, usually composed of optical effects and non-specific binding, \(\alpha\) is a gain factor, and \(S\) is the amount of measured specific binding. The signal \(S\) is considered a random variable as well and accounts for measurement error and probe effects. The measurement error is typically assumed to be multiplicative so we can write:
\[ \log(S) = \theta + \varphi + \varepsilon \]
Here \(\theta\) represents the logarithm of the true abundance, \(\varphi\) is a probe-specific effect, and \(\varepsilon\) accounts for the nonspecific error. This is the additive-multiplicative-error model for microarray data used by RMA and also the vsn algorithm (10). The algorithms differ in the way that \(B\) is removed and an estimate of \(\theta\) is obtained.
We now produce a clustering heatmap and another PCA plot using the calibrated data.
PCA analysis
First, we perform a PCA analysis of the calibrated data analogously to the one with the raw data:
exp_palmieri <- Biobase::exprs(palmieri_eset_norm)
PCA <- prcomp(t(exp_palmieri), scale = FALSE)
percentVar <- round(100*PCA$sdev^2/sum(PCA$sdev^2),1)
sd_ratio <- sqrt(percentVar[2] / percentVar[1])
dataGG <- data.frame(PC1 = PCA$x[,1], PC2 = PCA$x[,2],
                    Disease = 
                     Biobase::pData(palmieri_eset_norm)$Factor.Value.disease.,
                    Phenotype = 
                     Biobase::pData(palmieri_eset_norm)$Factor.Value.phenotype.)
ggplot(dataGG, aes(PC1, PC2)) +
      geom_point(aes(shape = Disease, colour = Phenotype)) +
  ggtitle("PCA plot of the calibrated, summarized data") +
  xlab(paste0("PC1, VarExp: ", percentVar[1], "%")) +
  ylab(paste0("PC2, VarExp: ", percentVar[2], "%")) +
  theme(plot.title = element_text(hjust = 0.5)) +
  coord_fixed(ratio = sd_ratio) +
  scale_shape_manual(values = c(4,15)) + 
  scale_color_manual(values = c("darkorange2", "dodgerblue4"))
Figure 6: PCA plot of the calibrated, summarized data
In comparison to the first PCA analysis before RMA (Figure 2), we see that now the first principal component separates between the tissues types (Figure 6). This indicates that now differential expression between the tissue types is the dominant source of variation. Note that the second principal component separates the diseases.
Heatmap clustering analysis
We want to plot a heatmap with the sample-to-sample distances with the sample names as row-names. Also, we want to see how well the samples cluster for phenotype (inflamed and non-inflamed tissue) and disease (UC and CD), respectively. We use row annotation for that: it means that these features get a colour code and will be displayed to the left of each row.
phenotype_names <- ifelse(str_detect(pData
                                    (palmieri_eset_norm)$Factor.Value.phenotype.,
                             "non"), "non_infl.", "infl.")
disease_names <- ifelse(str_detect(pData
                                    (palmieri_eset_norm)$Factor.Value.disease.,
                             "Crohn"), "CD", "UC")
annotation_for_heatmap <- 
  data.frame(Phenotype = phenotype_names,  Disease = disease_names)
row.names(annotation_for_heatmap) <- row.names(pData(palmieri_eset_norm))In order to map the sample-to-sample distances, we first compute the distances
using the dist function. We need to transpose the expression values since
the function computes the distances between the rows (i.e. genes in our case) by
default. The default distance is the Euclidean one. However this can be
changed and we choose the Manhattan distance here (it uses absolute
distances along rectangular paths instead of squared distances
of the direct path), as it is more robust.
We set the diagonal of the distance matrix to NA in order
to increase the contrast of the color coding. Those diagonal entries
do not contain information since the distance of a sample
to itself is always equal to zero.
dists <- as.matrix(dist(t(exp_palmieri), method = "manhattan"))
rownames(dists) <- row.names(pData(palmieri_eset_norm))
hmcol <- rev(colorRampPalette(RColorBrewer::brewer.pal(9, "YlOrRd"))(255))
colnames(dists) <- NULL
diag(dists) <- NA
ann_colors <- list(
  Phenotype = c(non_infl. = "chartreuse4", infl. = "burlywood3"),
  Disease = c(CD = "blue4", UC = "cadetblue2")
                   )
pheatmap(dists, col = (hmcol), 
         annotation_row = annotation_for_heatmap,
         annotation_colors = ann_colors,
         legend = TRUE, 
         treeheight_row = 0,
         legend_breaks = c(min(dists, na.rm = TRUE), 
                         max(dists, na.rm = TRUE)), 
         legend_labels = (c("small distance", "large distance")),
         main = "Clustering heatmap for the calibrated samples")
Figure 7: Heatmap of the sample-to-sample distances
On the heatmap plot (Figure 7) we also see that the samples do not cluster strongly by tissue, confirming the impression from the PCA plot (Figure 6) that the separation between the tissues is not perfect. The yellow stripes in the heatmap might correspond to outliers that could potentially be removed: the ones that could be flagged here are 2826_II, 3262_II, 3271_I, 2978_II and 3332_II. 2826_II, 3262_II and 3332_II were found to be outliers in both RLE and heatmap analysis and could therefore potentially be removed; however, in order to stay as close as possible to the original paper, we continue with the complete set of samples. Note again that more elaborate metrics to identify and remove outliers are provided by the arrayQualityMetrics package.
We now filter out lowly expressed genes. Microarray data commonly show a large numberof probes in the background intensity range. These probes also do not change much across arrays. Hence they combine a low variance with a low intensity. Thus, they could end up being detected as differentially expressed although they are barely above the “detection” limit and are not very informative in general.
We will perform a “soft” intensity based filtering here, since this is recommended by the limma (11,12) user guide (a package we will use below for the differential expression analysis).
However, note that a variance based filter might exclude a similar set of
probes in practice.
For intensity-based filtering, we calculate the row-wise medians from the
expression data, as they represent the transcript medians, and assign them to
palmieri_medians. From this we create a histogram:
palmieri_medians <- rowMedians(Biobase::exprs(palmieri_eset_norm))
hist_res <- hist(palmieri_medians, 100, col = "cornsilk1", freq = FALSE, 
            main = "Histogram of the median intensities", 
            border = "antiquewhite4",
            xlab = "Median intensities")
Figure 8: Histogram of the median intensities per gene
In the histogram of the gene-wise medians (Figure
8), we
can clearly see an enrichment of low medians on the left hand side.
These represent the genes we want to filter. In order to infer a cutoff from the
data, we inspect the histogram: We visually set a cutoff line man_threshold
to the left of the histogram peak in order not to exclude too many genes. In our
example, we choose a threshold of 4. We plot the same histogram as before
and add the threshold line with the abline() function:
man_threshold <- 4
hist_res <- hist(palmieri_medians, 100, col = "cornsilk", freq = FALSE, 
            main = "Histogram of the median intensities",
            border = "antiquewhite4",
            xlab = "Median intensities")
abline(v = man_threshold, col = "coral4", lwd = 2)
Figure 9: Histogram of the median intensities per gene with manual intensity filtering threshold (red line)
Transcripts that do not have intensities larger than the threshold in at least as many arrays as the smallest experimental group are excluded.
In order to do so, we first have to get a list with the number
of samples (=arrays) (no_of_samples) in the experimental groups:
no_of_samples <- 
  table(paste0(pData(palmieri_eset_norm)$Factor.Value.disease., "_", 
                  pData(palmieri_eset_norm)$Factor.Value.phenotype.))
no_of_samples    
          Crohn's disease_inflamed colonic mucosa 
                                               15 
      Crohn's disease_non-inflamed colonic mucosa 
                                               15 
       ulcerative colitis_inflamed colonic mucosa 
                                               14 
   ulcerative colitis_non-inflamed colonic mucosa 
                                               14We now filter out all transcripts
that do not have intensities greater than the threshold in at least
as many arrays as the smallest experimental group (14)
which we define as samples_cutoff.
A function idx_man_threshold is then applied to each row, i.e. 
to each transcript across all arrays.
It evaluates whether the number of arrays where the median intensity passes the
threshold (sum(x > man_threshold)) is greater than the samples_cutoff and
returns TRUE or FALSE for each row, i.e. each transcript.
We then create a table of idx_man_threshold to summarize the results and
get an overview over how many genes are filtered out.
In the last step, we subset our expression set to palmieri_manfiltered
and keep the TRUE elements of idx_man_threshold.
samples_cutoff <- min(no_of_samples)
idx_man_threshold <- apply(Biobase::exprs(palmieri_eset_norm), 1,
                           function(x){
                          sum(x > man_threshold) >= samples_cutoff})
                          table(idx_man_threshold)   idx_man_threshold
   FALSE  TRUE 
   10493 22804palmieri_manfiltered <- subset(palmieri_eset_norm, idx_man_threshold)Before we continue with the linear models for microarrays and differential expression, we first add “feature data”, i.e. annotation information to the transcript cluster identifiers stored in the featureData of our ExpressionSet:
anno_palmieri <- AnnotationDbi::select(hugene10sttranscriptcluster.db,
                                  keys = (featureNames(palmieri_manfiltered)),
                                  columns = c("SYMBOL", "GENENAME"),
                                  keytype = "PROBEID")
anno_palmieri <- subset(anno_palmieri, !is.na(SYMBOL))We used the function select from AnnotationDbi to query
the gene symbols and associated short descriptions for the transcript clusters.
For each cluster, we added the gene symbol (SYMBOL) and a short description
of the gene the cluster represents (GENENAME).
In a second step, we filtered out the probes that do not map to a gene, i.e. that do not have a gene symbol assigned.
Many transcript-cluster identifiers will map to multiple gene symbols, i.e. they can’t be unambigously assigned.
We compute a summary table in the code below to see how many there are:
anno_grouped <- group_by(anno_palmieri, PROBEID)
anno_summarized <- 
  dplyr::summarize(anno_grouped, no_of_matches = n_distinct(SYMBOL))
head(anno_summarized)   # A tibble: 6 × 2
     PROBEID no_of_matches
     <chr>           <int>
   1 7896742             2
   2 7896746             1
   3 7896754             1
   4 7896759             1
   5 7896761             1
   6 7896779             1anno_filtered <- filter(anno_summarized, no_of_matches > 1)
head(anno_filtered)   # A tibble: 6 × 2
     PROBEID no_of_matches
     <chr>           <int>
   1 7896742             2
   2 7896937             3
   3 7896961             3
   4 7897006             2
   5 7897632             3
   6 7897774             2probe_stats <- anno_filtered 
nrow(probe_stats)   [1] 1713First, we grouped anno_palmieri by their PROBEID; that way,
the subsequent operations are not carried through for each single row,
but for each group, i.e. each PROBEID.
We then summarized the groups and indicate the number of different genes
assigned to a transcript cluster in the column no_of_matches.
Finally, we filtered for PROBEIDs with multiple matches,
i.e. no_of_matches > 1.
With dim(probe_stats), we could see how many
probes have been mapped to multiple genes.
We have close to 2000 transcript clusters that map to multiple gene symbols. It is difficult to decide which mapping is “correct”. Therefore, we exclude these transcript clusters.
We want to remove those probe IDs that match the ones in probe_stats,
as those are the probes with multiple mappings. We assign these IDs
to the variable ids_to_exclude. Then, we generate palmieri_final,
an expression set without the ids_to_exclude.
ids_to_exlude <- (featureNames(palmieri_manfiltered) %in% probe_stats$PROBEID)
table(ids_to_exlude)   ids_to_exlude
   FALSE  TRUE 
   21091  1713palmieri_final <- subset(palmieri_manfiltered, !ids_to_exlude)
validObject(palmieri_final)   [1] TRUEAs we have just excluded probe IDs from the assay data, we now have
to also exclude them from the feature data anno_palmieri:
head(anno_palmieri)        PROBEID     SYMBOL
   2657 7896742     PCMTD2
   2658 7896742  LINC01881
   2659 7896746   MTND1P23
   2663 7896754 SEPTIN7P13
   2664 7896759  LINC01128
   2665 7896761     SAMD11
                                                                            GENENAME
   2657 protein-L-isoaspartate (D-aspartate) O-methyltransferase domain containing 2
   2658                                  long intergenic non-protein coding RNA 1881
   2659                                                         MT-ND1 pseudogene 23
   2663                                                       septin 7 pseudogene 13
   2664                                  long intergenic non-protein coding RNA 1128
   2665                                     sterile alpha motif domain containing 11Recall that fData enables us to access the feature data of an expression set.
Until now, no feature data whatsoever is stored in the fData(palmieri_final).
Only the row names are the row names of the assay data by default, which are the
PROBEIDs of the transcripts.
Therefore, we generate a column PROBEID in fData(palmieri_final) and
assign the row names of fData(palmieri_final) to it:
fData(palmieri_final)$PROBEID <- rownames(fData(palmieri_final))Then, we left-join fData(palmieri_final)with anno_palmieri, which already contains
the columns “SYMBOL” and “GENENAME”.
A left-join keeps the rows and columns of the first argument
and adds the corresponding column entries of the second argument:
fData(palmieri_final) <- left_join(fData(palmieri_final), anno_palmieri)   Joining with `by = join_by(PROBEID)`# restore rownames after left_join
rownames(fData(palmieri_final)) <- fData(palmieri_final)$PROBEID 
    
validObject(palmieri_final)   [1] TRUEBy left-joining with anno_palmieri, we thus add the “SYMBOL” and “GENENAME”
columns from anno_palmieri for only the PROBEIDs that are in
fData(palmieri_final) and thus get the feature Data for the filtered probes.
Alternatively, one can re-map the probes of the array to a current annotation. A workflow to do this for Illumina arrays is given in (13). Essentially, the individual probe sequences are re-aligned to an in-silico “exome” that consists of all annotated transcript exons.
In any case, the package pdInfoBuilder can be used to build customannotation packages for use with oligo. In order to do this, PGF / CLF files (called “Library files” on the Affymetrix website) as well as the probeset annotations are required. The probesets typically represent small stretches of the genome (such as a single exon) and multiple probesets are then used to form a transcript cluster.
The CLF file contains information about the location of individual probes on the array. The PGF file then contains the individual probe sequences and shows the probeset they belong to. Finally, the probeset annotation .csv then contains information about which probesets are used in which transcript cluster. Commonly, multiple probesets are used in one transcript cluster and some probesets are contained in multiple transcript clusters.
In order to analyse which genes are differentially expressed between inflamed and non-inflamed tissue, we have to fit a linear model to our expression data. Linear models are the “workhorse” for the analysis of experimental data. They can be used to analyse almost arbitrarily complex designs, however they also take a while to learn and understand and a thorough description is beyond the scope of this workflow.
Mike Love’s and Rafael Irizzary’s “Biomedical Data Science” (14) is a very good resource, especially the section on interactions and contrasts. It might also be helpful to learn some linear algebra to better understand the concepts here. The Khan Academy offers helpful (and free) online courses.
We now apply linear models to microarrays. Specifically, we discuss how to use the limma package for differential expression analysis. The package is designed to analyze complex experiments involving comparisons between many experimental groups simultaneously while remaining reasonably easy to use for simple experiments. The main idea is to fit a linear model to the expression data for each gene. Empirical Bayes and other methods are used to borrow information across genes for the residual variance estimation leading to “moderated” \(t\)-statistics, and stabilizing the analysis for experiments with just a small number of arrays (12). Conceptually, the final per gene variance is a mix of a prior variance and the per gene variance.
Typical experimental designs are disussed in chapter 9 of limma “User Guide”, which can be found on the Bioconductor landing page of limma.
In the following, we use appropriate design and contrast matrices for our linear models and fit a linear model to each gene separately.
For the subsequent linear modelling of the data, we introduce the abbreviations “UC” and “CD” for the disease types, and “non_infl.” and “infl.” for the phenotypes, respectively:
individual <- 
  as.character(Biobase::pData(palmieri_final)$Characteristics.individual.)
tissue <- str_replace_all(Biobase::pData(palmieri_final)$Factor.Value.phenotype.,
                  " ", "_")
tissue <- ifelse(tissue == "non-inflamed_colonic_mucosa",
                 "nI", "I")
disease <- 
  str_replace_all(Biobase::pData(palmieri_final)$Factor.Value.disease.,
                  " ", "_")
disease <- 
  ifelse(str_detect(Biobase::pData(palmieri_final)$Factor.Value.disease., 
                    "Crohn"), "CD", "UC")The original paper is interested in changes in transcription that occur between inflamed and adjacent non-inflamed mucosal areas of the colon. This is studied in both inflammatory bowel disease types.
For building our linear model, we have to think about which experimental variables we want to consider. As we want to find differential expression between non-inflamed and inflamed tissue, in principle, those are the only two variables we would have to consider.
However, since we have two arrays per individual patient, we have a “Paired Samples” design (see section 9.4 of the limma user guide). This means that the samples might be biased by the person they come from. Whenever a feature in an experimental setup is expected to have a systematic influence on the result, blocking factors on these features should be introduced.
Thus, the first factor we need is a blocking factor for the individuals that will absorb differences in expression between them. Therefore, we block on patients, which means that the patient IDs become variables of the linear model.
Then we create factors that give us the grouping for the tissue types (non-inflamed and inflamed).
Finally, we create two design matrices, one for each of the two diseases as we will analyze them separately in order to follow the analysis strategy of the original paper closely (one could also fit a joint model to the complete data set; however, the two diseases might show very different characteristics so that a joint fit might not be appropriate).
i_CD <- individual[disease == "CD"]
design_palmieri_CD <- model.matrix(~ 0 + tissue[disease == "CD"] + i_CD)
colnames(design_palmieri_CD)[1:2] <- c("I", "nI")
rownames(design_palmieri_CD) <- i_CD 
i_UC <- individual[disease == "UC"]
design_palmieri_UC <- model.matrix(~ 0 + tissue[disease == "UC"] + i_UC )
colnames(design_palmieri_UC)[1:2] <- c("I", "nI")
rownames(design_palmieri_UC) <- i_UC We can inspect the design matrices:
head(design_palmieri_CD[, 1:6])        I nI i_CD183 i_CD2114 i_CD2209 i_CD2255
   164  0  1       0        0        0        0
   164  1  0       0        0        0        0
   183  0  1       1        0        0        0
   183  1  0       1        0        0        0
   2114 0  1       0        1        0        0
   2114 1  0       0        1        0        0head(design_palmieri_UC[, 1:6])        I nI i_UC2424 i_UC2987 i_UC2992 i_UC2995
   2400 0  1        0        0        0        0
   2400 1  0        0        0        0        0
   2424 0  1        1        0        0        0
   2424 1  0        1        0        0        0
   2987 0  1        0        1        0        0
   2987 1  0        0        1        0        0In the design matrix, the rows represent the patient array,
and the columns are the variables we include in our linear model. The variables
correspond to our blocking factors: there are two for non-inflamed and inflamed
tissue, respectively, and one for each patient. “i_UC2424” for example is the
blocking variable of patient 2424; UC stands for the disease the patient is
suffering from.
For example, the first two rows of the design matrix design_palmieri_CD
correspond to the two arrays for individual “164”.
The design matrix entries are 0 or 1 and thus tell us which variables are “active” for which sample:
It can be seen as a switch that turns the single variables “on” (with a 1 at the corresponding position) or “off” (with a 0 at the corresponding position) for each row, i. e. each patient array. If we take a closer look at the single rows, we see that for each sample, there are two “ones” assigned: one to one of the variables “nI” or “I” corresponding to the tissue the array came from, and one to the corresponding patient-specific blocking variable.
Note that in the linear model, individual 164 serves as the baseline for all other individuals and thus isn’t included in the sample variables.
Note that instead of blocking on individuals, it would also be possible to use a
“mixed model” approach with the duplicateCorrelation()function from the
limma package.
It has advantages over the “fixed patient effect model” presented here
in terms of applicability to more complicated experimental designs,
where we want to perform comparisons both within and between the patients (e.g. 
comparing between the two diseases; “split-plot-designs”).
More information on it can be found in the limma User’s Guide (section 17.3.6). However, the above explained is more intuitive and is therefore used here.
Before heading on to find all differentially expressed genes for both diseases, we will first have a look at how this works in principle for one gene. We will fit the linear model for one gene and run a t-test in order to see whether the gene is differentially expressed or not.
For linear model fitting and subsequent testing for differential expression by t-test, we will pick the gene with the PROBEID 8164535. It has the gene symbol CRAT and will be named as such in the following code.
Illustration of the fitted linear model on the CRAT gene
Before fitting the linear model, we have a look at the expression intensities of this gene for each patient in non-inflamed and inflamed tissue, respectively:
tissue_CD <- tissue[disease == "CD"]
crat_expr <- Biobase::exprs(palmieri_final)["8164535", disease == "CD"]
crat_data <- as.data.frame(crat_expr)
colnames(crat_data)[1] <- "org_value"
crat_data <- mutate(crat_data, individual = i_CD, tissue_CD)
crat_data$tissue_CD <- factor(crat_data$tissue_CD, levels = c("nI", "I"))
ggplot(data = crat_data, aes(x = tissue_CD, y = org_value, 
                             group = individual, color = individual)) +
      geom_line() +
      ggtitle("Expression changes for the CRAT gene")
Figure 10: Visualization of expression changes
We see that overall, this gene is expressed less in inflamed tissue (Figure 10). We also see that the absolute expression intensities vary greatly between patients. However, we have already taken care of this problem by introducing blocking factors based on the individuals, which allows us to compare the tissues within each individual as represented by the single lines.
If we had not blocked for individuals, the linear model would treat them interchangably and a graphical depiction would only include a single line. Since the individuals haver very different baseline expression levels, this would lead to a very high variance of the estimated fold changes.
We now compute the variable coefficients by fitting a linear model.
We get a vector crat_coef with one entry for each variable.
crat_coef <- lmFit(palmieri_final[,disease == "CD"],
                design = design_palmieri_CD)$coefficients["8164535",]
crat_coef          I       nI  i_CD183 i_CD2114 i_CD2209 i_CD2255  i_CD255 i_CD2826 
    6.76669  7.19173  0.12382 -0.22145  0.55759 -0.39905  0.29204 -1.07285 
   i_CD2853 i_CD2978  i_CD321 i_CD3262 i_CD3266 i_CD3271 i_CD3302 i_CD3332 
   -0.78285 -0.11633  0.01692 -0.62480 -0.46209 -0.61732 -0.30257 -0.09709In order to now obtain the expression values fitted by the model,
we multiply the design matrix with this vector of coefficients crat_coef:
crat_fitted <- design_palmieri_CD %*% crat_coef
rownames(crat_fitted) <- names(crat_expr)
colnames(crat_fitted) <- "fitted_value"
crat_fitted               fitted_value
   164_I_.CEL         7.192
   164_II.CEL         6.767
   183_I.CEL          7.316
   183_II.CEL         6.891
   2114_I.CEL         6.970
   2114_II.CEL        6.545
   2209_A.CEL         7.749
   2209_B.CEL         7.324
   2255_I.CEL         6.793
   2255_II.CEL        6.368
   255_I.CEL          7.484
   255_II.CEL         7.059
   2826_I.CEL         6.119
   2826_II.CEL        5.694
   2853_I.CEL         6.409
   2853_II.CEL        5.984
   2978_I.CEL         7.075
   2978_II.CEL        6.650
   321_I.CEL          7.209
   321_II.CEL         6.784
   3262_I.CEL         6.567
   3262_II.CEL        6.142
   3266_I.CEL         6.730
   3266_II.CEL        6.305
   3271_I.CEL         6.574
   3271_II.CEL        6.149
   3302_I.CEL         6.889
   3302_II.CEL        6.464
   3332_I.CEL         7.095
   3332_II.CEL        6.670Recall that for every row in the design matrix (i.e. every patient sample) only the variables with a 1 in the design matrix are taken into account for calculating the fitted expression value.
This means that as output of the multiplication, we get
a vector crat_fitted whose entries are the sum of relevant
variable coefficients for each sample, respectively.
For example, the fitted value for patient sample 2114_I.CEL
is 6.9703: it is the sum of the respective activated variable coefficients “nI” (7.1917)
and “i_CD2114” (-0.2215).
Let’s visualize the difference between non-inflamed and inflamed tissue again after fitting:
crat_data$fitted_value <- crat_fitted
ggplot(data = crat_data, aes(x = tissue_CD, y = fitted_value, 
                             group = individual, color = individual)) +
      geom_line() +
      ggtitle("Fitted expression changes for the CRAT gene")
Figure 11: Expression changes for the CRAT gene
Note that the difference of the fitted expression values between inflamed and
non-inflamed samples of one patient is the same for all patients and
is determined by the difference between the variable coefficients of
I (6.7667) and nI (7.1917),
which is -0.425
(Figure 11).
This is the case because the same blocking variable is activated by the design matrix for both samples from a single patient, leading to a comparison within patients only. These blocking variables correct the fitted tissue specific expression values towards the expression levels of the individual patients. Therefore the final estimate is like an average of all the within-individual differences.
The “difference” between non-inflamed and inflamed tissue of -0.425 is actually a log2 fold change, as our expression data is on the log2 scale. -0.425 therefore is our log2 fold change for the CRAT gene.
Differential expression analysis of the CRAT gene
In order to test whether the gene is differentially expressed or not, a \(t\)-test with the null hypothesis that there is no difference in the expression between non-inflamed and inflamed tissue is carried out. Our blocking design is conceptually similar to a paired t-test for which the statistic is given by:
\[ t = \frac{\bar{d} }{s/\sqrt{n}} \]
Where, \(\bar{d}\) is the mean difference in expression values between the individuals. The paired t-test computes the variance \(s^2\) from the paired differences. This is lower than the variance of a standard t-test and thus the paired t-test has higher power as long as the expression values for the same individual are correlated (see e.g. the article on Wikipedia).
We thus have improved the power of the ordinary \(t\)-test by reducing the variance via blocking on individuals.
We now conduct the \(t\)-test on the linear model in order to find out whether the difference between non-inflamed and inflamed tissue differs significantly from 0:
crat_noninflamed <- na.exclude(crat_data$org_value[tissue == "nI"])
crat_inflamed <- na.exclude(crat_data$org_value[tissue == "I"])
res_t <- t.test(crat_noninflamed ,crat_inflamed , paired = TRUE)
res_t   
    Paired t-test
   
   data:  crat_noninflamed and crat_inflamed
   t = 6.8, df = 14, p-value = 8e-06
   alternative hypothesis: true mean difference is not equal to 0
   95 percent confidence interval:
    0.2919 0.5581
   sample estimates:
   mean difference 
             0.425We get a low p-value close to 0 and thus can conclude that the CRAT gene is differentially expressed between non-inflamed and inflamed tissue.
Note that the p-value isn’t exactly the same one as below when analyzing the differential expression of all genes. This is due to the variance moderation performed by limma .
We now fit the linear model for all genes and define appropriate contrasts to test hypotheses of interest.
We want to compare the inflamed to the non-inflamed tissue.
Thus, we create a contrast matrix consisting of only one contrast “I-nI”:
limma’s function
makeContrasts creates this matrix from
a symbolic description of the contrast of
interest.
We now fit a linear model to our data and apply the contrasts.fit() function
to it in order to find genes with significant differential expression between non-inflamed and inflamed tissue:
contrast_matrix_CD <- makeContrasts(I-nI, levels = design_palmieri_CD)
palmieri_fit_CD <- eBayes(contrasts.fit(lmFit(palmieri_final[,disease == "CD"],
                                design = design_palmieri_CD),
                                contrast_matrix_CD))
contrast_matrix_UC <- makeContrasts(I-nI, levels = design_palmieri_UC)
palmieri_fit_UC <- eBayes(contrasts.fit(lmFit(palmieri_final[,disease == "UC"],
                                design = design_palmieri_UC),
                                contrast_matrix_UC))We applied the empirical Bayes variance moderation method to the model
via the eBayes() function, which
computes moderated \(t\)-statistics. In microarray analysis,
the number of arrays often is quite small,
and thus variance estimation is difficult.
Using a combination of the per-gene-variance
and a prior variance we can improve the
variance estimate, hence the term “moderation”. “Empirical Bayes”
means that the prior is estimated from the data.
The result of the eBayes() step is that the individual variances are
shrunken towards the prior value.
Finally, we extract the number of differentially expressed genes.
Results can be extracted by use of the topTable function.
We extract the results for both Crohn’s disease and ulcerative colitis, and the
results are sorted by their absolute \(t\)-statistics.
As a diagnostic check, we also plot the p-value histogram: We expect a uniform
distribution for the p-values that correspond to true null hypotheses, while
a peak near zero shows an enrichment for low p-values corresponding to differentially
expressed (DE) genes.
Note that if the p-value distribution for a dataset is very different from the ones in the histograms below, this might lead to quality loss in the subsequent analysis. Reasons for a divergent p-value-distribution might be batch effects or a lack of consideration of other blocking factors in the design matrix. Thus, if the p-value is not as expected, try to include possible blocking factors and batches and rerun the analysis. If this does not help, empirical Bayes / null estimation methods for multiple testing are useful.
A good starting point to learn about these methods is the article on false discovery rate estimation by Korbininan Strimmer (15) and chapter 1–6 of Efron’s book on Large-Scale Inference (16), as well as the blog-post on “How to interpret a p-value histogram” by David Robinson (17).
table_CD <- topTable(palmieri_fit_CD, number = Inf)
head(table_CD)           PROBEID  SYMBOL                            GENENAME   logFC AveExpr
   7928695 7928695  PRXL2A               peroxiredoxin like 2A -0.5990   7.739
   8123695 8123695    ECI2         enoyl-CoA delta isomerase 2 -0.4855   6.876
   8164535 8164535    CRAT       carnitine O-acetyltransferase -0.4250   6.732
   8009746 8009746 SLC16A5   solute carrier family 16 member 5 -0.5182   5.562
   7952249 7952249    USP2      ubiquitin specific peptidase 2 -0.8484   5.606
   8105348 8105348    GPX8 glutathione peroxidase 8 (putative)  0.8312   5.301
                t   P.Value adj.P.Val     B
   7928695 -7.058 1.381e-06   0.02093 5.304
   8123695 -6.317 5.903e-06   0.02093 4.026
   8164535 -6.229 7.037e-06   0.02093 3.870
   8009746 -6.206 7.378e-06   0.02093 3.828
   7952249 -6.204 7.412e-06   0.02093 3.824
   8105348  6.075 9.635e-06   0.02093 3.590hist(table_CD$P.Value, col = brewer.pal(3, name = "Set2")[1],
     main = "inflamed vs non-inflamed - Crohn's disease", xlab = "p-values")
Figure 12: Histogram of the p–values for Crohn’s disease
table_UC <- topTable(palmieri_fit_UC, number = Inf)
head(table_UC)           PROBEID  SYMBOL                                              GENENAME
   8003875 8003875   SPNS2 SPNS lysolipid transporter 2, sphingosine-1-phosphate
   8082012 8082012 SLC15A2                     solute carrier family 15 member 2
   7952290 7952290  TRIM29                        tripartite motif containing 29
   7961693 7961693    LDHB                               lactate dehydrogenase B
   8072015 8072015    GRK3                   G protein-coupled receptor kinase 3
   8096070 8096070    BMP3                          bone morphogenetic protein 3
             logFC AveExpr      t   P.Value adj.P.Val     B
   8003875  0.7412   6.478  7.801 4.549e-07  0.003992 6.440
   8082012  0.8061   5.319  7.744 5.028e-07  0.003992 6.351
   7952290  1.0140   5.855  7.482 8.003e-07  0.003992 5.936
   7961693  0.3968   9.534  7.402 9.249e-07  0.003992 5.806
   8072015  0.4713   5.584  7.308 1.096e-06  0.003992 5.654
   8096070 -1.6961   6.420 -7.289 1.136e-06  0.003992 5.621hist(table_UC$P.Value, col = brewer.pal(3, name = "Set2")[2],
     main = "inflamed vs non-inflamed - Ulcerative colitis", xlab = "p-values")
Figure 13: Histogram of the p–values for ulcerative colitis
In the original paper, a p-value of 0.001 was used as a significance cutoff. Using this we get 944 genes identified as differentially expressed for UC:
nrow(subset(table_UC, P.Value < 0.001))   [1] 944However, it is impossible to determine a precise bound on the number of false positive genes in this list. All that we can say using p-values is that we have at most 21091 (total number of tests) * 0.001 = 21.091 false positive genes in our list. Therefore, by choosing a p-value cutoff of 0.001, as much as 2.23% of our genes identified as differentially expressed might be false positives.
Thus, we can see that the “raw” p-values are very “liberal” when looking at many tests simultaneously. We therefore need error rates adapted to the multiple testing situation. By far the most popular one in molecular biology is the false discovery rate or FDR for short. It is the percentage of false positives among all positives. As we have seen, the FDR of our genes list using a simple p-value cutoff might be quite high.
On the other hand, we can see a clear peak in the p-value histogram (Figure 13), caused by the differentially expressed genes. There we expect the actual FDR of our list to be lower.
The FDR at a given cutoff is given by the “adjusted” p-value in the results table.
tail(subset(table_UC, P.Value < 0.001))           PROBEID   SYMBOL
   7915640 7915640   EIF2B3
   7894577 7894577     <NA>
   7897877 7897877 TNFRSF1B
   8142671 8142671     WASL
   7941946 7941946   NDUFV1
   8140371 8140371 TMEM120A
                                                            GENENAME   logFC
   7915640 eukaryotic translation initiation factor 2B subunit gamma  0.2795
   7894577                                                      <NA>  0.4048
   7897877                        TNF receptor superfamily member 1B  0.3564
   8142671               WASP like actin nucleation promoting factor -0.4281
   7941946            NADH:ubiquinone oxidoreductase core subunit V1 -0.3832
   8140371                                transmembrane protein 120A -0.3657
           AveExpr      t   P.Value adj.P.Val       B
   7915640   6.271  3.957 0.0009861   0.02215 -0.6464
   7894577   5.598  3.956 0.0009888   0.02217 -0.6489
   7897877   6.656  3.956 0.0009894   0.02217 -0.6495
   8142671   7.907 -3.955 0.0009902   0.02217 -0.6502
   7941946   8.867 -3.952 0.0009971   0.02230 -0.6568
   8140371   6.519 -3.952 0.0009985   0.02231 -0.6581The adjusted p-value for a raw p-value of 0.001 in the table is 0.0223, which is an order of magnitude lower than the FDR we can infer from p-values alone.
So although this is not recommended in general, we also use a p-value cutoff at 0.001 in the following in order to be able to compare our workflow results to the paper results.
The paper results can be downloaded as excel files from
http://links.lww.com/IBD/A795 and
should be saved as an .xlsx file named palmieri_DE_res.xlsx in
your working directory.
Note that while the paper uses p-value cutoffs, it also reports the corresponding FDRs (just as we did for the UC data here).
For a p-value cutoff of 0.001, the corresponding FDRs are 0.05 in Crohn’s disease and 0.02 in ulcerative colitis. There are four tables in total, giving the list of up and downregulated genes in CD and UC, respectively. We calculate the overlap between our results and the ones from the paper as the ratio of the genes that were found in both analyses and the genes that were only found in the paper.
We also calculate the total number of diffentially expressed genes that we find in our workflow analysis.
fpath <- system.file("extdata", "palmieri_DE_res.xlsx", package = "maEndToEnd")
palmieri_DE_res <- sapply(1:4, function(i) read.xlsx(cols = 1, fpath, 
                                                     sheet = i, startRow = 4))
names(palmieri_DE_res) <- c("CD_UP", "CD_DOWN", "UC_UP", "UC_DOWN")
palmieri_DE_res <- lapply(palmieri_DE_res, as.character)
paper_DE_genes_CD <- Reduce("c", palmieri_DE_res[1:2])
paper_DE_genes_UC <- Reduce("c", palmieri_DE_res[3:4])
overlap_CD <- length(intersect(subset(table_CD, P.Value < 0.001)$SYMBOL,  
                               paper_DE_genes_CD)) / length(paper_DE_genes_CD)
overlap_UC <- length(intersect(subset(table_UC, P.Value < 0.001)$SYMBOL,
                               paper_DE_genes_UC)) / length(paper_DE_genes_UC)
overlap_CD   [1] 0.6376overlap_UC    [1] 0.6654total_genenumber_CD <- length(subset(table_CD, P.Value < 0.001)$SYMBOL)
total_genenumber_UC <- length(subset(table_UC, P.Value < 0.001)$SYMBOL)
total_genenumber_CD   [1] 577total_genenumber_UC   [1] 944We find 577 (CD) and 944 (UC) differentially expressed genes (“DE-genes”).
In the paper, 298 (CD) and 520 (UC) DE-genes were found for the two diseases at the same cutoff. This higher number of DE-genes identified is probably due to the increased power of the blocking according to the individuals and the moderated variance estimation that limma performs.
We see that we get a moderate overlap of 0.6376 for CD and 0.6654 for UC, showing that both analyses lead to somewhat comparable results.
For a visualization of the differentially expressed genes, we create a volcano plot, which is commonly used to summarize the results of a differential expression analysis in a single figure.
For a better overview, we only show gene symbols of genes with
a fold change greater than 1, which we define in the volcano_names object. The
highlight option in the volcanoplot function is set to 100 and thus only
labels the 100 genes with the lowest p-values.
volcano_names <- ifelse(abs(palmieri_fit_CD$coefficients)>=1, 
                        palmieri_fit_CD$genes$SYMBOL, NA)
             
             
volcanoplot(palmieri_fit_CD, coef = 1L, style = "p-value", highlight = 100, 
            names = volcano_names,
            xlab = "Log2 Fold Change", ylab = NULL, pch=16, cex=0.35)
Figure 14: Volcano plot of the DE-genes
We can now do a little research on the biological function of genes that show a high foldchange, for example the gene with the symbol S100A8 on the right side of the plot (Figure 14). If we search for this gene symbol on genecards.org, we find that it encodes for a protein that builds a pro-inflammatory complex in association with another protein.
As discussed above, it is recommended to use an FDR cutoff in differential expression analysis rather than a p-value cutoff, since this way you control an explicitly defined error rate and the results are easier to interpret and to compare. For the following enrichment analysis, we create tables with differentially expressed genes for CD and UC, respectively, and choose an FDR cutoff of 10%. Here, we focus on the CD subset of the data.
DE_genes_CD <- subset(table_CD, adj.P.Val < 0.1)$PROBEIDWe can now try to characterize the identified differentially expressed genes more in detail by performing a GO enrichment analysis. Essentially the gene ontology (http://www.geneontology.org/) is a hierarchically organized collection of functional gene sets (18–20).
The function genefinder from the genefilter
package (21)
will be used to find a background set of genes that are similar in expression
to the differentially expressed genes. We then check whether
the background has roughly the same distribution
of average expression strength as the foreground.
We do this in order not to select a biased background since the gene set testing is performed by a simple Fisher test on a 2x2 table. Note that this approach is very similar to commonly used web tools like GOrilla (22).
For every differentially expressed gene, we try to find genes with similar
expression with genefinder. The genefinder function returns a list with two
elements for each gene: one with the indices of the
background genes found and one with the distances to the DE-genes:
back_genes_idx <- genefilter::genefinder(palmieri_final, 
                                        as.character(DE_genes_CD), 
                                        method = "manhattan", scale = "none")We have to extract the PROBEIDs, which correspond to the indices.
We do that by using the sapply
function, which gives us a single matrix with the DE-genes as column names
and the PROBEIDs of the corresponding background genes in the cells below:
back_genes_idx <- sapply(back_genes_idx, function(x)x$indices)We then create a vector back_genes containing all background gene PROBEIDs:
In order to eliminate foreground genes, i.e. DE-genes, from the back_genes
set, we use the setdiff function. It returns all elements from the first
argument (back_genes) that are not part of the second argument
(DE_genes_CD).
With the intersect function, we verify that we were successful:
it should return 0, as there shouldn’t be any intersect anymore between
back_genesand DE_genes_CD:
back_genes <- featureNames(palmieri_final)[back_genes_idx]
back_genes <- setdiff(back_genes, DE_genes_CD)
    
intersect(back_genes, DE_genes_CD)   character(0)length(back_genes)   [1] 9772We create a multidensity plot with mean expression on the x-axis and curves for all genes, foreground genes and background genes, respectively (Figure 15). We want to see whether the background genes show a plot similar to the foreground genes so that the background is not biased for the gene enrichment analysis:
multidensity(list(
        all = table_CD[,"AveExpr"] ,
        fore = table_CD[DE_genes_CD , "AveExpr"],
        back = table_CD[rownames(table_CD) %in% back_genes, "AveExpr"]),
        col = c("#e46981", "#ae7ee2", "#a7ad4a"),
     xlab = "mean expression",
   main = "DE genes for CD-background-matching")
Figure 15: Selecting a background set of genes for the gene ontology analysis
When comparing the “background gene” curve to the “foreground gene” curve, we see a similar curve shape, indicating a sensible background matching (Figure 15). Note that the right-shift of the “foreground-gene” curve in comparison to the “background-gene” curve indicates that DE-genes are generally very highly expressed, so that it wasn’t possible to find background-genes with exactly equal overall expression distribution.
The “all gene” curve has the leftmost curve maximum; this can be explained by a high number of lowly expressed genes in all samples and shows that a background matching is sensible in order to avoid biases.
For the actual testing of which GO gene sets are enriched in inflamed tissue, we use the topGO package which implements a nice interface to Fisher testing and also has additional algorithms taking the GO structure into account, by e.g. only reporting the most specific gene set in the hierarchy (23).
The GO has three top ontologies: Cellular component (CC), biological processes (BP), and molecular function (MF). For illustrative purposes we limit ourselves to the BP category here.
topGO requires a topGOdata object containing the necessary information for the
analysis. We follow the steps described in the topGO vignettes: First,
we will create a named vector all_genes with all genes to be analyzed, i.e. 
DE-genes and background genes:
gene_IDs <- rownames(table_CD)
in_universe <- gene_IDs %in% c(DE_genes_CD, back_genes)
in_selection <- gene_IDs %in% DE_genes_CD 
all_genes <- in_selection[in_universe]
all_genes <- factor(as.integer(in_selection[in_universe]))
names(all_genes) <- gene_IDs[in_universe] The following steps were carried through:
we created an in_universe vector by using the %in% matching function.
We want to know which elements from gene_IDs are also contained in
DE_genes_CD and back_genes, as the latter two are our gene universe
we use for enrichment analysis.
We got a vector in_universe with the length of gene_IDs that has the
entry TRUE when the corresponding gene in gene_IDs could be also found in
DE_genes_CD or back_genes, and FALSE otherwise.
We did the same for our DE-genes and call this vector in_selection.
We created the all_genes vector:
a) First, we selected all the elements from in_selection
that are TRUE in in_universe by applying
all_genes <- in_selection[in_universe].
b) Then, we converted the elements in all_genes from TRUE and FALSE
to 0 and 1 by converting the vector to an integer vector.
This way, each element in the vector is a 0 if the corresponding gene is a
background gene and a 1 if the corresponding gene is a DE-gene. Also, we
converted the vector to a factor.
c) We named the vector elements with the corresponding gene_IDs.
We now initialize the topGO data set, using the GO annotations
contained in the annotation data base for the chip we are using. The nodeSize
parameter specifies a minimum size of a GO category we want to use: i.e. here,
categories with less than 10 genes are not included in the testing.
top_GO_data <- new("topGOdata", ontology = "BP", allGenes = all_genes,
 nodeSize = 10, annot = annFUN.db, affyLib = "hugene10sttranscriptcluster.db")Now the tests can be run. topGO offers a wide range of options, for details see the paper (23) or the package vignette.
We run two common tests: an ordinary Fisher test for every GO category, and the “elim” algorithm, which tries to incorporate the hierarchical structure of the GO and tries to “decorrelate” it in order to report the most specific significant term in the hierarchy.
The algorithm starts processing the nodes / GO categories from the highest (bottommost) level and then iteratively moves to nodes from a lower level. If a node is scored as significant, all of its genes are marked as removed in all ancestor nodes. This way, the “elim” algorithm aims at finding the most specific node for every gene.
The test uses a 0.01 p-value cutoff by default.
result_top_GO_elim <- 
  runTest(top_GO_data, algorithm = "elim", statistic = "Fisher")
result_top_GO_classic <- 
  runTest(top_GO_data, algorithm = "classic", statistic = "Fisher")We can now inspect the results. We look at the top 100 GO categories according
to the “Fisher elim” algorithm. The function GenTable produces
a table of significant GO categories, the function printGenes
gives genes annotated to them; the significant ones are denoted with a “2”
in the “raw p-value” column, the non-significant ones with a “1”.
We therefore select raw p-value == 2.
Note that we do not get the actual p-values here because our all_genes vector
doesn’t contain this information; it only tells us
whether a gene is differentially expressed or not.
res_top_GO <- GenTable(top_GO_data, Fisher.elim = result_top_GO_elim,
        Fisher.classic = result_top_GO_classic,
        orderBy = "Fisher.elim" , topNodes = 100)
genes_top_GO <- printGenes(top_GO_data, whichTerms = res_top_GO$GO.ID,
    chip = "hugene10sttranscriptcluster.db", geneCutOff = 1000)
res_top_GO$sig_genes <- sapply(genes_top_GO, function(x){
                str_c(paste0(x[x$'raw p-value' == 2, "Symbol.id"],";"), 
                      collapse = "")
    })
head(res_top_GO[,1:8], 20)           GO.ID                                        Term Annotated Significant
   1  GO:0030593                       neutrophil chemotaxis        69          42
   2  GO:0071222     cellular response to lipopolysaccharide       148          62
   3  GO:0006954                       inflammatory response       533         209
   4  GO:0051897 positive regulation of phosphatidylinosi...       110          49
   5  GO:0030335       positive regulation of cell migration       384         152
   6  GO:0001774                  microglial cell activation        41          24
   7  GO:0002250                    adaptive immune response       342         128
   8  GO:0042129          regulation of T cell proliferation       130          60
   9  GO:0001676     long-chain fatty acid metabolic process        54          29
   10 GO:0006067                   ethanol metabolic process        11          10
   11 GO:0006805                xenobiotic metabolic process        53          27
   12 GO:0042119                       neutrophil activation        34          20
   13 GO:0016525         negative regulation of angiogenesis        68          32
   14 GO:0070098        chemokine-mediated signaling pathway        54          27
   15 GO:0042102 positive regulation of T cell proliferat...        82          36
   16 GO:0009617                       response to bacterium       390         157
   17 GO:0045087                      innate immune response       589         189
   18 GO:2001238 positive regulation of extrinsic apoptot...        41          22
   19 GO:0071346     cellular response to type II interferon        83          36
   20 GO:0019221         cytokine-mediated signaling pathway       345         121
      Expected Rank in Fisher.classic Fisher.elim Fisher.classic
   1     14.96                     72     1.9e-12        1.9e-12
   2     32.08                    173     2.1e-08        2.1e-08
   3    115.53                      7     2.2e-08        2.6e-21
   4     23.84                    188     6.2e-08        6.2e-08
   5     83.23                     27     1.3e-07        3.9e-16
   6      8.89                    227     3.2e-07        3.2e-07
   7     74.13                     86     4.1e-07        1.0e-11
   8     28.18                    112     1.8e-06        3.7e-10
   9     11.70                    223     1.9e-06        2.4e-07
   10     2.38                    274     2.0e-06        2.0e-06
   11    11.49                    282     2.5e-06        2.5e-06
   12     7.37                    285     2.8e-06        2.8e-06
   13    14.74                    288     2.8e-06        2.8e-06
   14    11.70                    296     3.9e-06        3.9e-06
   15    17.77                    306     5.2e-06        5.2e-06
   16    84.53                     16     5.3e-06        2.0e-17
   17   127.67                    118     5.6e-06        9.1e-10
   18     8.89                    326     7.1e-06        7.1e-06
   19    17.99                    328     7.3e-06        7.3e-06
   20    74.78                    141     9.4e-06        4.2e-09A graph of the results can also be produced. Here we visualize the three most significant nodes according to the Fisher elim algorithm in the context of the GO hierarchy.
Note that the following code is hashed out so that it won’t be executed on the Bioconductor build, since it caused the build to fail. However, it can be run on a local machine without concerns.
# showSigOfNodes(top_GO_data, score(result_top_GO_elim), firstSigNodes = 3,
#                useInfo = 'def')We can see that indeed GO categories related to inflammation, signalling and immune response come up as significant (Figure ??) Gene set enrichment analysis has been a field of very extensive research in bioinformatics. For additional approaches see the topGO vignette and the references therein and also in the GeneSetEnrichment view.
The package ReactomePA offers the possibility to test enrichment
of specific pathways using the free, open-source, curated and peer reviewed
Reactome pathway database (24,25).
The package requires entrez identifiers, so we convert our PROBEIDs
(transcript cluster identifiers) to entrez identifiers
using the function mapIDs from the package AnnotationDbi.
This will create a named vector that maps the PROBEIDs to the entrez ones,
with the PROBEIDs as names and the entrez ids as vector elements.
entrez_ids <- mapIds(hugene10sttranscriptcluster.db, 
      keys = rownames(table_CD), 
      keytype = "PROBEID",
      column = "ENTREZID")We can now run the enrichment analysis that performs a statistical test based on the hypergeoemtric distribution that is the same as a one sided Fisher-test, which topGO calls “Fisher-classic”. Details can be found in the vignette of the DOSE package (26).
reactome_enrich <- enrichPathway(gene = entrez_ids[DE_genes_CD], 
                                universe = entrez_ids[c(DE_genes_CD, 
                                                        back_genes)],
                                organism = "human",
                                pvalueCutoff = 0.05,
                                qvalueCutoff = 0.9, 
                                readable = TRUE)
reactome_enrich@result$Description <- paste0(str_sub(
                                    reactome_enrich@result$Description, 1, 20),
                                    "...")
head(as.data.frame(reactome_enrich))[1:6]                            ID             Description GeneRatio  BgRatio
   R-HSA-8978868 R-HSA-8978868 Fatty acid metabolis...   54/1428 105/6265
   R-HSA-6785807 R-HSA-6785807 Interleukin-4 and In...   41/1428  77/6265
   R-HSA-556833   R-HSA-556833 Metabolism of lipids...  159/1428 467/6265
   R-HSA-380108   R-HSA-380108 Chemokine receptors ...   22/1428  31/6265
   R-HSA-6783783 R-HSA-6783783 Interleukin-10 signa...   22/1428  31/6265
   R-HSA-1474244 R-HSA-1474244 Extracellular matrix...   74/1428 187/6265
                    pvalue  p.adjust
   R-HSA-8978868 1.070e-10 1.313e-07
   R-HSA-6785807 5.118e-09 2.122e-06
   R-HSA-556833  5.188e-09 2.122e-06
   R-HSA-380108  1.497e-08 3.673e-06
   R-HSA-6783783 1.497e-08 3.673e-06
   R-HSA-1474244 1.320e-07 2.699e-05Note that we trimmed pathway names to 20 characters.
The top pathways can be displayed as a bar chart that displays all categories with a p-value below the specified cutoff (Figure 16).
barplot(reactome_enrich)
Figure 16: Enriched Reactome pathways and their p–values as a bar chart
The “enrichment map” from the package enrichplot displays the results of the enrichment analysis as a graph, where the color represents the p-value of the pathway and the edge-thickness (that is the line connecting two pathways) is proportional to the number of overlapping genes between two pathways.
reactome_enrich <- pairwise_termsim(reactome_enrich)
emapplot(reactome_enrich, showCategory = 10)
Figure 17: Enriched Reactome pathways enrichment results as a graph
Again, the graph in Figure 17 shows pathways related to signalling and immune response.
The package clusterProfiler (27) can also perform these analyses using downloaded KEGG data. Furthermore, the package EnrichmentBrowser (28) additionally offers network-based enrichment analysis of individual pathways. This allows the mapping of the expression data at hand to known regulatory interactions.
As the last part of this document, we call the function sessionInfo, which reports the version numbers of R and all the packages used in this session. It is good practice to always keep such a record of this as it will help to track down what has happened in case an R script ceases to work or gives different results because the functions have been changed in a newer version of one of your packages. By including it at the bottom of a script, your reports will become more reproducible.
The session information should also always be included in any emails to the Bioconductor support site along with all code used in the analysis.
gc()              used   (Mb) gc trigger (Mb)  max used (Mb)
   Ncells 1.38e+07  737.2   25258832 1349  25258832 1349
   Vcells 4.30e+08 3281.0  868370122 6625 868345513 6625length(getLoadedDLLs())   [1] 117sessionInfo()   R version 4.4.0 beta (2024-04-15 r86425)
   Platform: x86_64-pc-linux-gnu
   Running under: Ubuntu 22.04.4 LTS
   
   Matrix products: default
   BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
   LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
   
   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       
   
   time zone: America/New_York
   tzcode source: system (glibc)
   
   attached base packages:
   [1] grid      stats4    stats     graphics  grDevices utils     datasets 
   [8] methods   base     
   
   other attached packages:
    [1] maEndToEnd_2.24.0                    enrichplot_1.24.0                   
    [3] Rgraphviz_2.48.0                     openxlsx_4.2.5.2                    
    [5] genefilter_1.86.0                    matrixStats_1.3.0                   
    [7] stringr_1.5.1                        tidyr_1.3.1                         
    [9] dplyr_1.1.4                          RColorBrewer_1.1-3                  
   [11] pheatmap_1.0.12                      geneplotter_1.82.0                  
   [13] annotate_1.82.0                      XML_3.99-0.16.1                     
   [15] lattice_0.22-6                       ggplot2_3.5.1                       
   [17] gplots_3.1.3.1                       clusterProfiler_4.12.0              
   [19] ReactomePA_1.48.0                    topGO_2.56.0                        
   [21] SparseM_1.81                         GO.db_3.19.1                        
   [23] graph_1.82.0                         limma_3.60.0                        
   [25] arrayQualityMetrics_3.60.0           hugene10sttranscriptcluster.db_8.8.0
   [27] org.Hs.eg.db_3.19.1                  AnnotationDbi_1.66.0                
   [29] pd.hugene.1.0.st.v1_3.14.1           DBI_1.2.2                           
   [31] oligo_1.68.0                         RSQLite_2.3.6                       
   [33] Biostrings_2.72.0                    GenomeInfoDb_1.40.0                 
   [35] XVector_0.44.0                       IRanges_2.38.0                      
   [37] S4Vectors_0.42.0                     ArrayExpress_1.64.0                 
   [39] oligoClasses_1.66.0                  Biobase_2.64.0                      
   [41] BiocGenerics_0.50.0                  knitr_1.46                          
   [43] BiocStyle_2.32.0                    
   
   loaded via a namespace (and not attached):
     [1] fs_1.6.4                    bitops_1.0-7               
     [3] HDO.db_0.99.1               httr_1.4.7                 
     [5] tools_4.4.0                 gcrma_2.76.0               
     [7] backports_1.4.1             utf8_1.2.4                 
     [9] R6_2.5.1                    lazyeval_0.2.2             
    [11] withr_3.0.0                 graphite_1.50.0            
    [13] gridExtra_2.3               base64_2.0.1               
    [15] preprocessCore_1.66.0       cli_3.6.2                  
    [17] scatterpie_0.2.2            labeling_0.4.3             
    [19] sass_0.4.9                  askpass_1.2.0              
    [21] systemfonts_1.0.6           yulab.utils_0.1.4          
    [23] gson_0.1.0                  foreign_0.8-86             
    [25] illuminaio_0.46.0           DOSE_3.30.0                
    [27] svglite_2.1.3               affyPLM_1.80.0             
    [29] BeadDataPackR_1.56.0        rstudioapi_0.16.0          
    [31] generics_0.1.3              gridGraphics_0.5-1         
    [33] hwriter_1.3.2.1             gtools_3.9.5               
    [35] zip_2.3.1                   Matrix_1.7-0               
    [37] interp_1.1-6                fansi_1.0.6                
    [39] abind_1.4-5                 lifecycle_1.0.4            
    [41] yaml_2.3.8                  SummarizedExperiment_1.34.0
    [43] beadarray_2.54.0            qvalue_2.36.0              
    [45] SparseArray_1.4.0           blob_1.2.4                 
    [47] affxparser_1.76.0           crayon_1.5.2               
    [49] cowplot_1.1.3               KEGGREST_1.44.0            
    [51] magick_2.8.3                pillar_1.9.0               
    [53] fgsea_1.30.0                GenomicRanges_1.56.0       
    [55] codetools_0.2-20            fastmatch_1.1-4            
    [57] glue_1.7.0                  ggfun_0.1.4                
    [59] data.table_1.15.4           vctrs_0.6.5                
    [61] png_0.1-8                   treeio_1.28.0              
    [63] gtable_0.3.5                cachem_1.0.8               
    [65] xfun_0.43                   S4Arrays_1.4.0             
    [67] tidygraph_1.3.1             survival_3.6-4             
    [69] iterators_1.0.14            tinytex_0.50               
    [71] statmod_1.5.0               nlme_3.1-164               
    [73] ggtree_3.12.0               bit64_4.0.5                
    [75] bslib_0.7.0                 affyio_1.74.0              
    [77] KernSmooth_2.23-22          rpart_4.1.23               
    [79] colorspace_2.1-0            Hmisc_5.1-2                
    [81] nnet_7.3-19                 tidyselect_1.2.1           
    [83] curl_5.2.1                  bit_4.0.5                  
    [85] compiler_4.4.0              htmlTable_2.4.2            
    [87] DelayedArray_0.30.0         bookdown_0.39              
    [89] shadowtext_0.1.3            checkmate_2.3.1            
    [91] scales_1.3.0                caTools_1.18.2             
    [93] hexbin_1.28.3               affy_1.82.0                
    [95] rappdirs_0.3.3              digest_0.6.35              
    [97] rmarkdown_2.26              htmltools_0.5.8.1          
    [99] pkgconfig_2.0.3             jpeg_0.1-10                
   [101] base64enc_0.1-3             MatrixGenerics_1.16.0      
   [103] highr_0.10                  fastmap_1.1.1              
   [105] rlang_1.1.3                 htmlwidgets_1.6.4          
   [107] UCSC.utils_1.0.0            farver_2.1.1               
   [109] jquerylib_0.1.4             jsonlite_1.8.8             
   [111] BiocParallel_1.38.0         GOSemSim_2.30.0            
   [113] magrittr_2.0.3              Formula_1.2-5              
   [115] GenomeInfoDbData_1.2.12     ggplotify_0.1.2            
   [117] patchwork_1.2.0             munsell_0.5.1              
   [119] Rcpp_1.0.12                 ggnewscale_0.4.10          
   [121] ape_5.8                     viridis_0.6.5              
   [123] vsn_3.72.0                  stringi_1.8.3              
   [125] ggraph_2.2.1                zlibbioc_1.50.0            
   [127] MASS_7.3-60.2               plyr_1.8.9                 
   [129] parallel_4.4.0              ggrepel_0.9.5              
   [131] deldir_2.0-4                graphlayouts_1.1.1         
   [133] splines_4.4.0               igraph_2.0.3               
   [135] reshape2_1.4.4              evaluate_0.23              
   [137] latticeExtra_0.6-30         BiocManager_1.30.22        
   [139] foreach_1.5.2               tweenr_2.0.3               
   [141] openssl_2.1.2               purrr_1.0.2                
   [143] polyclip_1.10-6             ggforce_0.4.2              
   [145] xtable_1.8-4                ff_4.0.12                  
   [147] reactome.db_1.88.0          tidytree_0.4.6             
   [149] viridisLite_0.4.2           tibble_3.2.1               
   [151] aplot_0.2.2                 memoise_2.0.1              
   [153] setRNG_2024.2-1             cluster_2.1.6              
   [155] gridSVG_1.7-5The authors would like to thank Vladislava Milchevskaya, Julian Gehring and Mike Smith for helpful comments on and small contributions to the workflow. They also would like to thank Frederik Ziebell, Holly Giles and Jennifer Huellein for proof-reading.
This workflow draws a lot of inspiration from the Bioconductor books (29,30) as well as Love et. al.’s workflow for gene level analysis of RNA-Seq data (31). James W. MacDonald provided valuable information on the evolution of Affymetrix arrays in some of his posts of on the Biocondctor support site. BK would also like to thank him for some friendly personal correspondence about the annotation resources available for microarrays in Bioconductor.
1. Palmieri O, Creanza TM, Bossa F, Palumbo O, Maglietta R, Ancona N, et al. Genome-wide pathway analysis using gene expression data of colonic mucosa in patients with inflammatory bowel disease. Inflammatory Bowel Diseases [Internet]. 2015 Apr;1. Available from: http://dx.doi.org/10.1097/MIB.0000000000000370
2. Kolesnikov N, Hastings E, Keays M, Melnichuk O, Tang YA, Williams E, et al. ArrayExpress update–simplifying data submissions. Nucleic Acids Research [Internet]. 2014 Oct;43(D1):D1113–6. Available from: http://dx.doi.org/10.1093/nar/gku1057
3. Huber, W., Carey, J. V, Gentleman, R., et al. Orchestrating high-throughput genomic analysis with Bioconductor. Nature Methods [Internet]. 2015;12(2):115–21. Available from: http://www.nature.com/nmeth/journal/v12/n2/full/nmeth.3252.html
4. Carvalho BS, Irizarry RA. A framework for oligonucleotide microarray preprocessing. Bioinformatics. 2010;26(19):2363–7.
5. Kauffmann A, Gentleman R, Huber W. ArrayQualityMetrics–a bioconductor package for quality assessment of microarray data. Bioinformatics. 2009;25(3):415–6.
6. Irizarry RA. Summaries of affymetrix GeneChip probe level data. Nucleic Acids Research [Internet]. 2003 Feb;31(4):15e–15. Available from: http://dx.doi.org/10.1093/nar/gng015
7. Bolstad BM, Irizarry RA, Astrand M, Speed TP. A comparison of normalization methods for high density oligonucleotide array data based on variance and bias. Bioinformatics [Internet]. 2003 Jan;19(2):185–93. Available from: http://dx.doi.org/10.1093/bioinformatics/19.2.185
8. Irizarry RA. Exploration, normalization, and summaries of high density oligonucleotide array probe level data. Biostatistics [Internet]. 2003 Apr;4(2):249–64. Available from: http://dx.doi.org/10.1093/biostatistics/4.2.249
9. Gandolfo LC, Speed TP. RLE plots: Visualizing unwanted variation in high dimensional data. PLoS ONE [Internet]. 2018; Available from: http://www.ncbi.nlm.nih.gov/pmc/articles/PMC5798764/
10. Huber W, von Heydebreck A, Sueltmann H, Poustka A, Vingron M. Variance stabilization applied to microarray data calibration and to the quantification of differential expression. Bioinformatics. 2002;18 Suppl. 1:S96–S104.
11. Ritchie ME, Phipson B, Wu D, Hu Y, Law CW, Shi W, et al. limma powers differential expression analyses for RNA-sequencing and microarray studies. Nucleic Acids Research. 2015;43(7):e47.
12. Smyth GK. Linear models and empirical bayes methods for assessing differential expression in microarray experiments. Statistical Applications in Genetics and Molecular Biology [Internet]. 2004 Jan;3(1):1–25. Available from: http://dx.doi.org/10.2202/1544-6115.1027
13. Arloth J, Bader DM, Röh S, Altmann A. Re-annotator: Annotation pipeline for microarray probe sequences. Xue Y, editor. PLOS ONE [Internet]. 2015 Oct;10(10):e0139516. Available from: http://dx.doi.org/10.1371/journal.pone.0139516
14. Irizarry RA, Love MI. PH525x series - biomedical data science [Internet]. 2015. Available from: http://genomicsclass.github.io/book/
15. Strimmer K. A unified approach to false discovery rate estimation. BMC Bioinformatics. 2008;
16. Efron B, Cox D, Meng XL, Holmes S, Gillooly D, editors. Large-scale inference: Empirical bayes methods for estimation, testing and prediction [Internet]. 2010. Available from: https://statweb.stanford.edu/~ckirby/brad/LSI/monograph_CUP.pdf
17. Robinson D. How to interpret a p-value histogram. http://varianceexplained.org/statistics/interpreting-pvalue-histogram; 2014.
18. Ashburner M, Ball CA, Blake JA, Botstein D, Butler H, Cherry JM, et al. Gene ontology: Tool for the unification of biology. Nature Genetics [Internet]. 2000 May;25(1):25–9. Available from: http://dx.doi.org/10.1038/75556
19. Gene Ontology Consortium. Gene ontology consortium: Going forward. Nucleic acids research [Internet]. 2015;43(Database issue):D1049—56. Available from: http://europepmc.org/articles/PMC4383973
20. du Plessis L, Skunca N, Dessimoz C. The what, where, how and why of gene ontology–a primer for bioinformaticians. Briefings in Bioinformatics [Internet]. 2011 Feb;12(6):723–35. Available from: http://dx.doi.org/10.1093/bib/bbr002
21. Bourgon R, Gentleman R, Huber W. Independent filtering increases detection power for high-throughput experiments. Proceedings of the National Academy of Sciences [Internet]. 2010 May;107(21):9546–51. Available from: http://dx.doi.org/10.1073/pnas.0914005107
22. Eden E, Navon R, Steinfeld I, Lipson D, Yakhini Z. GOrilla: A tool for discovery and visualization of enriched GO terms in ranked gene lists. BMC Bioinformatics [Internet]. 2009;10(1):48. Available from: http://dx.doi.org/10.1186/1471-2105-10-48
23. Alexa A, Rahnenfuhrer J, Lengauer T. Improved scoring of functional groups from gene expression data by decorrelating GO graph structure. Bioinformatics [Internet]. 2006 Apr;22(13):1600–7. Available from: http://dx.doi.org/10.1093/bioinformatics/btl140
24. Croft D, Mundo AF, Haw R, Milacic M, Weiser J, Wu G, et al. The reactome pathway knowledgebase. Nucleic Acids Research [Internet]. 2013 Nov;42(D1):D472–7. Available from: http://dx.doi.org/10.1093/nar/gkt1102
25. Fabregat A, Jupe S, Matthews L, Sidiropoulos K, Gillespie M, Garapati P, et al. The reactome pathway knowledgebase. Nucleic Acids Research [Internet]. 2017 Nov;46(D1):D649–55. Available from: https://doi.org/10.1093/nar/gkx1132
26. Yu G, Wang LG, Yan GR, He QY. DOSE: An r/bioconductor package for disease ontology semantic and enrichment analysis. Bioinformatics [Internet]. 2014 Oct;31(4):608–9. Available from: http://dx.doi.org/10.1093/bioinformatics/btu684
27. Yu G, Wang LG, Han Y, He QY. clusterProfiler: An r package for comparing biological themes among gene clusters. OMICS: A Journal of Integrative Biology [Internet]. 2012 May;16(5):284–7. Available from: http://dx.doi.org/10.1089/omi.2011.0118
28. Geistlinger L, Csaba G, Zimmer R. Bioconductor’s EnrichmentBrowser: Seamless navigation through combined results of set- & network-based enrichment analysis. BMC Bioinformatics [Internet]. 2016 Jan;17(1). Available from: http://dx.doi.org/10.1186/s12859-016-0884-1
29. Gentleman R, Carey VJ, Huber W, Irizarry R, Dudoit S, editors. Bioinformatics and computational biology solutions using R and Bioconductor [Internet]. Springer; 2005. Available from: http://www.bioconductor.org/help/publications/books/bioinformatics-and-computational-biology-solutions
30. Hahne F, Huber W, Gentleman R, Falcon S. Bioconductor case studies [Internet]. Springer; 2008. (Use r). Available from: http://www.bioconductor.org/help/publications/books/bioconductor-case-studies
31. Love MI, Anders S, Kim V, Huber W. RNA-seq workflow: Gene-level exploratory analysis and differential expression. F1000Research [Internet]. 2016 Nov;4:1070. Available from: https://doi.org/10.12688/f1000research.7035.2