Contents

1 Presentation: RNA-seq work flow

Resources: Anders, CSAMA, 2015; Love, CSAMA, 2015; Huber, CSAMA, 2015; Pimentel, YouTube, 2015.

1.1 Experimental design

Keep it simple

Replicate

Avoid confounding experimental factors with other factors

Record co-variates

Be aware of batch effects

HapMap samples from one facility, ordered by date of processing.

1.2 Wet-lab

Confounding factors

Artifacts of your particular protocols

1.3 Sequencing

Axes of variation

Application-specific, e.g.,

1.4 Alignment

Alignment strategies

Splice-aware aligners (and Bioconductor wrappers)

1.5 Reduction to ‘count tables’

1.5.2 (kallisto / sailfish)

  • ‘Next generation’ differential expression tools; transcriptome alignment
  • E.g., kallisto takes a radically different approach: from FASTQ to count table without BAM files.
  • Very fast, almost as accurate.
  • Hints on how it works; arXiv
  • Integration with gene-level analyses – Soneson et al.

1.6 Analysis

Unique statistical aspects

Summarization

Normalization

Appropriate error model

Pre-filtering

Borrowing information

1.7 Statistical Issues In-depth: Normalization and Dispersion

1.7.1 Normalization

DESeq2 estimateSizeFactors(), Anders and Huber, 2010

  • For each gene: geometric mean of all samples.
  • For each sample: median ratio of the sample gene over the geometric mean of all samples
  • Functions other than the median can be used; control genes can be used instead

1.7.2 Dispersion

DESeq2 estimateDispersions()

  • Estimate per-gene dispersion
  • Fit a smoothed relationship between dispersion and abundance

1.8 Comprehension

Placing differentially expressed regions in context

Copy number / expression QC Correlation between genomic copy number and mRNA expression identified 38 mis-labeled samples in the TCGA ovarian cancer Affymetrix microarray dataset.

2 Lab: Gene-level RNA-seq differential expression

2.1 Background

This lab is derived from: RNA-Seq workflow: gene-level exploratory analysis and differential expression, by Michael Love, Simon Anders, Wolfgang Huber; modified by Martin Morgan, October 2015.

This lab will walk you through an end-to-end RNA-Seq differential expression workflow, using DESeq2 along with other Bioconductor packages. The complete work flow starts from the FASTQ files, but we will start after reads have been aligned to a reference genome and reads overlapping known genes have been counted. We will perform exploratory data analysis (EDA), differential gene expression analysis with DESeq2, and visually explore the results.

A number of other Bioconductor packages are important in statistical inference of differential expression at the gene level, including Rsubread, edgeR, limma, BaySeq, and others.

2.2 Experimental data

The data used in this workflow is an RNA-Seq experiment of airway smooth muscle cells treated with dexamethasone, a synthetic glucocorticoid steroid with anti-inflammatory effects. Glucocorticoids are used, for example, in asthma patients to prevent or reduce inflammation of the airways. In the experiment, four primary human airway smooth muscle cell lines were treated with 1 micromolar dexamethasone for 18 hours. For each of the four cell lines, we have a treated and an untreated sample. The reference for the experiment is:

Himes BE, Jiang X, Wagner P, Hu R, Wang Q, Klanderman B, Whitaker RM, Duan Q, Lasky-Su J, Nikolos C, Jester W, Johnson M, Panettieri R Jr, Tantisira KG, Weiss ST, Lu Q. “RNA-Seq Transcriptome Profiling Identifies CRISPLD2 as a Glucocorticoid Responsive Gene that Modulates Cytokine Function in Airway Smooth Muscle Cells.” PLoS One. 2014 Jun 13;9(6):e99625. PMID: 24926665. GEO: GSE52778.

2.3 Preparing count matrices

As input, DESeq2 package expects count data as obtained, e.g., from RNA-Seq or another high-throughput sequencing experiment, in the form of a matrix of integer values. The value in the i-th row and the j-th column of the matrix tells how many reads have been mapped to gene i in sample j. Analogously, for other types of assays, the rows of the matrix might correspond e.g., to binding regions (with ChIP-Seq) or peptide sequences (with quantitative mass spectrometry).

The count values must be raw counts of sequencing reads. This is important for DESeq2’s statistical model to hold, as only the actual counts allow assessing the measurement precision correctly. Hence, please do not supply other quantities, such as (rounded) normalized counts, or counts of covered base pairs – this will only lead to nonsensical results.

We will discuss how to summarize data from BAM files to a count table later in ther course. Here we’ll ‘jump right in’ and start with a prepared SummarizedExperiment.

2.4 Starting from SummarizedExperiment

We now use R’s data() command to load a prepared SummarizedExperiment that was generated from the publicly available sequencing data files associated with the Himes et al. paper, described above. The steps we used to produce this object were equivalent to those you worked through in the previous sections, except that we used all the reads and all the genes. For more details on the exact steps used to create this object type vignette("airway") into your R session.

library(airway)
data("airway")
se <- airway

The information in a SummarizedExperiment object can be accessed with accessor functions. For example, to see the actual data, i.e., here, the read counts, we use the assay() function. (The head() function restricts the output to the first few lines.)

head(assay(se))
##                 SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520
## ENSG00000000003        679        448        873        408       1138       1047        770
## ENSG00000000005          0          0          0          0          0          0          0
## ENSG00000000419        467        515        621        365        587        799        417
## ENSG00000000457        260        211        263        164        245        331        233
## ENSG00000000460         60         55         40         35         78         63         76
## ENSG00000000938          0          0          2          0          1          0          0
##                 SRR1039521
## ENSG00000000003        572
## ENSG00000000005          0
## ENSG00000000419        508
## ENSG00000000457        229
## ENSG00000000460         60
## ENSG00000000938          0

In this count matrix, each row represents an Ensembl gene, each column a sequenced RNA library, and the values give the raw numbers of sequencing reads that were mapped to the respective gene in each library. We also have metadata on each of the samples (the columns of the count matrix). If you’ve counted reads with some other software, you need to check at this step that the columns of the count matrix correspond to the rows of the column metadata.

We can quickly check the millions of fragments which uniquely aligned to the genes.

colSums(assay(se))
## SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520 SRR1039521 
##   20637971   18809481   25348649   15163415   24448408   30818215   19126151   21164133

Supposing we have constructed a SummarizedExperiment using one of the methods described in the previous section, we now need to make sure that the object contains all the necessary information about the samples, i.e., a table with metadata on the count matrix’s columns stored in the colData slot:

colData(se)
## DataFrame with 8 rows and 9 columns
##            SampleName     cell      dex    albut        Run avgLength Experiment    Sample
##              <factor> <factor> <factor> <factor>   <factor> <integer>   <factor>  <factor>
## SRR1039508 GSM1275862   N61311    untrt    untrt SRR1039508       126  SRX384345 SRS508568
## SRR1039509 GSM1275863   N61311      trt    untrt SRR1039509       126  SRX384346 SRS508567
## SRR1039512 GSM1275866  N052611    untrt    untrt SRR1039512       126  SRX384349 SRS508571
## SRR1039513 GSM1275867  N052611      trt    untrt SRR1039513        87  SRX384350 SRS508572
## SRR1039516 GSM1275870  N080611    untrt    untrt SRR1039516       120  SRX384353 SRS508575
## SRR1039517 GSM1275871  N080611      trt    untrt SRR1039517       126  SRX384354 SRS508576
## SRR1039520 GSM1275874  N061011    untrt    untrt SRR1039520       101  SRX384357 SRS508579
## SRR1039521 GSM1275875  N061011      trt    untrt SRR1039521        98  SRX384358 SRS508580
##               BioSample
##                <factor>
## SRR1039508 SAMN02422669
## SRR1039509 SAMN02422675
## SRR1039512 SAMN02422678
## SRR1039513 SAMN02422670
## SRR1039516 SAMN02422682
## SRR1039517 SAMN02422673
## SRR1039520 SAMN02422683
## SRR1039521 SAMN02422677

Here we see that this object already contains an informative colData slot – because we have already prepared it for you, as described in the airway vignette. However, when you work with your own data, you will have to add the pertinent sample / phenotypic information for the experiment at this stage. We highly recommend keeping this information in a comma-separated value (CSV) or tab-separated value (TSV) file, which can be exported from an Excel spreadsheet, and the assign this to the colData slot, making sure that the rows correspond to the columns of the SummarizedExperiment. We made sure of this correspondence by specifying the BAM files using a column of the sample table.

Check out the rowRanges() of the summarized experiment; these are the genomic ranges over which counting occurred.

rowRanges(se)
## GRangesList object of length 64102:
## $ENSG00000000003 
## GRanges object with 17 ranges and 2 metadata columns:
##        seqnames               ranges strand |   exon_id       exon_name
##           <Rle>            <IRanges>  <Rle> | <integer>     <character>
##    [1]        X [99883667, 99884983]      - |    667145 ENSE00001459322
##    [2]        X [99885756, 99885863]      - |    667146 ENSE00000868868
##    [3]        X [99887482, 99887565]      - |    667147 ENSE00000401072
##    [4]        X [99887538, 99887565]      - |    667148 ENSE00001849132
##    [5]        X [99888402, 99888536]      - |    667149 ENSE00003554016
##    ...      ...                  ...    ... .       ...             ...
##   [13]        X [99890555, 99890743]      - |    667156 ENSE00003512331
##   [14]        X [99891188, 99891686]      - |    667158 ENSE00001886883
##   [15]        X [99891605, 99891803]      - |    667159 ENSE00001855382
##   [16]        X [99891790, 99892101]      - |    667160 ENSE00001863395
##   [17]        X [99894942, 99894988]      - |    667161 ENSE00001828996
## 
## ...
## <64101 more elements>
## -------
## seqinfo: 722 sequences (1 circular) from an unspecified genome

2.5 From SummarizedExperiment to DESeqDataSet

We will use the DESeq2 package for assessing differential expression. The package uses an extended version of the SummarizedExperiment class, called DESeqDataSet. It’s easy to go from a SummarizedExperiment to DESeqDataSet:

library("DESeq2")
dds <- DESeqDataSet(se, design = ~ cell + dex)

The ‘design’ argument is a formula which expresses how the counts for each gene depend on the variables in colData. Remember you can always get information on method arguments with ?, e.g ?DESeqDataSet.

2.6 Differential expression analysis

It will be convenient to make sure that untrt is the first level in the dex factor, so that the default log2 fold changes are calculated as treated over untreated (by default R will chose the first alphabetical level, remember: computers don’t know what to do unless you tell them). The function relevel() achieves this:

dds$dex <- relevel(dds$dex, "untrt")

In addition, if you have at any point subset the columns of the DESeqDataSet you should similarly call droplevels() on the factors if the subsetting has resulted in some levels having 0 samples.

2.6.1 Running the pipeline

Finally, we are ready to run the differential expression pipeline. With the data object prepared, the DESeq2 analysis can now be run with a single call to the function DESeq():

dds <- DESeq(dds)
## estimating size factors
## estimating dispersions
## gene-wise dispersion estimates
## mean-dispersion relationship
## final dispersion estimates
## fitting model and testing

This function will print out a message for the various steps it performs. These are described in more detail in the manual page ?DESeq. Briefly these are: the estimation of size factors (which control for differences in the library size of the sequencing experiments), the estimation of dispersion for each gene, and fitting a generalized linear model.

A DESeqDataSet is returned which contains all the fitted information within it, and the following section describes how to extract out result tables of interest from this object.

2.6.2 Building the results table

Calling results() without any arguments will extract the estimated log2 fold changes and p values for the last variable in the design formula. If there are more than 2 levels for this variable, results() will extract the results table for a comparison of the last level over the first level.

(res <- results(dds))
## log2 fold change (MAP): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 64102 rows and 6 columns
##                  baseMean log2FoldChange      lfcSE       stat       pvalue        padj
##                 <numeric>      <numeric>  <numeric>  <numeric>    <numeric>   <numeric>
## ENSG00000000003 708.60217    -0.37415246 0.09884435 -3.7852692 0.0001535422 0.001289269
## ENSG00000000005   0.00000             NA         NA         NA           NA          NA
## ENSG00000000419 520.29790     0.20206175 0.10974241  1.8412367 0.0655868755 0.197066699
## ENSG00000000457 237.16304     0.03616686 0.13834540  0.2614244 0.7937652378 0.913855995
## ENSG00000000460  57.93263    -0.08445399 0.24990710 -0.3379415 0.7354072485 0.884141561
## ...                   ...            ...        ...        ...          ...         ...
## LRG_94                  0             NA         NA         NA           NA          NA
## LRG_96                  0             NA         NA         NA           NA          NA
## LRG_97                  0             NA         NA         NA           NA          NA
## LRG_98                  0             NA         NA         NA           NA          NA
## LRG_99                  0             NA         NA         NA           NA          NA

As res is a DataFrame object, it carries metadata with information on the meaning of the columns:

mcols(res, use.names=TRUE)
## DataFrame with 6 rows and 2 columns
##                        type                               description
##                 <character>                               <character>
## baseMean       intermediate mean of normalized counts for all samples
## log2FoldChange      results  log2 fold change (MAP): dex trt vs untrt
## lfcSE               results          standard error: dex trt vs untrt
## stat                results          Wald statistic: dex trt vs untrt
## pvalue              results       Wald test p-value: dex trt vs untrt
## padj                results                      BH adjusted p-values

The first column, baseMean, is a just the average of the normalized count values, dividing by size factors, taken over all samples. The remaining four columns refer to a specific contrast, namely the comparison of the trt level over the untrt level for the factor variable dex. See the help page for results() (by typing ?results) for information on how to obtain other contrasts.

The column log2FoldChange is the effect size estimate. It tells us how much the gene’s expression seems to have changed due to treatment with dexamethasone in comparison to untreated samples. This value is reported on a logarithmic scale to base 2: for example, a log2 fold change of 1.5 means that the gene’s expression is increased by a multiplicative factor of \(2^{1.5} \approx 2.82\).

Of course, this estimate has an uncertainty associated with it, which is available in the column lfcSE, the standard error estimate for the log2 fold change estimate. We can also express the uncertainty of a particular effect size estimate as the result of a statistical test. The purpose of a test for differential expression is to test whether the data provides sufficient evidence to conclude that this value is really different from zero. DESeq2 performs for each gene a hypothesis test to see whether evidence is sufficient to decide against the null hypothesis that there is no effect of the treatment on the gene and that the observed difference between treatment and control was merely caused by experimental variability (i.e., the type of variability that you can just as well expect between different samples in the same treatment group). As usual in statistics, the result of this test is reported as a p value, and it is found in the column pvalue. (Remember that a p value indicates the probability that a fold change as strong as the observed one, or even stronger, would be seen under the situation described by the null hypothesis.)

We can also summarize the results with the following line of code, which reports some additional information.

summary(res)
## 
## out of 33469 with nonzero total read count
## adjusted p-value < 0.1
## LFC > 0 (up)     : 2617, 7.8% 
## LFC < 0 (down)   : 2203, 6.6% 
## outliers [1]     : 0, 0% 
## low counts [2]   : 15441, 46% 
## (mean count < 5)
## [1] see 'cooksCutoff' argument of ?results
## [2] see 'independentFiltering' argument of ?results

Note that there are many genes with differential expression due to dexamethasone treatment at the FDR level of 10%. This makes sense, as the smooth muscle cells of the airway are known to react to glucocorticoid steroids. However, there are two ways to be more strict about which set of genes are considered significant:

  • lower the false discovery rate threshold (the threshold on padj in the results table)
  • raise the log2 fold change threshold from 0 using the lfcThreshold argument of results(). See the DESeq2 vignette for a demonstration of the use of this argument.

Sometimes a subset of the p values in res will be NA (“not available”). This is DESeq()’s way of reporting that all counts for this gene were zero, and hence not test was applied. In addition, p values can be assigned NA if the gene was excluded from analysis because it contained an extreme count outlier. For more information, see the outlier detection section of the vignette.

2.6.3 Multiple testing

Novices in high-throughput biology often assume that thresholding these p values at a low value, say 0.05, as is often done in other settings, would be appropriate – but it is not. We briefly explain why:

There are 5648 genes with a p value below 0.05 among the 33469 genes, for which the test succeeded in reporting a p value:

sum(res$pvalue < 0.05, na.rm=TRUE)
## [1] 5648
sum(!is.na(res$pvalue))
## [1] 33469

Now, assume for a moment that the null hypothesis is true for all genes, i.e., no gene is affected by the treatment with dexamethasone. Then, by the definition of p value, we expect up to 5% of the genes to have a p value below 0.05. This amounts to 1673 genes. If we just considered the list of genes with a p value below 0.05 as differentially expressed, this list should therefore be expected to contain up to 1673 / 5648 = 30% false positives.

DESeq2 uses the Benjamini-Hochberg (BH) adjustment as described in the base R p.adjust function; in brief, this method calculates for each gene an adjusted p value which answers the following question: if one called significant all genes with a p value less than or equal to this gene’s p value threshold, what would be the fraction of false positives (the false discovery rate, FDR) among them (in the sense of the calculation outlined above)? These values, called the BH-adjusted p values, are given in the column padj of the res object.

Hence, if we consider a fraction of 10% false positives acceptable, we can consider all genes with an adjusted p value below \(10% = 0.1\) as significant. How many such genes are there?

sum(res$padj < 0.1, na.rm=TRUE)
## [1] 4820

We subset the results table to these genes and then sort it by the log2 fold change estimate to get the significant genes with the strongest down-regulation.

resSig <- subset(res, padj < 0.1)
head(resSig[ order( resSig$log2FoldChange ), ])
## log2 fold change (MAP): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 6 columns
##                  baseMean log2FoldChange     lfcSE       stat       pvalue         padj
##                 <numeric>      <numeric> <numeric>  <numeric>    <numeric>    <numeric>
## ENSG00000162692 508.17023      -3.449475 0.1767138 -19.520120 7.406383e-85 9.537305e-82
## ENSG00000105989 333.21469      -2.847364 0.1763098 -16.149771 1.139834e-58 5.871121e-56
## ENSG00000146006  46.80760      -2.828093 0.3377003  -8.374564 5.542861e-17 2.708041e-15
## ENSG00000214814 243.27698      -2.753559 0.2235537 -12.317215 7.317772e-35 1.232942e-32
## ENSG00000267339  26.23357      -2.704476 0.3519722  -7.683776 1.544666e-14 5.924944e-13
## ENSG00000013293 244.49733      -2.641050 0.1992872 -13.252481 4.365309e-40 8.842448e-38

…and with the strongest upregulation. The order() function gives the indices in increasing order, so a simple way to ask for decreasing order is to add a - sign. Alternatively, you can use the argument decreasing=TRUE.

head(resSig[ order( -resSig$log2FoldChange ), ])
## log2 fold change (MAP): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 6 columns
##                  baseMean log2FoldChange     lfcSE      stat        pvalue          padj
##                 <numeric>      <numeric> <numeric> <numeric>     <numeric>     <numeric>
## ENSG00000109906 385.07103       4.847164 0.3313657  14.62784  1.866158e-48  5.902298e-46
## ENSG00000179593  67.24305       4.830815 0.3314192  14.57615  3.983670e-48  1.196960e-45
## ENSG00000152583 997.43977       4.313961 0.1721375  25.06114 1.320082e-138 2.379844e-134
## ENSG00000163884 561.10717       4.074297 0.2104708  19.35802  1.744626e-83  1.965757e-80
## ENSG00000250978  56.31819       4.054731 0.3294746  12.30666  8.340703e-35  1.392280e-32
## ENSG00000168309 159.52692       3.977191 0.2558532  15.54481  1.725211e-54  7.775524e-52

2.7 Diagnostic plots

A quick way to visualize the counts for a particular gene is to use the plotCounts() function, which takes as arguments the DESeqDataSet, a gene name, and the group over which to plot the counts.

topGene <- rownames(res)[which.min(res$padj)]
data <- plotCounts(dds, gene=topGene, intgroup=c("dex"), returnData=TRUE)

We can also make more customizable plots using the ggplot() function from the ggplot2 package:

library(ggplot2)
ggplot(data, aes(x=dex, y=count, fill=dex)) +
  scale_y_log10() + 
  geom_dotplot(binaxis="y", stackdir="center")
## `stat_bindot()` using `bins = 30`. Pick better value with `binwidth`.

An “MA-plot” provides a useful overview for an experiment with a two-group comparison. The log2 fold change for a particular comparison is plotted on the y-axis and the average of the counts normalized by size factor is shown on the x-axis (“M” for minus, because a log ratio is equal to log minus log, and “A” for average).

plotMA(res, ylim=c(-5,5))

Each gene is represented with a dot. Genes with an adjusted \(p\) value below a threshold (here 0.1, the default) are shown in red. The DESeq2 package incorporates a prior on log2 fold changes, resulting in moderated log2 fold changes from genes with low counts and highly variable counts, as can be seen by the narrowing of spread of points on the left side of the plot. This plot demonstrates that only genes with a large average normalized count contain sufficient information to yield a significant call.

We can label individual points on the MA plot as well. Here we use the with() R function to plot a circle and text for a selected row of the results object. Within the with() function, only the baseMean and log2FoldChange values for the selected rows of res are used.

plotMA(res, ylim=c(-5,5))
with(res[topGene, ], {
  points(baseMean, log2FoldChange, col="dodgerblue", cex=2, lwd=2)
  text(baseMean, log2FoldChange, topGene, pos=2, col="dodgerblue")
})

Whether a gene is called significant depends not only on its LFC but also on its within-group variability, which DESeq2 quantifies as the dispersion. For strongly expressed genes, the dispersion can be understood as a squared coefficient of variation: a dispersion value of 0.01 means that the gene’s expression tends to differ by typically \(\sqrt{0.01} = 10\%\) between samples of the same treatment group. For weak genes, the Poisson noise is an additional source of noise.

The function plotDispEsts() visualizes DESeq2’s dispersion estimates:

plotDispEsts(dds)

The black points are the dispersion estimates for each gene as obtained by considering the information from each gene separately. Unless one has many samples, these values fluctuate strongly around their true values. Therefore, we fit the red trend line, which shows the dispersions’ dependence on the mean, and then shrink each gene’s estimate towards the red line to obtain the final estimates (blue points) that are then used in the hypothesis test. The blue circles above the main “cloud” of points are genes which have high gene-wise dispersion estimates which are labelled as dispersion outliers. These estimates are therefore not shrunk toward the fitted trend line.

Another useful diagnostic plot is the histogram of the p values.

hist(res$pvalue, breaks=20, col="grey50", border="white")

This plot becomes a bit smoother by excluding genes with very small counts:

hist(res$pvalue[res$baseMean > 1], breaks=20, col="grey50", border="white")

2.8 Independent filtering

The MA plot highlights an important property of RNA-Seq data. For weakly expressed genes, we have no chance of seeing differential expression, because the low read counts suffer from so high Poisson noise that any biological effect is drowned in the uncertainties from the read counting. We can also show this by examining the ratio of small p values (say, less than, 0.01) for genes binned by mean normalized count:

# create bins using the quantile function
qs <- c(0, quantile(res$baseMean[res$baseMean > 0], 0:7/7))
# cut the genes into the bins
bins <- cut(res$baseMean, qs)
# rename the levels of the bins using the middle point
levels(bins) <- paste0("~",round(.5*qs[-1] + .5*qs[-length(qs)]))
# calculate the ratio of $p$ values less than .01 for each bin
ratios <- tapply(res$pvalue, bins, function(p) mean(p < .01, na.rm=TRUE))
# plot these ratios
barplot(ratios, xlab="mean normalized count", ylab="ratio of small p values")

At first sight, there may seem to be little benefit in filtering out these genes. After all, the test found them to be non-significant anyway. However, these genes have an influence on the multiple testing adjustment, whose performance improves if such genes are removed. By removing the weakly-expressed genes from the input to the FDR procedure, we can find more genes to be significant among those which we keep, and so improved the power of our test. This approach is known as independent filtering.

The term independent highlights an important caveat. Such filtering is permissible only if the filter criterion is independent of the actual test statistic. Otherwise, the filtering would invalidate the test and consequently the assumptions of the BH procedure. This is why we filtered on the average over all samples: this filter is blind to the assignment of samples to the treatment and control group and hence independent. The independent filtering software used inside DESeq2 comes from the genefilter package, which contains a reference to a paper describing the statistical foundation for independent filtering.

2.9 Annotation: adding gene names

Our result table only uses Ensembl gene IDs, but gene names may be more informative. Bioconductor’s annotation packages help with mapping various ID schemes to each other.

We load the AnnotationDbi package and the annotation package org.Hs.eg.db:

library(org.Hs.eg.db)

This is the organism annotation package (“org”) for Homo sapiens (“Hs”), organized as an AnnotationDbi database package (“db”), using Entrez Gene IDs (“eg”) as primary key. To get a list of all available key types, use:

columns(org.Hs.eg.db)
##  [1] "ACCNUM"       "ALIAS"        "ENSEMBL"      "ENSEMBLPROT"  "ENSEMBLTRANS" "ENTREZID"    
##  [7] "ENZYME"       "EVIDENCE"     "EVIDENCEALL"  "GENENAME"     "GO"           "GOALL"       
## [13] "IPI"          "MAP"          "OMIM"         "ONTOLOGY"     "ONTOLOGYALL"  "PATH"        
## [19] "PFAM"         "PMID"         "PROSITE"      "REFSEQ"       "SYMBOL"       "UCSCKG"      
## [25] "UNIGENE"      "UNIPROT"
res$hgnc_symbol <- 
    unname(mapIds(org.Hs.eg.db, rownames(res), "SYMBOL", "ENSEMBL"))
## 'select()' returned 1:many mapping between keys and columns
res$entrezgene <- 
    unname(mapIds(org.Hs.eg.db, rownames(res), "ENTREZID", "ENSEMBL"))
## 'select()' returned 1:many mapping between keys and columns

Now the results have the desired external gene ids:

resOrdered <- res[order(res$pvalue),]
head(resOrdered)
## log2 fold change (MAP): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 8 columns
##                   baseMean log2FoldChange     lfcSE      stat        pvalue          padj
##                  <numeric>      <numeric> <numeric> <numeric>     <numeric>     <numeric>
## ENSG00000152583   997.4398       4.313961 0.1721375  25.06114 1.320082e-138 2.379844e-134
## ENSG00000165995   495.0929       3.186823 0.1281565  24.86665 1.708459e-136 1.540005e-132
## ENSG00000101347 12703.3871       3.618735 0.1489428  24.29613 2.153705e-130 1.294233e-126
## ENSG00000120129  3409.0294       2.871488 0.1182491  24.28337 2.937761e-130 1.324049e-126
## ENSG00000189221  2341.7673       3.230395 0.1366746  23.63567 1.657244e-123 5.975359e-120
## ENSG00000211445 12285.6151       3.553360 0.1579821  22.49216 4.952520e-112 1.488067e-108
##                 hgnc_symbol  entrezgene
##                 <character> <character>
## ENSG00000152583     SPARCL1        8404
## ENSG00000165995      CACNB2         783
## ENSG00000101347      SAMHD1       25939
## ENSG00000120129       DUSP1        1843
## ENSG00000189221        MAOA        4128
## ENSG00000211445        GPX3        2878

2.10 Exporting results

You can easily save the results table in a CSV file, which you can then load with a spreadsheet program such as Excel. The call to as.data.frame is necessary to convert the DataFrame object (IRanges package) to a data.frame object which can be processed by write.csv.

write.csv(as.data.frame(resOrdered), file="results.csv")

2.11 Session information

As 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 as it will help to trace down what has happened in case that an R script ceases to work because the functions have been changed in a newer version of a package. The session information should also always be included in any emails to the Bioconductor support site along with all code used in the analysis.

sessionInfo()
## R version 3.3.0 Patched (2016-05-17 r70630)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 14.04.4 LTS
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C               LC_TIME=en_US.UTF-8       
##  [4] LC_COLLATE=C               LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                  LC_ADDRESS=C              
## [10] LC_TELEPHONE=C             LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] stats4    parallel  stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
##  [1] gplots_3.0.1                            DESeq2_1.12.2                          
##  [3] rtracklayer_1.32.0                      RNAseqData.HNRNPC.bam.chr14_0.10.0     
##  [5] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2 GenomicFeatures_1.24.2                 
##  [7] org.Hs.eg.db_3.3.0                      AnnotationDbi_1.34.2                   
##  [9] airway_0.106.0                          biomaRt_2.28.0                         
## [11] AnnotationHub_2.4.2                     GenomicAlignments_1.8.0                
## [13] Rsamtools_1.24.0                        SummarizedExperiment_1.2.2             
## [15] Biobase_2.32.0                          GenomicRanges_1.24.0                   
## [17] GenomeInfoDb_1.8.1                      Biostrings_2.40.0                      
## [19] XVector_0.12.0                          IRanges_2.6.0                          
## [21] S4Vectors_0.10.0                        BiocGenerics_0.18.0                    
## [23] ggplot2_2.1.0                           class_7.3-14                           
## [25] RColorBrewer_1.1-2                      BiocStyle_2.0.2                        
## 
## loaded via a namespace (and not attached):
##  [1] httr_1.1.0                    splines_3.3.0                 gtools_3.5.0                 
##  [4] Formula_1.2-1                 shiny_0.13.2                  interactiveDisplayBase_1.10.3
##  [7] latticeExtra_0.6-28           yaml_2.1.13                   RSQLite_1.0.0                
## [10] lattice_0.20-33               chron_2.3-47                  digest_0.6.9                 
## [13] colorspace_1.2-6              htmltools_0.3.5               httpuv_1.3.3                 
## [16] Matrix_1.2-6                  plyr_1.8.3                    XML_3.98-1.4                 
## [19] genefilter_1.54.2             zlibbioc_1.18.0               xtable_1.8-2                 
## [22] scales_0.4.0                  gdata_2.17.0                  BiocParallel_1.6.2           
## [25] annotate_1.50.0               nnet_7.3-12                   survival_2.39-4              
## [28] magrittr_1.5                  mime_0.4                      evaluate_0.9                 
## [31] foreign_0.8-66                BiocInstaller_1.22.2          tools_3.3.0                  
## [34] data.table_1.9.6              formatR_1.4                   stringr_1.0.0                
## [37] locfit_1.5-9.1                munsell_0.4.3                 cluster_2.0.4                
## [40] caTools_1.17.1                grid_3.3.0                    RCurl_1.95-4.8               
## [43] bitops_1.0-6                  labeling_0.3                  rmarkdown_0.9.6              
## [46] gtable_0.2.0                  codetools_0.2-14              DBI_0.4-1                    
## [49] curl_0.9.7                    reshape2_1.4.1                R6_2.1.2                     
## [52] gridExtra_2.2.1               knitr_1.13                    Hmisc_3.17-4                 
## [55] KernSmooth_2.23-15            stringi_1.0-1                 Rcpp_0.12.5                  
## [58] geneplotter_1.50.0            rpart_4.1-10                  acepack_1.3-3.3