Contents

1 Abstract

We walk through an end-to-end gene-level RNA-seq differential expression workflow using Bioconductor packages. We will start from the FASTQ files, show how these were aligned to the reference genome, and prepare a count matrix which tallies the number of RNA-seq reads/fragments within each gene for each sample. We will perform exploratory data analysis (EDA) for quality assessment and to explore the relationship between samples, perform differential gene expression analysis and visually explore the results.

2 Introduction

Bioconductor has many packages which support analysis of high-throughput sequence data, including RNA sequencing (RNA-seq). The packages which we will use in this workflow include core packages maintained by the Bioconductor core team for importing and processing raw sequencing data and loading gene annotations. We will also use contributed packages for statistical analysis and visualization of sequencing data. Through scheduled releases every 6 months, the Bioconductor project ensures that all the packages within a release will work together in harmony (hence the “conductor” metaphor). The packages used in this workflow are loaded with the library function and can be installed by following the Bioconductor package installation instructions.

A published (but essentially similar) version of this workflow, including reviewer reports and comments is available at F1000Research.

If you have questions about this workflow or any Bioconductor software, please post these to the Bioconductor support site. If the questions concern a specific package, you can tag the post with the name of the package, or for general questions about the workflow, tag the post with rnaseqgene. Note the posting guide for crafting an optimal question for the support site.

2.1 Experimental data

The data used in this workflow is stored in the airway package that summarizes an RNA-seq experiment wherein airway smooth muscle cells were treated with dexamethasone, a synthetic glucocorticoid steroid with anti-inflammatory effects (Himes et al. 2014). Glucocorticoids are used, for example, by people with asthma to 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. For more description of the experiment see the PubMed entry 24926665 and for raw data see the GEO entry GSE52778.

3 Preparing count matrices

As input, the count-based statistical methods, such as DESeq2 (M. I. Love, Huber, and Anders 2014), edgeR (M. D. Robinson, McCarthy, and Smyth 2009), limma with the voom method (Law et al. 2014), DSS (H. Wu, Wang, and Wu 2013), EBSeq (Leng et al. 2013) and baySeq (Hardcastle and Kelly 2010), expect input 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 (or fragments, for paired-end RNA-seq) have been assigned 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), species of bacteria (with metagenomic datasets), or peptide sequences (with quantitative mass spectrometry).

The values in the matrix should be counts of sequencing reads/fragments. This is important for the statistical models used by DESeq2 and edgeR to hold, as only counts allow assessing the measurement precision correctly. It is important to never provide counts that were pre-normalized for sequencing depth/library size, as the statistical model is most powerful when applied to un-normalized counts and is designed to account for library size differences internally.

3.2 Aligning reads to a reference genome

The computational analysis of an RNA-seq experiment begins from the FASTQ files that contain the nucleotide sequence of each read and a quality score at each position. These reads must first be aligned to a reference genome or transcriptome, or the abundances and estimated counts per transcript can be estimated without alignment, as described above. In either case, it is important to know if the sequencing experiment was single-end or paired-end, as the alignment software will require the user to specify both FASTQ files for a paired-end experiment. The output of this alignment step is commonly stored in a file format called SAM/BAM.

A number of software programs exist to align reads to a reference genome, and we encourage you to check out some of the benchmarking papers that discuss the advantages and disadvantages of each software, which include accuracy, sensitivity in aligning reads over splice junctions, speed, memory required, usability and many other features. Here, we used the STAR read aligner (Dobin et al. 2013) to align the reads for our current experiment to the Ensembl release 75 (Flicek et al. 2014) human reference genome. In the following code example, it is assumed that there is a file in the current directory called files with each line containing an identifier for each experiment, and we have all the FASTQ files in a subdirectory fastq. If you have downloaded the FASTQ files from the Sequence Read Archive, the identifiers would be SRA run IDs, e.g. SRR1039520. You should have two files for a paired-end experiment for each ID, fastq/SRR1039520_1.fastq1 and fastq/SRR1039520_2.fastq, which give the first and second read for the paired-end fragments. If you have performed a single-end experiment, you would only have one file per ID. We have also created a subdirectory, aligned, where STAR will output its alignment files.

for f in `cat files`; do STAR \
--genomeDir ../STAR/ENSEMBL.homo_sapiens.release-75 \
--readFilesIn fastq/$f\_1.fastq fastq/$f\_2.fastq \
--runThreadN 12 --outFileNamePrefix aligned/$f.; done

SAMtools (H. Li et al. 2009) was used to generate BAM files. The -@ flag can be used to allocate additional threads.

for f in `cat files`; do samtools view -bS aligned/$f.Aligned.out.sam \
-o aligned/$f.bam; done

The BAM files for a number of sequencing runs can then be used to generate count matrices, as described in the following section.

3.3 Locating alignment files

Besides the count matrix that we will use later, the airway package also contains eight files with a small subset of the total number of reads in the experiment. The reads were selected which aligned to a small region of chromosome 1. Here, for demonstration, we chose a subset of reads because the full alignment files are large (a few gigabytes each), and because it takes between 10-30 minutes to count the fragments for each sample. We will use these files to demonstrate how a count matrix can be constructed from BAM files. Afterwards, we will load the full count matrix corresponding to all samples and all data, which is already provided in the same package, and will continue the analysis with that full matrix.

We load the data package with the example data:

library("airway")

The R function system.file can be used to find out where on your computer the files from a package have been installed. Here we ask for the full path to the extdata directory, where R packages store external data, that is part of the airway package.

indir <- system.file("extdata", package = "airway", mustWork = TRUE)

In this directory, we find the eight BAM files (and some other files):

list.files(indir)
##  [1] "GSE52778_series_matrix.txt"        "Homo_sapiens.GRCh37.75_subset.gtf"
##  [3] "SRR1039508_subset.bam"             "SRR1039509_subset.bam"            
##  [5] "SRR1039512_subset.bam"             "SRR1039513_subset.bam"            
##  [7] "SRR1039516_subset.bam"             "SRR1039517_subset.bam"            
##  [9] "SRR1039520_subset.bam"             "SRR1039521_subset.bam"            
## [11] "SraRunInfo_SRP033351.csv"          "sample_table.csv"

Typically, we have a table with detailed information for each of our samples that links samples to the associated FASTQ and BAM files. For your own project, you might create such a comma-separated value (CSV) file using a text editor or spreadsheet software such as Excel.

We load such a CSV file with read.csv:

csvfile <- file.path(indir, "sample_table.csv")
sampleTable <- read.csv(csvfile, row.names = 1)
sampleTable
##            SampleName    cell   dex albut        Run avgLength Experiment
## SRR1039508 GSM1275862  N61311 untrt untrt SRR1039508       126  SRX384345
## SRR1039509 GSM1275863  N61311   trt untrt SRR1039509       126  SRX384346
## SRR1039512 GSM1275866 N052611 untrt untrt SRR1039512       126  SRX384349
## SRR1039513 GSM1275867 N052611   trt untrt SRR1039513        87  SRX384350
## SRR1039516 GSM1275870 N080611 untrt untrt SRR1039516       120  SRX384353
## SRR1039517 GSM1275871 N080611   trt untrt SRR1039517       126  SRX384354
## SRR1039520 GSM1275874 N061011 untrt untrt SRR1039520       101  SRX384357
## SRR1039521 GSM1275875 N061011   trt untrt SRR1039521        98  SRX384358
##               Sample    BioSample
## SRR1039508 SRS508568 SAMN02422669
## SRR1039509 SRS508567 SAMN02422675
## SRR1039512 SRS508571 SAMN02422678
## SRR1039513 SRS508572 SAMN02422670
## SRR1039516 SRS508575 SAMN02422682
## SRR1039517 SRS508576 SAMN02422673
## SRR1039520 SRS508579 SAMN02422683
## SRR1039521 SRS508580 SAMN02422677

Once the reads have been aligned, there are a number of tools that can be used to count the number of reads/fragments that can be assigned to genomic features for each sample. These often take as input SAM/BAM alignment files and a file specifying the genomic features, e.g. a GFF3 or GTF file specifying the gene models.

3.4 DESeq2 import functions

The following tools can be used generate count matrices: summarizeOverlaps (Lawrence et al. 2013), featureCounts (Liao, Smyth, and Shi 2014), tximport (Soneson, Love, and Robinson 2015), htseq-count (Anders, Pyl, and Huber 2015).

function package framework output DESeq2 input function
summarizeOverlaps GenomicAlignments R/Bioconductor SummarizedExperiment DESeqDataSet
featureCounts Rsubread R/Bioconductor matrix DESeqDataSetFromMatrix
tximport tximport R/Bioconductor list of matrices DESeqDataSetFromTximport
htseq-count HTSeq Python files DESeqDataSetFromHTSeq

We now proceed using summarizeOverlaps. Using the Run column in the sample table, we construct the full paths to the files we want to perform the counting operation on:

filenames <- file.path(indir, paste0(sampleTable$Run, "_subset.bam"))
file.exists(filenames)
## [1] TRUE TRUE TRUE TRUE TRUE TRUE TRUE TRUE

We indicate in Bioconductor that these files are BAM files using the BamFileList function from the Rsamtools package that provides an R interface to BAM files. Here we also specify details about how the BAM files should be treated, e.g., only process 2 million reads at a time. See ?BamFileList for more information.

library("Rsamtools")
bamfiles <- BamFileList(filenames, yieldSize = 2000000)

Note: make sure that the chromosome names of the genomic features in the annotation you use are consistent with the chromosome names of the reference used for read alignment. Otherwise, the scripts might fail to count any reads to features due to the mismatching names. For example, a common mistake is when the alignment files contain chromosome names in the style of 1 and the gene annotation in the style of chr1, or the other way around. See the seqlevelsStyle function in the GenomeInfoDb package for solutions. We can check the chromosome names (here called “seqnames”) in the alignment files like so:

seqinfo(bamfiles[1])
## Seqinfo object with 84 sequences from an unspecified genome:
##   seqnames   seqlengths isCircular genome
##   1           249250621       <NA>   <NA>
##   10          135534747       <NA>   <NA>
##   11          135006516       <NA>   <NA>
##   12          133851895       <NA>   <NA>
##   13          115169878       <NA>   <NA>
##   ...               ...        ...    ...
##   GL000210.1      27682       <NA>   <NA>
##   GL000231.1      27386       <NA>   <NA>
##   GL000229.1      19913       <NA>   <NA>
##   GL000226.1      15008       <NA>   <NA>
##   GL000207.1       4262       <NA>   <NA>

3.5 Defining gene models

Next, we need to read in the gene model that will be used for counting reads/fragments. We will read the gene model from an Ensembl GTF file (Flicek et al. 2014), using makeTxDbFromGFF from the GenomicFeatures package. GTF files can be downloaded from Ensembl’s FTP site or other gene model repositories. A TxDb object is a database that can be used to generate a variety of range-based objects, such as exons, transcripts, and genes. We want to make a list of exons grouped by gene for counting read/fragments.

There are other options for constructing a TxDb. For the known genes track from the UCSC Genome Browser (Kent et al. 2002), one can use the pre-built Transcript DataBase: TxDb.Hsapiens.UCSC.hg19.knownGene. If the annotation file is accessible from AnnotationHub (as is the case for the Ensembl genes), a pre-scanned GTF file can be imported using makeTxDbFromGRanges.

Here we will demonstrate loading from a GTF file:

library("GenomicFeatures")

We indicate that none of our sequences (chromosomes) are circular using a 0-length character vector.

gtffile <- file.path(indir,"Homo_sapiens.GRCh37.75_subset.gtf")
txdb <- makeTxDbFromGFF(gtffile, format = "gtf", circ_seqs = character())
txdb
## TxDb object:
## # Db type: TxDb
## # Supporting package: GenomicFeatures
## # Data source: /Library/Frameworks/R.framework/Versions/3.4/Resources/library/airway/extdata/Homo_sapiens.GRCh37.75_subset.gtf
## # Organism: NA
## # Taxonomy ID: NA
## # miRBase build ID: NA
## # Genome: NA
## # transcript_nrow: 65
## # exon_nrow: 279
## # cds_nrow: 158
## # Db created by: GenomicFeatures package from Bioconductor
## # Creation time: 2017-06-11 09:38:46 +0200 (Sun, 11 Jun 2017)
## # GenomicFeatures version at creation time: 1.28.3
## # RSQLite version at creation time: 1.1-2
## # DBSCHEMAVERSION: 1.1

The following line produces a GRangesList of all the exons grouped by gene (Lawrence et al. 2013). Each element of the list is a GRanges object of the exons for a gene.

ebg <- exonsBy(txdb, by = "gene")
ebg
## GRangesList object of length 20:
## $ENSG00000009724 
## GRanges object with 18 ranges and 2 metadata columns:
##        seqnames               ranges strand |   exon_id       exon_name
##           <Rle>            <IRanges>  <Rle> | <integer>     <character>
##    [1]        1 [11086580, 11087705]      - |        98 ENSE00000818830
##    [2]        1 [11090233, 11090307]      - |        99 ENSE00000472123
##    [3]        1 [11090805, 11090939]      - |       100 ENSE00000743084
##    [4]        1 [11094885, 11094963]      - |       101 ENSE00000743085
##    [5]        1 [11097750, 11097868]      - |       102 ENSE00003482788
##    ...      ...                  ...    ... .       ...             ...
##   [14]        1 [11106948, 11107176]      - |       111 ENSE00003467404
##   [15]        1 [11106948, 11107176]      - |       112 ENSE00003489217
##   [16]        1 [11107260, 11107280]      - |       113 ENSE00001833377
##   [17]        1 [11107260, 11107284]      - |       114 ENSE00001472289
##   [18]        1 [11107260, 11107290]      - |       115 ENSE00001881401
## 
## ...
## <19 more elements>
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths

3.6 Read counting step

After these preparations, the actual counting is easy. The function summarizeOverlaps from the GenomicAlignments package will do this. This produces a SummarizedExperiment object that contains a variety of information about the experiment and will be described in more detail below.

Note: If it is desired to perform counting using multiple cores, one can use the register and MulticoreParam or SnowParam functions from the BiocParallel package before the counting call below. Expect that the summarizeOverlaps call will take at least 30 minutes per file for a human RNA-seq file with 30 million aligned reads. By sending the files to separate cores, one can speed up the entire counting process.

library("GenomicAlignments")
library("BiocParallel")

Here we specify to use one core, not multiple cores. We could have also skipped this line and the counting step would run in serial.

register(SerialParam())

The following call creates the SummarizedExperiment object with counts:

se <- summarizeOverlaps(features = ebg, reads = bamfiles,
                        mode = "Union",
                        singleEnd = FALSE,
                        ignore.strand = TRUE,
                        fragments = TRUE)

We specify a number of arguments besides the features and the reads. The mode argument describes what kind of read overlaps will be counted. These modes are shown in Figure 1 of the Counting reads with summarizeOverlaps vignette for the GenomicAlignments package. Note that fragments will be counted only once to each gene, even if they overlap multiple exons of a gene which may themselves be overlapping. Setting singleEnd to FALSE indicates that the experiment produced paired-end reads, and we want to count a pair of reads (a fragment) only once toward the count for a gene. The fragments argument can be used when singleEnd = FALSE to specify if unpaired reads should be counted (yes if fragments = TRUE).

In order to produce correct counts, it is important to know if the RNA-seq experiment was strand-specific or not. This experiment was not strand-specific so we set ignore.strand to TRUE. However, certain strand-specific protocols could have the reads align only to the opposite strand of the genes. The user must check if the experiment was strand-specific and if so, whether the reads should align to the forward or reverse strand of the genes. For various counting/quantifying tools, one specifies counting on the forward or reverse strand in different ways, although this task is currently easiest with htseq-count, featureCounts, or the transcript abundance quantifiers mentioned previously. It is always a good idea to check the column sums of the count matrix (see below) to make sure these totals match the expected of the number of reads or fragments aligning to genes. Additionally, one can visually check the read alignments using a genome visualization tool.

3.7 SummarizedExperiment

The component parts of a SummarizedExperiment object. The assay (pink block) contains the matrix of counts, the rowRanges (blue block) contains information about the genomic ranges and the colData (green block) contains information about the samples. The highlighted line in each block represents the first row (note that the first row of colData lines up with the first column of the assay).

The SummarizedExperiment container is diagrammed in the Figure above and discussed in the latest Bioconductor paper (Huber et al. 2015). In our case we have created a single matrix named “counts” that contains the fragment counts for each gene and sample, which is stored in assay. It is also possible to store multiple matrices, accessed with assays. The rowRanges for our object is the GRangesList we used for counting (one GRanges of exons for each row of the count matrix). The component parts of the SummarizedExperiment are accessed with an R function of the same name: assay (or assays), rowRanges and colData.

This example code above actually only counted a small subset of fragments from the original experiment. Nevertheless, we can still investigate the resulting SummarizedExperiment by looking at the counts in the assay slot, the phenotypic data about the samples in colData slot (in this case an empty DataFrame) and the data about the genes in the rowRanges slot.

se
## class: RangedSummarizedExperiment 
## dim: 20 8 
## metadata(0):
## assays(1): counts
## rownames(20): ENSG00000009724 ENSG00000116649 ... ENSG00000271794
##   ENSG00000271895
## rowData names(0):
## colnames(8): SRR1039508_subset.bam SRR1039509_subset.bam ...
##   SRR1039520_subset.bam SRR1039521_subset.bam
## colData names(0):
dim(se)
## [1] 20  8
assayNames(se)
## [1] "counts"
head(assay(se), 3)
##                 SRR1039508_subset.bam SRR1039509_subset.bam
## ENSG00000009724                    38                    28
## ENSG00000116649                  1004                  1255
## ENSG00000120942                   218                   256
##                 SRR1039512_subset.bam SRR1039513_subset.bam
## ENSG00000009724                    66                    24
## ENSG00000116649                  1122                  1313
## ENSG00000120942                   233                   252
##                 SRR1039516_subset.bam SRR1039517_subset.bam
## ENSG00000009724                    42                    41
## ENSG00000116649                  1100                  1879
## ENSG00000120942                   269                   465
##                 SRR1039520_subset.bam SRR1039521_subset.bam
## ENSG00000009724                    47                    36
## ENSG00000116649                   745                  1536
## ENSG00000120942                   207                   400
colSums(assay(se))
## SRR1039508_subset.bam SRR1039509_subset.bam SRR1039512_subset.bam 
##                  6478                  6501                  7699 
## SRR1039513_subset.bam SRR1039516_subset.bam SRR1039517_subset.bam 
##                  6801                  8009                 10849 
## SRR1039520_subset.bam SRR1039521_subset.bam 
##                  5254                  9168

The rowRanges, when printed, only shows the first GRanges, and tells us there are 19 more elements:

rowRanges(se)
## GRangesList object of length 20:
## $ENSG00000009724 
## GRanges object with 18 ranges and 2 metadata columns:
##        seqnames               ranges strand |   exon_id       exon_name
##           <Rle>            <IRanges>  <Rle> | <integer>     <character>
##    [1]        1 [11086580, 11087705]      - |        98 ENSE00000818830
##    [2]        1 [11090233, 11090307]      - |        99 ENSE00000472123
##    [3]        1 [11090805, 11090939]      - |       100 ENSE00000743084
##    [4]        1 [11094885, 11094963]      - |       101 ENSE00000743085
##    [5]        1 [11097750, 11097868]      - |       102 ENSE00003482788
##    ...      ...                  ...    ... .       ...             ...
##   [14]        1 [11106948, 11107176]      - |       111 ENSE00003467404
##   [15]        1 [11106948, 11107176]      - |       112 ENSE00003489217
##   [16]        1 [11107260, 11107280]      - |       113 ENSE00001833377
##   [17]        1 [11107260, 11107284]      - |       114 ENSE00001472289
##   [18]        1 [11107260, 11107290]      - |       115 ENSE00001881401
## 
## ...
## <19 more elements>
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths

The rowRanges also contains metadata about the construction of the gene model in the metadata slot. Here we use a helpful R function, str, to display the metadata compactly:

str(metadata(rowRanges(se)))
## List of 1
##  $ genomeInfo:List of 15
##   ..$ Db type                                 : chr "TxDb"
##   ..$ Supporting package                      : chr "GenomicFeatures"
##   ..$ Data source                             : chr "/Library/Frameworks/R.framework/Versions/3.4/Resources/library/airway/extdata/Homo_sapiens.GRCh37.75_subset.gtf"
##   ..$ Organism                                : chr NA
##   ..$ Taxonomy ID                             : chr NA
##   ..$ miRBase build ID                        : chr NA
##   ..$ Genome                                  : chr NA
##   ..$ transcript_nrow                         : chr "65"
##   ..$ exon_nrow                               : chr "279"
##   ..$ cds_nrow                                : chr "158"
##   ..$ Db created by                           : chr "GenomicFeatures package from Bioconductor"
##   ..$ Creation time                           : chr "2017-06-11 09:38:46 +0200 (Sun, 11 Jun 2017)"
##   ..$ GenomicFeatures version at creation time: chr "1.28.3"
##   ..$ RSQLite version at creation time        : chr "1.1-2"
##   ..$ DBSCHEMAVERSION                         : chr "1.1"

The colData:

colData(se)
## DataFrame with 8 rows and 0 columns

The colData slot, so far empty, should contain all the metadata. Because we used a column of sampleTable to produce the bamfiles vector, we know the columns of se are in the same order as the rows of sampleTable. We can assign the sampleTable as the colData of the summarized experiment, by converting it into a DataFrame and using the assignment function:

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

3.8 Branching point

At this point, we have counted the fragments which overlap the genes in the gene model we specified. This is a branching point where we could use a variety of Bioconductor packages for exploration and differential expression of the count data, including edgeR (M. D. Robinson, McCarthy, and Smyth 2009), limma with the voom method (Law et al. 2014), DSS (H. Wu, Wang, and Wu 2013), EBSeq (Leng et al. 2013) and baySeq (Hardcastle and Kelly 2010). Schurch et al. (2016) compared performance of different statistical methods for RNA-seq using a large number of biological replicates and can help users to decide which tools make sense to use, and how many biological replicates are necessary to obtain a certain sensitivity. We will continue using DESeq2 (M. I. Love, Huber, and Anders 2014) and edgeR (M. D. Robinson, McCarthy, and Smyth 2009) below. The SummarizedExperiment object is all we need to start our analysis. In the following section we will show how to use it to create the data object used by DESeq2.

4 The DESeqDataSet object, sample information and the design formula

Bioconductor software packages often define and use a custom class for storing data that makes sure that all the needed data slots are consistently provided and fulfill the requirements. In addition, Bioconductor has general data classes (such as the SummarizedExperiment) that can be used to move data between packages. Additionally, the core Bioconductor classes provide useful functionality: for example, subsetting or reordering the rows or columns of a SummarizedExperiment automatically subsets or reorders the associated rowRanges and colData, which can help to prevent accidental sample swaps that would otherwise lead to spurious results. With SummarizedExperiment this is all taken care of behind the scenes.

In DESeq2, the custom class is called DESeqDataSet. It is built on top of the SummarizedExperiment class, and it is easy to convert SummarizedExperiment objects into DESeqDataSet objects, which we show below. One of the two main differences is that the assay slot is instead accessed using the counts accessor function, and the DESeqDataSet class enforces that the values in this matrix are non-negative integers.

A second difference is that the DESeqDataSet has an associated design formula. The experimental design is specified at the beginning of the analysis, as it will inform many of the DESeq2 functions how to treat the samples in the analysis (one exception is the size factor estimation, i.e., the adjustment for differing library sizes, which does not depend on the design formula). The design formula tells which columns in the sample information table (colData) specify the experimental design and how these factors should be used in the analysis.

The simplest design formula for differential expression would be ~ condition, where condition is a column in colData(dds) that specifies which of two (or more groups) the samples belong to. For the airway experiment, we will specify ~ cell + dex meaning that we want to test for the effect of dexamethasone (dex) controlling for the effect of different cell line (cell). We can see each of the columns just using the $ directly on the SummarizedExperiment or DESeqDataSet:

se$cell
## [1] N61311  N61311  N052611 N052611 N080611 N080611 N061011 N061011
## Levels: N052611 N061011 N080611 N61311
se$dex
## [1] untrt trt   untrt trt   untrt trt   untrt trt  
## Levels: trt untrt

Note: it is prefered in R that the first level of a factor be the reference level (e.g. control, or untreated samples), so we can relevel the dex factor like so:

library("magrittr")
se$dex %<>% relevel("untrt")
se$dex
## [1] untrt trt   untrt trt   untrt trt   untrt trt  
## Levels: untrt trt

%<>% is the compound assignment pipe-operator from the magrittr package, the above line of code is a concise way of saying

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

For running DESeq2 or edgeR models, you can use R’s formula notation to express any fixed-effects experimental design. Note that DESeq2 and edgeR use the same formula notation as, for instance, the lm function of base R. If the research aim is to determine for which genes the effect of treatment is different across groups, then interaction terms can be included and tested using a design such as ~ group + treatment + group:treatment. See the manual page for ?results for more examples. We will show how to use an interaction term to test for condition-specific changes over time in a time course example below.

In the following sections, we will demonstrate the construction of the DESeqDataSet from two starting points:

For a full example of using the HTSeq Python package for read counting, please see the pasilla vignette. For an example of generating the DESeqDataSet from files produced by htseq-count, please see the DESeq2 vignette.

4.1 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 Himes et al. (2014), 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.

data("airway")
se <- airway

Again, we want to specify that untrt is the reference level for the dex variable:

se$dex %<>% relevel("untrt")
se$dex
## [1] untrt trt   untrt trt   untrt trt   untrt trt  
## Levels: untrt trt

We can quickly check the millions of fragments that uniquely aligned to the genes (the second argument of round tells how many decimal points to keep).

round(colSums(assay(se)) / 1e6, 1)
## SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520 
##       20.6       18.8       25.3       15.2       24.4       30.8       19.1 
## SRR1039521 
##       21.2

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
##              <factor> <factor> <factor> <factor>   <factor> <integer>
## SRR1039508 GSM1275862   N61311    untrt    untrt SRR1039508       126
## SRR1039509 GSM1275863   N61311      trt    untrt SRR1039509       126
## SRR1039512 GSM1275866  N052611    untrt    untrt SRR1039512       126
## SRR1039513 GSM1275867  N052611      trt    untrt SRR1039513        87
## SRR1039516 GSM1275870  N080611    untrt    untrt SRR1039516       120
## SRR1039517 GSM1275871  N080611      trt    untrt SRR1039517       126
## SRR1039520 GSM1275874  N061011    untrt    untrt SRR1039520       101
## SRR1039521 GSM1275875  N061011      trt    untrt SRR1039521        98
##            Experiment    Sample    BioSample
##              <factor>  <factor>     <factor>
## SRR1039508  SRX384345 SRS508568 SAMN02422669
## SRR1039509  SRX384346 SRS508567 SAMN02422675
## SRR1039512  SRX384349 SRS508571 SAMN02422678
## SRR1039513  SRX384350 SRS508572 SAMN02422670
## SRR1039516  SRX384353 SRS508575 SAMN02422682
## SRR1039517  SRX384354 SRS508576 SAMN02422673
## SRR1039520  SRX384357 SRS508579 SAMN02422683
## SRR1039521  SRX384358 SRS508580 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 earlier by specifying the BAM files using a column of the sample table.

Once we have our fully annotated SummarizedExperiment object, we can construct a DESeqDataSet object from it that will then form the starting point of the analysis. We add an appropriate design for the analysis:

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

4.2 Starting from count matrices

In this section, we will show how to build an DESeqDataSet supposing we only have a count matrix and a table of sample information.

Note: if you have prepared a SummarizedExperiment you should skip this section. While the previous section would be used to construct a DESeqDataSet from a SummarizedExperiment, here we first extract the individual object (count matrix and sample info) from the SummarizedExperiment in order to build it back up into a new object – only for demonstration purposes. In practice, the count matrix would either be read in from a file or perhaps generated by an R function like featureCounts from the Rsubread package (Liao, Smyth, and Shi 2014).

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

countdata <- assay(se)
head(countdata, 3)
##                 SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516
## ENSG00000000003        679        448        873        408       1138
## ENSG00000000005          0          0          0          0          0
## ENSG00000000419        467        515        621        365        587
##                 SRR1039517 SRR1039520 SRR1039521
## ENSG00000000003       1047        770        572
## ENSG00000000005          0          0          0
## ENSG00000000419        799        417        508

In this count matrix, each row represents an Ensembl gene, each column a sequenced RNA library, and the values give the raw numbers of fragments that were uniquely assigned to the respective gene in each library. We also have information on each of the samples (the columns of the count matrix). If you’ve counted reads with some other software, it is very important to check that the columns of the count matrix correspond to the rows of the sample information table.

coldata <- colData(se)

We now have all the ingredients to prepare our data object in a form that is suitable for analysis, namely:

  • countdata: a table with the fragment counts
  • coldata: a table with information about the samples

To now construct the DESeqDataSet object from the matrix of counts and the sample information table, we use:

ddsMat <- DESeqDataSetFromMatrix(countData = countdata,
                                 colData = coldata,
                                 design = ~cell + dex)

We will continue with the object generated from the SummarizedExperiment section.

4.3 Creating a DGEList for use with edgeR

The edgeR package uses another type of data container, namely a DGEList object. It is just as easy to create a DGEList object using the count matrix and information about samples. We can additionally add information about the genes:

library("edgeR")
genetable <- data.frame(gene.id = rownames(se))
dge <- DGEList(counts = countdata, 
               samples = coldata, 
               genes = genetable)
names(dge)
## [1] "counts"  "samples" "genes"

Just like the SummarizedExperiment and the DESeqDataSet the DGEList contains all the information we need: the count matrix, information about the samples (columns of the count matrix), and information about the genes (rows of the count matrix).

5 Exploratory analysis and visualization

There are two separate paths in this workflow; the one we will see first involves transformations of the counts in order to visually explore sample relationships. In the second part, we will go back to the original raw counts for statistical testing. This is critical because the statistical testing methods rely on original count data (not scaled or transformed) for calculating the precision of measurements.

5.1 Pre-filtering the dataset

Our count matrix with our DESeqDataSet contains many rows with only zeros, and additionally many rows with only a few fragments total. In order to reduce the size of the object, and to increase the speed of our functions, we can remove the rows that have no or nearly no information about the amount of gene expression. Here we apply the most minimal filtering rule: removing rows of the DESeqDataSet that have no counts, or only a single count across all samples. Additional weighting/filtering to improve power is applied at a later step in the workflow.

nrow(dds)
## [1] 64102
dds <- dds[ rowSums(counts(dds)) > 1, ]
nrow(dds)
## [1] 29391

Here, we also filter the DGEList for edgeR in the same way. Note, however, that edgeR does not apply filtering internally, and for this reason, typically a stronger prefiltering criterion is used at this stage for edgeR.

dge <- dge[rowSums(dge$counts) > 1, ]
all(rownames(dge) == rownames(dds))
## [1] TRUE

5.2 The rlog and variance stabilizing transformations

Many common statistical methods for exploratory analysis of multidimensional data, for example clustering and principal components analysis (PCA), work best for data that generally has the same range of variance at different ranges of the mean values. When the expected amount of variance is approximately the same across different mean values, the data is said to be homoskedastic. For RNA-seq counts, however, the expected variance grows with the mean. For example, if one performs PCA directly on a matrix of counts or normalized counts (e.g. correcting for differences in sequencing depth), the resulting plot typically depends mostly on the genes with highest counts because they show the largest absolute differences between samples. A simple and often used strategy to avoid this is to take the logarithm of the normalized count values plus a pseudocount of 1; however, depending on the choice of pseudocount, now the genes with the very lowest counts will contribute a great deal of noise to the resulting plot, because taking the logarithm of small counts actually inflates their variance. We can quickly show this property of counts with some simulated data (here, Poisson counts with a range of lambda from 0.1 to 100). We plot the standard deviation of each row (genes) against the mean:

lambda <- 10^seq(from = -1, to = 2, length = 1000)
cts <- matrix(rpois(1000*100, lambda), ncol = 100)
library("vsn")
meanSdPlot(cts, ranks = FALSE)

And for logarithm-transformed counts after adding a pseudocount of 1:

log.cts.one <- log2(cts + 1)
meanSdPlot(log.cts.one, ranks = FALSE)

The logarithm with a small pseudocount amplifies differences when the values are close to 0. The low count genes with low signal-to-noise ratio will overly contribute to sample-sample distances and PCA plots.

As a solution, DESeq2 offers two transformations for count data that stabilize the variance across the mean: the regularized-logarithm transformation or rlog (M. I. Love, Huber, and Anders 2014), and the variance stabilizing transformation (VST) for negative binomial data with a dispersion-mean trend (Anders and Huber 2010), implemented in the vst function.

For genes with high counts, the rlog and VST will give similar result to the ordinary log2 transformation of normalized counts. For genes with lower counts, however, the values are shrunken towards the genes’ averages across all samples. The rlog-transformed or VST data then becomes approximately homoskedastic, and can be used directly for computing distances between samples, making PCA plots, or as input to downstream methods which perform best with homoskedastic data.

Which transformation to choose? The rlog tends to work well on small datasets (n < 30), sometimes outperforming the VST when there is a large range of sequencing depth across samples (an order of magnitude difference). The VST is much faster to compute and is less sensitive to high count outliers than the rlog. We therefore recommend the VST for large datasets (hundreds of samples). You can perform both transformations and compare the meanSdPlot or PCA plots generated, as described below.

Note that the two transformations offered by DESeq2 are provided for applications other than differential testing. For differential testing we recommend the DESeq function applied to raw counts, as described later in this workflow, which also takes into account the dependence of the variance of counts on the mean value during the dispersion estimation step.

The function rlog returns an object based on the SummarizedExperiment class that contains the rlog-transformed values in its assay slot.

rld <- rlog(dds, blind = FALSE)
head(assay(rld), 3)
##                 SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516
## ENSG00000000003   9.385683   9.052608   9.516875   9.285338   9.839085
## ENSG00000000419   8.869616   9.138271   9.036116   9.075295   8.972126
## ENSG00000000457   7.961369   7.881385   7.824079   7.921490   7.751699
##                 SRR1039517 SRR1039520 SRR1039521
## ENSG00000000003   9.530311   9.663255   9.277699
## ENSG00000000419   9.131824   8.861534   9.060905
## ENSG00000000457   7.886441   7.957121   7.912123
meanSdPlot(assay(rld), ranks = FALSE)

The function vst returns a similar object:

vsd <- vst(dds, blind = FALSE)
head(assay(vsd), 3)
##                 SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516
## ENSG00000000003   9.742340   9.430742   9.867872   9.646127  10.183344
## ENSG00000000419   9.334009   9.582000   9.486456   9.523397   9.427605
## ENSG00000000457   8.765748   8.698941   8.651978   8.732909   8.593308
##                 SRR1039517 SRR1039520 SRR1039521
## ENSG00000000003   9.880660  10.010591   9.640065
## ENSG00000000419   9.575154   9.326341   9.509553
## ENSG00000000457   8.703164   8.762420   8.724586
meanSdPlot(assay(vsd), ranks = FALSE)

In the above function calls, we specified blind = FALSE, which means that differences between cell lines and treatment (the variables in the design) will not contribute to the expected variance-mean trend of the experiment. The experimental design is not used directly in the transformation, only in estimating the global amount of variability in the counts. For a fully unsupervised transformation, one can set blind = TRUE (which is the default).

To show the effect of the transformation, in the figure below we plot the first sample against the second, first simply using the log2 function (after adding 1, to avoid taking the log of zero), and then using the rlog- and VST-transformed values. For the log2 approach, we need to first estimate size factors to account for sequencing depth, and then specify normalized = TRUE. Sequencing depth correction is done automatically for the rlog and the vst.

library("dplyr")
library("ggplot2")

dds <- estimateSizeFactors(dds)

df <- bind_rows(
  as_data_frame(log2(counts(dds, normalized = TRUE)[, 1:2]+1)) %>%
         mutate(transformation = "log2(x + 1)"),
  as_data_frame(assay(rld)[, 1:2]) %>% mutate(transformation = "rlog"),
  as_data_frame(assay(vsd)[, 1:2]) %>% mutate(transformation = "vst"))
  
colnames(df)[1:2] <- c("x", "y")  

ggplot(df, aes(x = x, y = y)) + geom_hex(bins = 80) +
  coord_fixed() + facet_grid(. ~ transformation)  

Scatterplot of transformed counts from two samples. Shown are scatterplots using the log2 transform of normalized counts (left), using the rlog (middle), and using the VST (right). While the rlog is on roughly the same scale as the log2 counts, the VST has a upward shift for the smaller values. It is the differences between samples (deviation from y = x in these scatterplots) which will contribute to the distance calculations and the PCA plot.

We can see how genes with low counts (bottom left-hand corner) seem to be excessively variable on the ordinary logarithmic scale, while the rlog transform and VST compress differences for the low count genes for which the data provide little information about differential expression.

5.3 Sample distances

A useful first step in an RNA-seq analysis is often to assess overall similarity between samples: Which samples are similar to each other, which are different? Does this fit to the expectation from the experiment’s design?

We use the R function dist to calculate the Euclidean distance between samples. To ensure we have a roughly equal contribution from all genes, we use it on the rlog-transformed data. We need to transpose the matrix of values using t, because the dist function expects the different samples to be rows of its argument, and different dimensions (here, genes) to be columns.

sampleDists <- dist(t(assay(rld)))
sampleDists
##            SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517
## SRR1039509   45.69859                                                       
## SRR1039512   39.25239   54.90828                                            
## SRR1039513   62.63201   44.52740   48.72579                                 
## SRR1039516   44.50557   59.06364   43.57856   63.74275                      
## SRR1039517   64.49410   51.44882   59.22962   49.87992   47.48200           
## SRR1039520   39.57693   57.46259   36.74434   58.49014   46.40786   63.59942
## SRR1039521   63.36124   45.05732   57.87616   36.49484   65.54600   52.31695
##            SRR1039520
## SRR1039509           
## SRR1039512           
## SRR1039513           
## SRR1039516           
## SRR1039517           
## SRR1039520           
## SRR1039521   50.13430

We visualize the distances in a heatmap in a figure below, using the function pheatmap from the pheatmap package.

library("pheatmap")
library("RColorBrewer")

In order to plot the sample distance matrix with the rows/columns arranged by the distances in our distance matrix, we manually provide sampleDists to the clustering_distance argument of the pheatmap function. Otherwise the pheatmap function would assume that the matrix contains the data values themselves, and would calculate distances between the rows/columns of the distance matrix, which is not desired. We also manually specify a blue color palette using the colorRampPalette function from the RColorBrewer package.

sampleDistMatrix <- as.matrix(sampleDists)
rownames(sampleDistMatrix) <- paste(rld$dex, rld$cell, sep = " - ")
colnames(sampleDistMatrix) <- NULL
colors <- colorRampPalette(rev(brewer.pal(9, "Blues")))(255)
pheatmap(sampleDistMatrix,
         clustering_distance_rows = sampleDists,
         clustering_distance_cols = sampleDists,
         col = colors)

Heatmap of sample-to-sample distances using the rlog-transformed values.

Note that we have changed the row names of the distance matrix to contain treatment type and patient number instead of sample ID, so that we have all this information in view when looking at the heatmap.

Another option for calculating sample distances is to use the Poisson Distance (Witten 2011), implemented in the PoiClaClu package. This measure of dissimilarity between counts also takes the inherent variance structure of counts into consideration when calculating the distances between samples. The PoissonDistance function takes the original count matrix (not normalized) with samples as rows instead of columns, so we need to transpose the counts in dds.

library("PoiClaClu")
poisd <- PoissonDistance(t(counts(dds)))

We plot the heatmap in a Figure below.

samplePoisDistMatrix <- as.matrix(poisd$dd)
rownames(samplePoisDistMatrix) <- paste(rld$dex, rld$cell, sep = " - ")
colnames(samplePoisDistMatrix) <- NULL
pheatmap(samplePoisDistMatrix,
         clustering_distance_rows = poisd$dd,
         clustering_distance_cols = poisd$dd,
         col = colors)

Heatmap of sample-to-sample distances using the Poisson Distance.

5.4 PCA plot

Another way to visualize sample-to-sample distances is a principal components analysis (PCA). In this ordination method, the data points (here, the samples) are projected onto the 2D plane such that they spread out in the two directions that explain most of the differences (figure below). The x-axis is the direction that separates the data points the most. The values of the samples in this direction are written PC1. The y-axis is a direction (it must be orthogonal to the first direction) that separates the data the second most. The values of the samples in this direction are written PC2. The percent of the total variance that is contained in the direction is printed in the axis label. Note that these percentages do not add to 100%, because there are more dimensions that contain the remaining variance (although each of these remaining dimensions will explain less than the two that we see).

plotPCA(rld, intgroup = c("dex", "cell"))

PCA plot using the rlog-transformed values. Each unique combination of treatment and cell line is given its own color.

Here, we have used the function plotPCA that comes with DESeq2. The two terms specified by intgroup are the interesting groups for labeling the samples; they tell the function to use them to choose colors. We can also build the PCA plot from scratch using the ggplot2 package (Wickham 2009). This is done by asking the plotPCA function to return the data used for plotting rather than building the plot. See the ggplot2 documentation for more details on using ggplot.

pcaData <- plotPCA(rld, intgroup = c("dex", "cell"), returnData = TRUE)
pcaData
##                  PC1        PC2           group   dex    cell       name
## SRR1039508 -17.81773  -4.020836  untrt : N61311 untrt  N61311 SRR1039508
## SRR1039509   8.38790  -1.490805    trt : N61311   trt  N61311 SRR1039509
## SRR1039512 -10.22735  -5.004069 untrt : N052611 untrt N052611 SRR1039512
## SRR1039513  17.53277  -3.909890   trt : N052611   trt N052611 SRR1039513
## SRR1039516 -14.67169  15.873239 untrt : N080611 untrt N080611 SRR1039516
## SRR1039517  10.98782  20.598625   trt : N080611   trt N080611 SRR1039517
## SRR1039520 -12.06035 -11.985876 untrt : N061011 untrt N061011 SRR1039520
## SRR1039521  17.86863 -10.060389   trt : N061011   trt N061011 SRR1039521
percentVar <- round(100 * attr(pcaData, "percentVar"))

We can then use these data to build up a second plot in a figure below, specifying that the color of the points should reflect dexamethasone treatment and the shape should reflect the cell line.

ggplot(pcaData, aes(x = PC1, y = PC2, color = dex, shape = cell)) +
  geom_point(size = 3) +
  xlab(paste0("PC1: ", percentVar[1], "% variance")) +
  ylab(paste0("PC2: ", percentVar[2], "% variance")) +
  coord_fixed()

PCA plot using the rlog-transformed values with custom ggplot2 code. Here we specify cell line (plotting symbol) and dexamethasone treatment (color).

From the PCA plot, we see that the differences between cells (the different plotting shapes) are considerable, though not stronger than the differences due to treatment with dexamethasone (red vs blue color). This shows why it will be important to account for this in differential testing by using a paired design (“paired”, because each dex treated sample is paired with one untreated sample from the same cell line). We are already set up for this design by assigning the formula ~ cell + dex earlier.

5.5 MDS plot

Another plot, very similar to the PCA plot, can be made using the multidimensional scaling (MDS) function in base R. This is useful when we don’t have a matrix of data, but only a matrix of distances. Here we compute the MDS for the distances calculated from the rlog transformed counts and plot these in a figure below.

mds <- as.data.frame(colData(rld))  %>%
         cbind(cmdscale(sampleDistMatrix))
ggplot(mds, aes(x = `1`, y = `2`, color = dex, shape = cell)) +
  geom_point(size = 3) + coord_fixed()

MDS plot using rlog-transformed values.

In a figure below we show the same plot for the PoissonDistance:

mdsPois <- as.data.frame(colData(dds)) %>%
   cbind(cmdscale(samplePoisDistMatrix))
ggplot(mdsPois, aes(x = `1`, y = `2`, color = dex, shape = cell)) +
  geom_point(size = 3) + coord_fixed()

MDS plot using the Poisson Distance.

6 Differential expression analysis

6.1 Performing differential expression testing with DESeq2

As we have already specified an experimental design when we created the DESeqDataSet, we can run the differential expression pipeline on the raw counts with a single call to the function DESeq:

dds <- DESeq(dds)

This function will print out a message for the various steps it performs. These are described in more detail in the manual page for DESeq, which can be accessed by typing ?DESeq. Briefly these are: the estimation of size factors (controlling for differences in the sequencing depth of the samples), the estimation of dispersion values for each gene, and fitting a generalized linear model.

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

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. The comparison is printed at the top of the output: dex trt vs untrt.

res <- results(dds)
res
## log2 fold change (MLE): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 29391 rows and 6 columns
##                    baseMean log2FoldChange     lfcSE        stat
##                   <numeric>      <numeric> <numeric>   <numeric>
## ENSG00000000003 708.6021697    -0.38125388 0.1006560  -3.7876928
## ENSG00000000419 520.2979006     0.20681271 0.1122218   1.8428925
## ENSG00000000457 237.1630368     0.03792050 0.1434532   0.2643405
## ENSG00000000460  57.9326331    -0.08816322 0.2871677  -0.3070095
## ENSG00000000938   0.3180984    -1.37823397 3.4998753  -0.3937952
## ...                     ...            ...       ...         ...
## ENSG00000273485   1.2864477     -0.1271363 1.6005571 -0.07943253
## ENSG00000273486  15.4525365     -0.1509944 0.4865490 -0.31033758
## ENSG00000273487   8.1632350      1.0464169 0.6990336  1.49694792
## ENSG00000273488   8.5844790      0.1078190 0.6381646  0.16895172
## ENSG00000273489   0.2758994      1.4837258 3.5139452  0.42223933
##                       pvalue        padj
##                    <numeric>   <numeric>
## ENSG00000000003 0.0001520527 0.001281516
## ENSG00000000419 0.0653447023 0.196230403
## ENSG00000000457 0.7915175170 0.911219224
## ENSG00000000460 0.7588361136 0.894672849
## ENSG00000000938 0.6937322727          NA
## ...                      ...         ...
## ENSG00000273485    0.9366886          NA
## ENSG00000273486    0.7563043   0.8936787
## ENSG00000273487    0.1344068   0.3291460
## ENSG00000273488    0.8658346   0.9451750
## ENSG00000273489    0.6728503          NA

We could have equivalently produced this results table with the following more specific command. Because dex is the last variable in the design, we could optionally leave off the contrast argument to extract the comparison of the two levels of dex.

res <- results(dds, contrast = c("dex", "trt", "untrt"))

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 (MLE): 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, divided by the size factors, taken over all samples in the DESeqDataSet. 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. We will find out below 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 zero 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 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, that will be covered in later sections.

summary(res)
## 
## out of 29391 with nonzero total read count
## adjusted p-value < 0.1
## LFC > 0 (up)     : 2607, 8.9% 
## LFC < 0 (down)   : 2218, 7.5% 
## outliers [1]     : 0, 0% 
## low counts [2]   : 11397, 39% 
## (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

If we lower the false discovery rate threshold, we should also inform the results() function about it, so that the function can use this threshold for the optimal independent filtering that it performs:

res.05 <- results(dds, alpha = 0.05)
table(res.05$padj < 0.05)
## 
## FALSE  TRUE 
## 12831  4024

If we want to raise the log2 fold change threshold, so that we test for genes that show more substantial changes due to treatment, we simply supply a value on the log2 scale. For example, by specifying lfcThreshold = 1, we test for genes that show significant effects of treatment on gene counts more than doubling or less than halving, because \(2^1=2\).

resLFC1 <- results(dds, lfcThreshold = 1)
table(resLFC1$padj < 0.1)
## 
## FALSE  TRUE 
## 20034   240

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 no 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 DESeq2 vignette.

If you use the results from an R analysis package in published research, you can find the proper citation for the software by typing citation("pkgName"), where you would substitute the name of the package for pkgName. Citing methods papers helps to support and reward the individuals who put time into open source software for genomic data analysis.

6.3 Other comparisons

In general, the results for a comparison of any two levels of a variable can be extracted using the contrast argument to results. The user should specify three values: the name of the variable, the name of the level for the numerator, and the name of the level for the denominator. Here we extract results for the log2 of the fold change of one cell line over another:

results(dds, contrast = c("cell", "N061011", "N61311"))
## log2 fold change (MLE): cell N061011 vs N61311 
## Wald test p-value: cell N061011 vs N61311 
## DataFrame with 29391 rows and 6 columns
##                    baseMean log2FoldChange     lfcSE         stat     pvalue
##                   <numeric>      <numeric> <numeric>    <numeric>  <numeric>
## ENSG00000000003 708.6021697     0.30632635 0.1435330  2.134187335 0.03282745
## ENSG00000000419 520.2979006    -0.05404672 0.1597205 -0.338383153 0.73507447
## ENSG00000000457 237.1630368     0.01630854 0.2030380  0.080322608 0.93598068
## ENSG00000000460  57.9326331     0.27912688 0.4007042  0.696590931 0.48605883
## ENSG00000000938   0.3180984     0.03849337 4.9975798  0.007702403 0.99385443
## ...                     ...            ...       ...          ...        ...
## ENSG00000273485   1.2864477    -2.32142337 2.3530918 -0.986541762  0.3238673
## ENSG00000273486  15.4525365    -0.07271569 0.7119582 -0.102134773  0.9186497
## ENSG00000273487   8.1632350    -0.02693764 1.0051679 -0.026799147  0.9786199
## ENSG00000273488   8.5844790     1.09584871 0.8990390  1.218911152  0.2228779
## ENSG00000273489   0.2758994     0.03849337 4.9975798  0.007702403  0.9938544
##                      padj
##                 <numeric>
## ENSG00000000003 0.2206411
## ENSG00000000419 0.9395460
## ENSG00000000457 0.9904710
## ENSG00000000460 0.8550869
## ENSG00000000938        NA
## ...                   ...
## ENSG00000273485        NA
## ENSG00000273486 0.9865440
## ENSG00000273487 0.9958619
## ENSG00000273488 0.6482708
## ENSG00000273489        NA

There are additional ways to build results tables for certain comparisons after running DESeq once. If results for an interaction term are desired, the name argument of results should be used. Please see the help page for the results function for details on the additional ways to build results tables. In particular, the Examples section of the help page for results gives some pertinent examples.

6.4 Performing differential expression testing with edgeR

Next we will show how to perform differential expression analysis with edgeR. Recall that we have a DGEList dge, containing three objects:

names(dge)
## [1] "counts"  "samples" "genes"

We first define a design matrix, using the same formula syntax as for DESeq2 above.

design <- model.matrix(~ cell + dex, dge$samples)

Then, we calculate normalization factors and estimate the dispersion for each gene. Note that we need to specify the design in the dispersion calculation.

dge <- calcNormFactors(dge)
dge <- estimateDisp(dge, design)

Finally, we fit the generalized linear model and perform the test. In the glmLRT function, we indicate which coefficient (which column in the design matrix) that we would like to test for. It is possible to test more general contrasts as well, and the user guide contains many examples on how to do this. The topTags function extracts the top-ranked genes. You can indicate the adjusted p-value cutoff, and/or the number of genes to keep.

fit <- glmFit(dge, design)
lrt <- glmLRT(fit, coef = ncol(design))
tt <- topTags(lrt, n = nrow(dge), p.value = 0.1)
tt10 <- topTags(lrt) # just the top 10 by default
tt10
## Coefficient:  dextrt 
##               gene.id     logFC   logCPM        LR        PValue
## 3751  ENSG00000109906  7.141823 4.143740 1293.2297 3.345882e-283
## 9658  ENSG00000152583  4.552568 5.528678  901.7880 4.009940e-198
## 11864 ENSG00000165995  3.274517 4.504348  749.6856 4.696274e-165
## 11292 ENSG00000163884  4.458891 4.685058  703.8823 4.280321e-155
## 17245 ENSG00000189221  3.335100 6.766777  647.7861 6.770158e-143
## 6020  ENSG00000127954  5.197020 3.652850  638.6802 6.471125e-141
## 10923 ENSG00000162692 -3.677380 4.599336  634.3727 5.595117e-140
## 2529  ENSG00000101347  3.751932 9.206504  631.1283 2.840687e-139
## 13316 ENSG00000171819  5.668260 3.500283  615.0520 8.910996e-136
## 5054  ENSG00000120129  2.932273 7.308257  603.5211 2.870368e-133
##                 FDR
## 3751  9.833881e-279
## 9658  5.892808e-194
## 11864 4.600940e-161
## 11292 3.145073e-151
## 17245 3.979634e-139
## 6020  3.169881e-137
## 10923 2.349230e-136
## 2529  1.043633e-135
## 13316 2.910034e-132
## 5054  8.436298e-130

We can compare to see how the results from the two software overlap:

tt.all <- topTags(lrt, n = nrow(dge), sort.by = "none")
table(DESeq2 = res$padj < 0.1, edgeR = tt.all$table$FDR < 0.1)
##        edgeR
## DESeq2  FALSE  TRUE
##   FALSE 11826  1343
##   TRUE    148  4677

We can also compare the two lists by the ranks:

common <- !is.na(res$padj)
plot(rank(res$padj[common]), 
     rank(tt.all$table$FDR[common]), cex = .1,
     xlab = "DESeq2", ylab = "edgeR")

Also with edgeR we can test for significance relative to a fold-change threshold, using the function glmTreat instead of glmLRT. Below we set the log fold-change threshold to 1 (i.e., fold change threshold equal to 2), as for DESeq2 above. We also compare the results to those obtained with DESeq2.

treatres <- glmTreat(fit, coef = ncol(design), lfc = 1)
tt.treat <- topTags(treatres, n = nrow(dge), sort.by = "none")
table(DESeq2 = resLFC1$padj < 0.1, edgeR = tt.treat$table$FDR < 0.1)
##        edgeR
## DESeq2  FALSE  TRUE
##   FALSE 19759   275
##   TRUE      0   240

6.5 Multiple testing

In high-throughput biology, we are careful to not use the p values directly as evidence against the null, but to correct for multiple testing. What would happen if we were to simply threshold the p values at a low value, say 0.05? With DESeq2, there are 5676 genes with a p value below 0.05 among the 29391 genes for which the test succeeded in reporting a p value:

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

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 the p value, we expect up to 5% of the genes to have a p value below 0.05. This amounts to 1470 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 1470 / 5676= 26% false positives.

DESeq2 uses the Benjamini-Hochberg (BH) adjustment (Benjamini and Hochberg 1995) as implemented in the base R p.adjust function; in brief, this method calculates for each gene an adjusted p value that answers the following question: if one called significant all genes with an adjusted p value less than or equal to this gene’s adjusted 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.

The FDR is a useful statistic for many high-throughput experiments, as we are often interested in reporting or focusing on a set of interesting genes, and we would like to put an upper bound on the percent of false positives in this set.

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] 4825

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 (MLE): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 6 columns
##                  baseMean log2FoldChange     lfcSE      stat       pvalue
##                 <numeric>      <numeric> <numeric> <numeric>    <numeric>
## ENSG00000128285  6.624741      -5.325912 1.2578863 -4.234017 2.295537e-05
## ENSG00000267339 26.233573      -4.611553 0.6731316 -6.850894 7.338997e-12
## ENSG00000019186 14.087605      -4.325920 0.8578247 -5.042895 4.585398e-07
## ENSG00000183454  5.804171      -4.264087 1.1669498 -3.654045 2.581412e-04
## ENSG00000146006 46.807597      -4.211875 0.5288797 -7.963767 1.668799e-15
## ENSG00000141469 53.436528      -4.124784 1.1297977 -3.650905 2.613174e-04
##                         padj
##                    <numeric>
## ENSG00000128285 2.378002e-04
## ENSG00000267339 2.053778e-10
## ENSG00000019186 6.611351e-06
## ENSG00000183454 2.048057e-03
## ENSG00000146006 7.166675e-14
## ENSG00000141469 2.069606e-03

…and with the strongest up-regulation:

head(resSig[ order(resSig$log2FoldChange, decreasing = TRUE), ])
## log2 fold change (MLE): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 6 columns
##                   baseMean log2FoldChange     lfcSE      stat       pvalue
##                  <numeric>      <numeric> <numeric> <numeric>    <numeric>
## ENSG00000179593  67.243048       9.505972 1.0545111  9.014578 1.976299e-19
## ENSG00000109906 385.071029       7.352628 0.5363902 13.707610 9.141988e-43
## ENSG00000250978  56.318194       6.327393 0.6778153  9.334981 1.010098e-20
## ENSG00000132518   5.654654       5.885113 1.3241367  4.444491 8.810031e-06
## ENSG00000127954 286.384119       5.207160 0.4930828 10.560419 4.546302e-26
## ENSG00000249364   8.839061       5.098168 1.1596852  4.396166 1.101798e-05
##                         padj
##                    <numeric>
## ENSG00000179593 1.252166e-17
## ENSG00000109906 2.253437e-40
## ENSG00000250978 7.212582e-19
## ENSG00000132518 1.000175e-04
## ENSG00000127954 5.049763e-24
## ENSG00000249364 1.223812e-04

7 Plotting results

7.1 Counts plot

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

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

Normalized counts for a single gene over treatment group.

We can also make custom plots using the ggplot function from the ggplot2 package (figures below).

library("ggbeeswarm")
geneCounts <- plotCounts(dds, gene = topGene, intgroup = c("dex","cell"),
                         returnData = TRUE)
ggplot(geneCounts, aes(x = dex, y = count, color = cell)) +
  scale_y_log10() +  geom_beeswarm(cex = 3)

ggplot(geneCounts, aes(x = dex, y = count, color = cell, group = cell)) +
  scale_y_log10() + geom_point(size = 3) + geom_line()

Normalized counts with lines connecting cell lines. Note that the DESeq test actually takes into account the cell line effect, so this figure more closely depicts the difference being tested.

7.2 MA plot with DESeq2

An MA plot (R. Dudoit et al. 2002) provides a useful overview for the distribution of the estimated coefficients in the model, e.g. the comparisons of interest, across all genes. On the y-axis, the “M” stands for “minus” – subtraction of log values is equivalent to the log of the ratio – and on the x-axis, the “A” stands for “average”. You may hear this plot also referred to as a mean-difference plot, or a Bland-Altman plot.

Before making the MA plot, we use the lfcShrink function to shrink the log2 fold changes for the comparison of dex treated vs untreated samples:

res <- lfcShrink(dds, contrast = c("dex","trt","untrt"), res = res)
DESeq2::plotMA(res, ylim = c(-5, 5))

An MA plot of changes induced by treatment. 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. 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 uses a Bayesian procedure to moderate (or “shrink”) log2 fold changes from genes with very low counts and highly variable counts, as can be seen by the narrowing of the vertical spread of points on the left side of the MA plot. As shown above, the lfcShrink function performs this operation. For a detailed explanation of the rationale of moderated fold changes, please see the DESeq2 paper (M. I. Love, Huber, and Anders 2014).

If we had not used statistical moderation to shrink the noisy log2 fold changes, we would have instead seen the following plot:

res.noshr <- results(dds)
DESeq2::plotMA(res.noshr, ylim = c(-5, 5))

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.

DESeq2::plotMA(res, ylim = c(-5,5))
topGene <- rownames(res)[which.min(res$padj)]
with(res[topGene, ], {
  points(baseMean, log2FoldChange, col = "dodgerblue", cex = 2, lwd = 2)
  text(baseMean, log2FoldChange, topGene, pos = 2, col = "dodgerblue")
})

Another useful diagnostic plot is the histogram of the p values (figure below). This plot is best formed by excluding genes with very small counts, which otherwise generate spikes in the histogram.

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

Histogram of p values for genes with mean normalized count larger than 1.

7.3 MA / Smear plot with edgeR

In edgeR, the MA plot is obtained via the plotSmear function.

plotSmear(lrt, de.tags = tt$table$gene.id)

7.4 Gene clustering

In the sample distance heatmap made previously, the dendrogram at the side shows us a hierarchical clustering of the samples. Such a clustering can also be performed for the genes. Since the clustering is only relevant for genes that actually carry a signal, one usually would only cluster a subset of the most highly variable genes. Here, for demonstration, let us select the 20 genes with the highest variance across samples. We will work with the rlog transformed counts:

library("genefilter")
topVarGenes <- head(order(rowVars(assay(rld)), decreasing = TRUE), 20)

The heatmap becomes more interesting if we do not look at absolute expression strength but rather at the amount by which each gene deviates in a specific sample from the gene’s average across all samples. Hence, we center each genes’ values across samples, and plot a heatmap (figure below). We provide a data.frame that instructs the pheatmap function how to label the columns.

mat  <- assay(rld)[ topVarGenes, ]
mat  <- mat - rowMeans(mat)
anno <- as.data.frame(colData(rld)[, c("cell","dex")])
pheatmap(mat, annotation_col = anno)

Heatmap of relative rlog-transformed values across samples. Treatment status and cell line information are shown with colored bars at the top of the heatmap. Blocks of genes that covary across patients. Note that a set of genes at the top of the heatmap are separating the N061011 cell line from the others. In the center of the heatmap, we see a set of genes for which the dexamethasone treated samples have higher gene expression.

7.5 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 such high Poisson noise that any biological effect is drowned in the uncertainties from the sampling at a low rate. We can also show this by examining the ratio of small p values (say, less than 0.05) for genes binned by mean normalized count. We will use the results table subjected to the threshold to show what this looks like in a case when there are few tests with small p value.

In the following code chunk, we create bins using the quantile function, bin the genes by base mean using cut, rename the levels of the bins using the middle point, calculate the ratio of p values less than 0.05 for each bin, and finally plot these ratios (figure below).

qs <- c(0, quantile(resLFC1$baseMean[resLFC1$baseMean > 0], 0:6/6))
bins <- cut(resLFC1$baseMean, qs)
levels(bins) <- paste0("~", round(signif((qs[-1] + qs[-length(qs)])/2, 2)))
fractionSig <- tapply(resLFC1$pvalue, bins, function(p)
                          mean(p < .05, na.rm = TRUE))
barplot(fractionSig, xlab = "mean normalized count",
                     ylab = "fraction of small p values")

The ratio of small p values for genes binned by mean normalized count. The p values are from a test of log2 fold change greater than 1 or less than -1. This plot demonstrates that genes with very low mean count have little or no power, and are best excluded from testing.

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 low count genes from the input to the FDR procedure, we can find more genes to be significant among those that we keep, and so improved the power of our test. This approach is known as independent filtering.

The DESeq2 software automatically performs independent filtering that maximizes the number of genes with adjusted p value less than a critical value (by default, alpha is set to 0.1). This automatic independent filtering is performed by, and can be controlled by, the results function.

The term independent highlights an important caveat. Such filtering is permissible only if the statistic that we filter on (here the mean of normalized counts across all samples) is independent of the actual test statistic (the p value) under the null hypothesis. Otherwise, the filtering would invalidate the test and consequently the assumptions of the BH procedure. The independent filtering software used inside DESeq2 comes from the genefilter package, that contains a reference to a paper describing the statistical foundation for independent filtering (Bourgon, Gentleman, and Huber 2010).

8 Annotating and exporting results

Our result table so far only contains the Ensembl gene IDs, but alternative gene names may be more informative for interpretation. 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("AnnotationDbi")
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" 
##  [5] "ENSEMBLTRANS" "ENTREZID"     "ENZYME"       "EVIDENCE"    
##  [9] "EVIDENCEALL"  "GENENAME"     "GO"           "GOALL"       
## [13] "IPI"          "MAP"          "OMIM"         "ONTOLOGY"    
## [17] "ONTOLOGYALL"  "PATH"         "PFAM"         "PMID"        
## [21] "PROSITE"      "REFSEQ"       "SYMBOL"       "UCSCKG"      
## [25] "UNIGENE"      "UNIPROT"

We can use the mapIds function to add individual columns to our results table. We provide the row names of our results table as a key, and specify that keytype = ENSEMBL. The column argument tells the mapIds function which information we want, and the multiVals argument tells the function what to do if there are multiple possible values for a single input value. Here we ask to just give us back the first one that occurs in the database. To add the gene symbol and Entrez ID, we call mapIds twice.

res$symbol <- mapIds(org.Hs.eg.db,
                     keys = row.names(res),
                     column = "SYMBOL",
                     keytype = "ENSEMBL",
                     multiVals = "first")
res$entrez <- mapIds(org.Hs.eg.db,
                     keys = row.names(res),
                     column = "ENTREZID",
                     keytype = "ENSEMBL",
                     multiVals = "first")

Now the results have the desired external gene IDs:

resOrdered <- res[order(res$padj),]
head(resOrdered)
## log2 fold change (MAP): dex trt vs untrt 
## Wald test p-value: dex trt vs untrt 
## DataFrame with 6 rows and 7 columns
##                   baseMean log2FoldChange      stat        pvalue
##                  <numeric>      <numeric> <numeric>     <numeric>
## ENSG00000152583   997.4398       4.313968  24.85607 2.223017e-136
## ENSG00000165995   495.0929       3.186818  24.71198 7.951622e-135
## ENSG00000120129  3409.0294       2.871488  24.27480 3.619142e-130
## ENSG00000101347 12703.3871       3.618752  24.23541 9.423844e-130
## ENSG00000189221  2341.7673       3.230386  23.65337 1.089779e-123
## ENSG00000211445 12285.6151       3.553363  22.49579 4.563626e-112
##                          padj      symbol      entrez
##                     <numeric> <character> <character>
## ENSG00000152583 4.000096e-132     SPARCL1        8404
## ENSG00000165995 7.154074e-131      CACNB2         783
## ENSG00000120129 2.170761e-126       DUSP1        1843
## ENSG00000101347 4.239316e-126      SAMHD1       25939
## ENSG00000189221 3.921897e-120        MAOA        4128
## ENSG00000211445 1.368631e-108        GPX3        2878

8.1 Exporting results

You can easily save the results table in a CSV file that you can then share or 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 that can be processed by write.csv. Here, we take just the top 100 genes for demonstration.

resOrderedDF <- as.data.frame(resOrdered)[1:100, ]
write.csv(resOrderedDF, file = "results.csv")

A more sophisticated way for exporting results the Bioconductor package ReportingTools (Huntley et al. 2013). ReportingTools will automatically generate dynamic HTML documents, including links to external databases using gene identifiers and boxplots summarizing the normalized counts across groups. See the ReportingTools vignettes for full details. The simplest version of creating a dynamic ReportingTools report is performed with the following code:

library("ReportingTools")
htmlRep <- HTMLReport(shortName = "report", title = "My report",
                      reportDirectory = "./report")
publish(resOrderedDF, htmlRep)
url <- finish(htmlRep)
browseURL(url)

8.2 Plotting fold changes in genomic space

If we have used the summarizeOverlaps function to count the reads, then our DESeqDataSet object is built on top of ready-to-use Bioconductor objects specifying the genomic coordinates of the genes. We can therefore easily plot our differential expression results in genomic space. While the results function by default returns a DataFrame, using the format argument, we can ask for GRanges or GRangesList output.

resGR <- results(dds, lfcThreshold = 1, format = "GRanges")
resGR
## GRanges object with 29391 ranges and 6 metadata columns:
##                   seqnames                 ranges strand |          baseMean
##                      <Rle>              <IRanges>  <Rle> |         <numeric>
##   ENSG00000000003        X [ 99883667,  99894988]      - |  708.602169691234
##   ENSG00000000419       20 [ 49551404,  49575092]      - |  520.297900552084
##   ENSG00000000457        1 [169818772, 169863408]      - |  237.163036796015
##   ENSG00000000460        1 [169631245, 169823221]      + |  57.9326331250967
##   ENSG00000000938        1 [ 27938575,  27961788]      - | 0.318098378392895
##               ...      ...                    ...    ... .               ...
##   ENSG00000273485       10 [105209953, 105210609]      + |  1.28644765243289
##   ENSG00000273486        3 [136556180, 136557863]      - |  15.4525365439045
##   ENSG00000273487        1 [ 92654794,  92656264]      + |   8.1632349843654
##   ENSG00000273488        3 [100080031, 100080481]      + |  8.58447903624707
##   ENSG00000273489        7 [131178723, 131182453]      - | 0.275899382507797
##                        log2FoldChange             lfcSE               stat
##                             <numeric>         <numeric>          <numeric>
##   ENSG00000000003  -0.381253879698533 0.100655966721595                  0
##   ENSG00000000419   0.206812707891892 0.112221794734617                  0
##   ENSG00000000457  0.0379205016307726  0.14345321468656                  0
##   ENSG00000000460 -0.0881632238896495 0.287167704842381                  0
##   ENSG00000000938     -1.378233965407  3.49987526565862 -0.108070698724122
##               ...                 ...               ...                ...
##   ENSG00000273485  -0.127136302037982  1.60055705381885                  0
##   ENSG00000273486  -0.150994443836803 0.486549011011151                  0
##   ENSG00000273487    1.04641689862923 0.699033601456549 0.0664015270975729
##   ENSG00000273488   0.107819010750248 0.638164612522881                  0
##   ENSG00000273489    1.48372584344306  3.51394515550545   0.13765890531478
##                              pvalue      padj
##                           <numeric> <numeric>
##   ENSG00000000003                 1         1
##   ENSG00000000419                 1         1
##   ENSG00000000457                 1         1
##   ENSG00000000460                 1         1
##   ENSG00000000938 0.913939611026568      <NA>
##               ...               ...       ...
##   ENSG00000273485                 1      <NA>
##   ENSG00000273486                 1         1
##   ENSG00000273487  0.94705815444012         1
##   ENSG00000273488                 1         1
##   ENSG00000273489 0.890509998916736      <NA>
##   -------
##   seqinfo: 722 sequences (1 circular) from an unspecified genome

We need to add the symbol again for labeling the genes on the plot:

resGR$symbol <- mapIds(org.Hs.eg.db, names(resGR), "SYMBOL", "ENSEMBL")

We will use the Gviz package for plotting the GRanges and associated metadata: the log fold changes due to dexamethasone treatment.

library("Gviz")

The following code chunk specifies a window of 1 million base pairs upstream and downstream from the gene with the smallest p value. We create a subset of our full results, for genes within the window. We add the gene symbol as a name if the symbol exists and is not duplicated in our subset.

window <- resGR[topGene] + 1e6
strand(window) <- "*"
resGRsub <- resGR[resGR %over% window]
naOrDup <- is.na(resGRsub$symbol) | duplicated(resGRsub$symbol)
resGRsub$group <- ifelse(naOrDup, names(resGRsub), resGRsub$symbol)

We create a vector specifying if the genes in this subset had a low value of padj.

status <- factor(ifelse(resGRsub$padj < 0.1 & !is.na(resGRsub$padj),
                     "sig", "notsig"))

We can then plot the results using Gviz functions (figure below). We create an axis track specifying our location in the genome, a track that will show the genes and their names, colored by significance, and a data track that will draw vertical bars showing the moderated log fold change produced by DESeq2, which we know are only large when the effect is well supported by the information in the counts.

options(ucscChromosomeNames = FALSE)
g <- GenomeAxisTrack()
a <- AnnotationTrack(resGRsub, name = "gene ranges", feature = status)
d <- DataTrack(resGRsub, data = "log2FoldChange", baseline = 0,
               type = "h", name = "log2 fold change", strand = "+")
plotTracks(list(g, d, a), groupAnnotation = "group",
           notsig = "grey", sig = "hotpink")

log2 fold changes in genomic region surrounding the gene with smallest adjusted p value. Genes highlighted in pink have adjusted p value less than 0.1.

9 Removing hidden batch effects

Suppose we did not know that there were different cell lines involved in the experiment, only that there was treatment with dexamethasone. The cell line effect on the counts then would represent some hidden and unwanted variation that might be affecting many or all of the genes in the dataset. We can use statistical methods designed for RNA-seq from the sva package (Leek 2014) to detect such groupings of the samples, and then we can add these to the DESeqDataSet design, in order to account for them. The SVA package uses the term surrogate variables for the estimated variables that we want to account for in our analysis. Another package for detecting hidden batches is the RUVSeq package (Risso et al. 2014), with the acronym “Remove Unwanted Variation”.

library("sva")

Below we obtain a matrix of normalized counts for which the average count across samples is larger than 1. As we described above, we are trying to recover any hidden batch effects, supposing that we do not know the cell line information. So we use a full model matrix with the dex variable, and a reduced, or null, model matrix with only an intercept term. Finally we specify that we want to estimate 2 surrogate variables. For more information read the manual page for the svaseq function by typing ?svaseq.

dat  <- counts(dds, normalized = TRUE)
idx  <- rowMeans(dat) > 1
dat  <- dat[idx, ]
mod  <- model.matrix(~ dex, colData(dds))
mod0 <- model.matrix(~   1, colData(dds))
svseq <- svaseq(dat, mod, mod0, n.sv = 2)
## Number of significant surrogate variables is:  2 
## Iteration (out of 5 ):1  2  3  4  5
svseq$sv
##            [,1]        [,2]
## [1,]  0.2481108 -0.52600157
## [2,]  0.2629867 -0.58115433
## [3,]  0.1502704  0.27428267
## [4,]  0.2023800  0.38419545
## [5,] -0.6086586 -0.07854931
## [6,] -0.6101210 -0.02923693
## [7,]  0.1788509  0.25708985
## [8,]  0.1761807  0.29937417

Because we actually do know the cell lines, we can see how well the SVA method did at recovering these variables (figure below).

par(mfrow = c(2, 1), mar = c(3,5,3,1))
for (i in 1:2) {
  stripchart(svseq$sv[, i] ~ dds$cell, vertical = TRUE, main = paste0("SV", i))
  abline(h = 0)
 }

Surrogate variables 1 and 2 plotted over cell line. Here, we know the hidden source of variation (cell line), and therefore can see how the SVA procedure is able to identify a source of variation which is correlated with cell line.

Finally, in order to use SVA to remove any effect on the counts from our surrogate variables, we simply add these two surrogate variables as columns to the DESeqDataSet and then add them to the design:

ddssva <- dds
ddssva$SV1 <- svseq$sv[,1]
ddssva$SV2 <- svseq$sv[,2]
design(ddssva) <- ~ SV1 + SV2 + dex

We could then produce results controlling for surrogate variables by running DESeq with the new design:

ddssva %<>% DESeq

10 Time course experiments

DESeq2 can be used to analyze time course experiments, for example to find those genes that react in a condition-specific manner over time, compared to a set of baseline samples. Here we demonstrate a basic time course analysis with the fission data package, which contains gene counts for an RNA-seq time course of fission yeast (Leong et al. 2014). The yeast were exposed to oxidative stress, and half of the samples contained a deletion of the gene atf21. We use a design formula that models the strain difference at time 0, the difference over time, and any strain-specific differences over time (the interaction term strain:minute).

library("fission")
data("fission")
ddsTC <- DESeqDataSet(fission, ~ strain + minute + strain:minute)

The following chunk of code performs a likelihood ratio test, where we remove the strain-specific differences over time. Genes with small p values from this test are those which at one or more time points after time 0 showed a strain-specific effect. Note therefore that this will not give small p values to genes that moved up or down over time in the same way in both strains.

ddsTC <- DESeq(ddsTC, test = "LRT", reduced = ~ strain + minute)
resTC <- results(ddsTC)
resTC$symbol <- mcols(ddsTC)$symbol
head(resTC[order(resTC$padj),], 4)
## log2 fold change (MLE): strainmut.minute180 
## LRT p-value: '~ strain + minute + strain:minute' vs '~ strain + minute' 
## DataFrame with 4 rows and 7 columns
##               baseMean log2FoldChange     lfcSE      stat       pvalue
##              <numeric>      <numeric> <numeric> <numeric>    <numeric>
## SPBC2F12.09c  174.6712    -2.65671953 0.7522613  97.28339 1.974151e-19
## SPAC1002.18   444.5050    -0.05093214 0.2042995  56.95360 5.169552e-11
## SPAC1002.19   336.3732    -0.39274898 0.5734940  43.53391 2.879804e-08
## SPAC1002.17c  261.7731    -1.13876477 0.6061288  39.31584 2.051371e-07
##                      padj      symbol
##                 <numeric> <character>
## SPBC2F12.09c 1.334526e-15       atf21
## SPAC1002.18  1.747308e-07        urg3
## SPAC1002.19  6.489157e-05        urg1
## SPAC1002.17c 3.466817e-04        urg2

This is just one of the tests that can be applied to time series data. Another option would be to model the counts as a smooth function of time, and to include an interaction term of the condition with the smooth function. It is possible to build such a model using spline basis functions within R, and another, more modern approach is using Gaussian processes (Tonner et al. 2017).

We can plot the counts for the groups over time using ggplot2, for the gene with the smallest adjusted p value, testing for condition-dependent time profile and accounting for differences at time 0 (figure below). Keep in mind that the interaction terms are the difference between the two groups at a given time after accounting for the difference at time 0.

fiss <- plotCounts(ddsTC, which.min(resTC$padj), 
                   intgroup = c("minute","strain"), returnData = TRUE)
ggplot(fiss,
  aes(x = as.numeric(minute), y = count, color = strain, group = strain)) + 
  geom_point() + geom_smooth(se = FALSE, method = "loess") + scale_y_log10()

Normalized counts for a gene with condition-specific changes over time.

Wald tests for the log2 fold changes at individual time points can be investigated using the test argument to results:

resultsNames(ddsTC)
##  [1] "Intercept"           "strain_mut_vs_wt"    "minute_15_vs_0"     
##  [4] "minute_30_vs_0"      "minute_60_vs_0"      "minute_120_vs_0"    
##  [7] "minute_180_vs_0"     "strainmut.minute15"  "strainmut.minute30" 
## [10] "strainmut.minute60"  "strainmut.minute120" "strainmut.minute180"
res30 <- results(ddsTC, name = "strainmut.minute30", test = "Wald")
res30[which.min(resTC$padj),]
## log2 fold change (MLE): strainmut.minute30 
## Wald test p-value: strainmut.minute30 
## DataFrame with 1 row and 6 columns
##               baseMean log2FoldChange     lfcSE      stat       pvalue
##              <numeric>      <numeric> <numeric> <numeric>    <numeric>
## SPBC2F12.09c  174.6712      -2.600469 0.6343429 -4.099469 4.140994e-05
##                   padj
##              <numeric>
## SPBC2F12.09c 0.2799312

We can furthermore cluster significant genes by their profiles. We extract a matrix of the shrunken log2 fold changes using the coef function:

betas <- coef(ddsTC)
colnames(betas)
##  [1] "Intercept"           "strain_mut_vs_wt"    "minute_15_vs_0"     
##  [4] "minute_30_vs_0"      "minute_60_vs_0"      "minute_120_vs_0"    
##  [7] "minute_180_vs_0"     "strainmut.minute15"  "strainmut.minute30" 
## [10] "strainmut.minute60"  "strainmut.minute120" "strainmut.minute180"

We can now plot the log2 fold changes in a heatmap (figure below).

topGenes <- head(order(resTC$padj),20)
mat <- betas[topGenes, -c(1,2)]
thr <- 3 
mat[mat < -thr] <- -thr
mat[mat > thr] <- thr
pheatmap(mat, breaks = seq(from = -thr, to = thr, length = 101),
         cluster_col = FALSE)

Heatmap of log2 fold changes for genes with smallest adjusted p value. The bottom set of genes show strong induction of expression for the baseline samples in minutes 15-60 (red boxes in the bottom left corner), but then have slight differences for the mutant strain (shown in the boxes in the bottom right corner).

11 Session information

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.

devtools::session_info()
##  setting  value                       
##  version  R version 3.4.0 (2017-04-21)
##  system   x86_64, darwin15.6.0        
##  ui       X11                         
##  language (EN)                        
##  collate  C                           
##  tz       Europe/Berlin               
##  date     2017-06-11                  
## 
##  package                * version  date       source        
##  acepack                  1.4.1    2016-10-29 CRAN (R 3.4.0)
##  affy                     1.54.0   2017-04-25 Bioconductor  
##  affyio                   1.46.0   2017-04-25 Bioconductor  
##  airway                 * 0.110.0  2017-04-28 Bioconductor  
##  annotate                 1.54.0   2017-04-25 Bioconductor  
##  AnnotationDbi          * 1.38.1   2017-06-01 Bioconductor  
##  AnnotationFilter         1.0.0    2017-04-25 Bioconductor  
##  AnnotationHub            2.8.1    2017-05-04 Bioconductor  
##  assertthat               0.2.0    2017-04-11 CRAN (R 3.4.0)
##  backports                1.1.0    2017-05-22 CRAN (R 3.4.0)
##  base                   * 3.4.0    2017-04-21 local         
##  base64enc                0.1-3    2015-07-28 CRAN (R 3.4.0)
##  beeswarm                 0.2.3    2016-04-25 CRAN (R 3.4.0)
##  Biobase                * 2.36.2   2017-05-04 Bioconductor  
##  BiocGenerics           * 0.22.0   2017-04-25 Bioconductor  
##  BiocInstaller            1.26.0   2017-04-25 Bioconductor  
##  BiocParallel           * 1.10.1   2017-05-03 Bioconductor  
##  BiocStyle              * 2.4.0    2017-04-25 Bioconductor  
##  biomaRt                  2.32.1   2017-06-09 Bioconductor  
##  Biostrings             * 2.44.1   2017-06-01 Bioconductor  
##  biovizBase               1.24.0   2017-04-25 Bioconductor  
##  bitops                   1.0-6    2013-08-17 CRAN (R 3.4.0)
##  bookdown                 0.4      2017-05-20 CRAN (R 3.4.0)
##  BSgenome                 1.44.0   2017-04-25 Bioconductor  
##  checkmate                1.8.2    2016-11-02 CRAN (R 3.4.0)
##  cluster                  2.0.6    2017-03-10 CRAN (R 3.4.0)
##  colorspace               1.3-2    2016-12-14 CRAN (R 3.4.0)
##  compiler                 3.4.0    2017-04-21 local         
##  curl                     2.6      2017-04-27 CRAN (R 3.4.0)
##  data.table               1.10.4   2017-02-01 CRAN (R 3.4.0)
##  datasets               * 3.4.0    2017-04-21 local         
##  DBI                      0.6-1    2017-04-01 CRAN (R 3.4.0)
##  DelayedArray           * 0.2.7    2017-06-03 Bioconductor  
##  DESeq2                 * 1.16.1   2017-05-06 Bioconductor  
##  devtools                 1.13.2   2017-06-02 CRAN (R 3.4.0)
##  dichromat                2.0-0    2013-01-24 CRAN (R 3.4.0)
##  digest                   0.6.12   2017-01-27 CRAN (R 3.4.0)
##  dplyr                  * 0.7.0    2017-06-09 CRAN (R 3.4.0)
##  edgeR                  * 3.18.1   2017-05-06 Bioconductor  
##  ensembldb                2.0.3    2017-05-27 Bioconductor  
##  evaluate                 0.10     2016-10-11 CRAN (R 3.4.0)
##  fission                * 0.110.0  2017-06-11 Bioconductor  
##  foreign                  0.8-68   2017-04-24 CRAN (R 3.4.0)
##  Formula                  1.2-1    2015-04-07 CRAN (R 3.4.0)
##  genefilter             * 1.58.1   2017-05-06 Bioconductor  
##  geneplotter              1.54.0   2017-04-25 Bioconductor  
##  GenomeInfoDb           * 1.12.2   2017-06-09 Bioconductor  
##  GenomeInfoDbData         0.99.0   2017-04-28 Bioconductor  
##  GenomicAlignments      * 1.12.1   2017-05-12 Bioconductor  
##  GenomicFeatures        * 1.28.3   2017-06-09 Bioconductor  
##  GenomicRanges          * 1.28.3   2017-05-25 Bioconductor  
##  ggbeeswarm             * 0.5.3    2016-12-01 CRAN (R 3.4.0)
##  ggplot2                * 2.2.1    2016-12-30 CRAN (R 3.4.0)
##  glue                     1.0.0    2017-04-17 CRAN (R 3.4.0)
##  graphics               * 3.4.0    2017-04-21 local         
##  grDevices              * 3.4.0    2017-04-21 local         
##  grid                   * 3.4.0    2017-04-21 local         
##  gridExtra                2.2.1    2016-02-29 CRAN (R 3.4.0)
##  gtable                   0.2.0    2016-02-26 CRAN (R 3.4.0)
##  Gviz                   * 1.20.0   2017-04-25 Bioconductor  
##  hexbin                 * 1.27.1   2016-12-05 CRAN (R 3.4.0)
##  Hmisc                    4.0-3    2017-05-02 CRAN (R 3.4.0)
##  htmlTable                1.9      2017-01-26 CRAN (R 3.4.0)
##  htmltools                0.3.6    2017-04-28 CRAN (R 3.4.0)
##  htmlwidgets              0.8      2016-11-09 CRAN (R 3.4.0)
##  httpuv                   1.3.3    2015-08-04 CRAN (R 3.4.0)
##  httr                     1.2.1    2016-07-03 CRAN (R 3.4.0)
##  interactiveDisplayBase   1.14.0   2017-04-25 Bioconductor  
##  IRanges                * 2.10.2   2017-05-25 Bioconductor  
##  knitr                  * 1.16     2017-05-18 CRAN (R 3.4.0)
##  labeling                 0.3      2014-08-23 CRAN (R 3.4.0)
##  lattice                  0.20-35  2017-03-25 CRAN (R 3.4.0)
##  latticeExtra             0.6-28   2016-02-09 CRAN (R 3.4.0)
##  lazyeval                 0.2.0    2016-06-12 CRAN (R 3.4.0)
##  limma                  * 3.32.2   2017-05-02 Bioconductor  
##  locfit                   1.5-9.1  2013-04-20 CRAN (R 3.4.0)
##  magrittr               * 1.5      2014-11-22 CRAN (R 3.4.0)
##  Matrix                   1.2-10   2017-04-28 CRAN (R 3.4.0)
##  matrixStats            * 0.52.2   2017-04-14 CRAN (R 3.4.0)
##  memoise                  1.1.0    2017-04-21 CRAN (R 3.4.0)
##  methods                * 3.4.0    2017-04-21 local         
##  mgcv                   * 1.8-17   2017-02-08 CRAN (R 3.4.0)
##  mime                     0.5      2016-07-07 CRAN (R 3.4.0)
##  munsell                  0.4.3    2016-02-13 CRAN (R 3.4.0)
##  nlme                   * 3.1-131  2017-02-06 CRAN (R 3.4.0)
##  nnet                     7.3-12   2016-02-02 CRAN (R 3.4.0)
##  org.Hs.eg.db           * 3.4.1    2017-04-28 Bioconductor  
##  parallel               * 3.4.0    2017-04-21 local         
##  pheatmap               * 1.0.8    2015-12-11 CRAN (R 3.4.0)
##  plyr                     1.8.4    2016-06-08 CRAN (R 3.4.0)
##  PoiClaClu              * 1.0.2    2013-12-02 CRAN (R 3.4.0)
##  preprocessCore           1.38.1   2017-05-06 Bioconductor  
##  ProtGenerics             1.8.0    2017-04-25 Bioconductor  
##  R6                       2.2.1    2017-05-10 CRAN (R 3.4.0)
##  RColorBrewer           * 1.1-2    2014-12-07 CRAN (R 3.4.0)
##  Rcpp                     0.12.11  2017-05-22 CRAN (R 3.4.0)
##  RCurl                    1.95-4.8 2016-03-01 CRAN (R 3.4.0)
##  reshape2                 1.4.2    2016-10-22 CRAN (R 3.4.0)
##  rlang                    0.1.1    2017-05-18 CRAN (R 3.4.0)
##  rmarkdown              * 1.5      2017-04-26 CRAN (R 3.4.0)
##  rpart                    4.1-11   2017-03-13 CRAN (R 3.4.0)
##  rprojroot                1.2      2017-01-16 CRAN (R 3.4.0)
##  Rsamtools              * 1.28.0   2017-04-25 Bioconductor  
##  RSQLite                  1.1-2    2017-01-08 CRAN (R 3.4.0)
##  rtracklayer              1.36.3   2017-05-25 Bioconductor  
##  S4Vectors              * 0.14.3   2017-06-03 Bioconductor  
##  scales                   0.4.1    2016-11-09 CRAN (R 3.4.0)
##  shiny                    1.0.3    2017-04-26 CRAN (R 3.4.0)
##  splines                  3.4.0    2017-04-21 local         
##  stats                  * 3.4.0    2017-04-21 local         
##  stats4                 * 3.4.0    2017-04-21 local         
##  stringi                  1.1.5    2017-04-07 CRAN (R 3.4.0)
##  stringr                  1.2.0    2017-02-18 CRAN (R 3.4.0)
##  SummarizedExperiment   * 1.6.3    2017-05-29 Bioconductor  
##  survival                 2.41-3   2017-04-04 CRAN (R 3.4.0)
##  sva                    * 3.24.0   2017-04-25 Bioconductor  
##  tibble                   1.3.3    2017-05-28 CRAN (R 3.4.0)
##  tools                    3.4.0    2017-04-21 local         
##  utils                  * 3.4.0    2017-04-21 local         
##  VariantAnnotation        1.22.1   2017-05-25 Bioconductor  
##  vipor                    0.4.5    2017-03-22 CRAN (R 3.4.0)
##  vsn                    * 3.44.0   2017-04-25 Bioconductor  
##  withr                    1.0.2    2016-06-20 CRAN (R 3.4.0)
##  XML                      3.98-1.7 2017-05-03 CRAN (R 3.4.0)
##  xtable                   1.8-2    2016-02-05 CRAN (R 3.4.0)
##  XVector                * 0.16.0   2017-04-25 Bioconductor  
##  yaml                     2.1.14   2016-11-12 CRAN (R 3.4.0)
##  zlibbioc                 1.22.0   2017-04-25 Bioconductor

References

Anders, Simon, and Wolfgang Huber. 2010. “Differential expression analysis for sequence count data.” Genome Biology 11 (10). BioMed Central Ltd: R106+. doi:10.1186/gb-2010-11-10-r106.

Anders, Simon, Paul T. Pyl, and Wolfgang Huber. 2015. “HTSeq – a Python framework to work with high-throughput sequencing data.” Bioinformatics 31 (2). Oxford University Press: 166–69. doi:10.1093/bioinformatics/btu638.

Benjamini, Yoav, and Yosef Hochberg. 1995. “Controlling the False Discovery Rate: A Practical and Powerful Approach to Multiple Testing.” Journal of the Royal Statistical Society. Series B (Methodological) 57 (1): 289–300. http://www.jstor.org/stable/2346101.

Bourgon, R., R. Gentleman, and W. Huber. 2010. “Independent filtering increases detection power for high-throughput experiments.” Proceedings of the National Academy of Sciences 107 (21). National Academy of Sciences: 9546–51. doi:10.1073/pnas.0914005107.

Bray, Nicolas L, Harold Pimentel, Páll Melsted, and Lior Pachter. 2016. “Near-Optimal RNA-Seq Quantification.” Nat. Biotechnol.

Dobin, Alexander, Carrie A. Davis, Felix Schlesinger, Jorg Drenkow, Chris Zaleski, Sonali Jha, Philippe Batut, Mark Chaisson, and Thomas R. Gingeras. 2013. “STAR: ultrafast universal RNA-seq aligner.” Bioinformatics 29 (1). Oxford University Press: 15–21. doi:10.1093/bioinformatics/bts635.

Dudoit, Rine, Yee H. Yang, Matthew J. Callow, and Terence P. Speed. 2002. “Statistical methods for identifying differentially expressed genes in replicated cDNA microarray experiments.” In Statistica Sinica, 111–39. http://citeseerx.ist.psu.edu/viewdoc/summary?doi=10.1.1.117.9702.

Flicek, Paul, M. Ridwan Amode, Daniel Barrell, Kathryn Beal, Konstantinos Billis, Simon Brent, Denise Carvalho-Silva, et al. 2014. “Ensembl 2014.” Nucleic Acids Research 42 (D1). Oxford University Press: D749–D755. doi:10.1093/nar/gkt1196.

Hardcastle, Thomas, and Krystyna Kelly. 2010. “baySeq: Empirical Bayesian methods for identifying differential expression in sequence count data.” BMC Bioinformatics 11 (1): 422+. doi:10.1186/1471-2105-11-422.

Himes, Blanca E., Xiaofeng Jiang, Peter Wagner, Ruoxi Hu, Qiyu Wang, Barbara Klanderman, Reid M. Whitaker, et al. 2014. “RNA-Seq transcriptome profiling identifies CRISPLD2 as a glucocorticoid responsive gene that modulates cytokine function in airway smooth muscle cells.” PloS One 9 (6). doi:10.1371/journal.pone.0099625.

Huber, Wolfgang, Vincent J. Carey, Robert Gentleman, Simon Anders, Marc Carlson, Benilton S. Carvalho, Hector Corrada C. Bravo, et al. 2015. “Orchestrating high-throughput genomic analysis with Bioconductor.” Nature Methods 12 (2). Nature Publishing Group: 115–21. doi:10.1038/nmeth.3252.

Huntley, Melanie A., Jessica L. Larson, Christina Chaivorapol, Gabriel Becker, Michael Lawrence, Jason A. Hackney, and Joshua S. Kaminker. 2013. “ReportingTools: an automated result processing and presentation toolkit for high-throughput genomic analyses.” Bioinformatics 29 (24). Oxford University Press: 3220–1. doi:10.1093/bioinformatics/btt551.

Kent, W. James, Charles W. Sugnet, Terrence S. Furey, Krishna M. Roskin, Tom H. Pringle, Alan M. Zahler, and David Haussler. 2002. “The human genome browser at UCSC.” Genome Research 12 (6). Cold Spring Harbor Laboratory Press: 996–1006. doi:10.1101/gr.229102.

Law, Charity W., Yunshun Chen, Wei Shi, and Gordon K. Smyth. 2014. “Voom: precision weights unlock linear model analysis tools for RNA-seq read counts.” Genome Biology 15 (2). BioMed Central Ltd: R29+. doi:10.1186/gb-2014-15-2-r29.

Lawrence, Michael, Wolfgang Huber, Hervé Pagès, Patrick Aboyoun, Marc Carlson, Robert Gentleman, Martin T. Morgan, and Vincent J. Carey. 2013. “Software for Computing and Annotating Genomic Ranges.” Edited by Andreas Prlic. PLoS Computational Biology 9 (8). Public Library of Science: e1003118+. doi:10.1371/journal.pcbi.1003118.

Leek, Jeffrey T. 2014. “svaseq: removing batch effects and other unwanted noise from sequencing data.” Nucleic Acids Research 42 (21). Oxford University Press: 000. doi:10.1093/nar/gku864.

Leng, N., J. A. Dawson, J. A. Thomson, V. Ruotti, A. I. Rissman, B. M. G. Smits, J. D. Haag, M. N. Gould, R. M. Stewart, and C. Kendziorski. 2013. “EBSeq: an empirical Bayes hierarchical model for inference in RNA-seq experiments.” Bioinformatics 29 (8). Oxford University Press: 1035–43. doi:10.1093/bioinformatics/btt087.

Leong, Hui S., Keren Dawson, Chris Wirth, Yaoyong Li, Yvonne Connolly, Duncan L. Smith, Caroline R. Wilkinson, and Crispin J. Miller. 2014. “A global non-coding RNA system modulates fission yeast protein levels in response to stress.” Nature Communications 5. doi:10.1038/ncomms4947.

Li, Bo, and Colin N. Dewey. 2011. “RSEM: accurate transcript quantification from RNA-Seq data with or without a reference genome.” BMC Bioinformatics 12: 323+. doi:10.1186/1471-2105-12-3231.

Li, Heng, Bob Handsaker, Alec Wysoker, Tim Fennell, Jue Ruan, Nils Homer, Gabor Marth, Goncalo Abecasis, Richard Durbin, and 1000 Genome Project Data Processing Subgroup. 2009. “The Sequence Alignment/Map format and SAMtools.” Bioinformatics (Oxford, England) 25 (16). Oxford University Press: 2078–9. doi:10.1093/bioinformatics/btp352.

Liao, Y., G. K. Smyth, and W. Shi. 2014. “featureCounts: an efficient general purpose program for assigning sequence reads to genomic features.” Bioinformatics 30 (7). Oxford University Press: 923–30. doi:10.1093/bioinformatics/btt656.

Love, Michael I., Wolfgang Huber, and Simon Anders. 2014. “Moderated estimation of fold change and dispersion for RNA-seq data with DESeq2.” Genome Biology 15 (12). BioMed Central Ltd: 550+. doi:10.1186/s13059-014-0550-8.

Patro, Rob, Geet Duggal, Michael I Love, Rafael A Irizarry, and Carl Kingsford. 2017. “Salmon Provides Fast and Bias-Aware Quantification of Transcript Expression.” Nat. Methods 14 (6~mar): 417–19.

Patro, Rob, Stephen M. Mount, and Carl Kingsford. 2014. “Sailfish enables alignment-free isoform quantification from RNA-seq reads using lightweight algorithms.” Nature Biotechnology 32: 462–64. doi:10.1038/nbt.2862.

Risso, Davide, John Ngai, Terence P. Speed, and Sandrine Dudoit. 2014. “Normalization of RNA-seq data using factor analysis of control genes or samples.” Nature Biotechnology 32 (9). Nature Publishing Group: 896–902. doi:10.1038/nbt.2931.

Robert, Christelle, and Mick Watson. 2015. “Errors in RNA-Seq quantification affect genes of relevance to human disease.” Genome Biology. doi:10.1186/s13059-015-0734-x.

Robinson, M. D., D. J. McCarthy, and G. K. Smyth. 2009. “edgeR: a Bioconductor package for differential expression analysis of digital gene expression data.” Bioinformatics 26 (1). Oxford University Press: 139–40. doi:10.1093/bioinformatics/btp616.

Schurch, Nicholas J, Pietá Schofield, Marek Gierliński, Christian Cole, Alexander Sherstnev, Vijender Singh, Nicola Wrobel, et al. 2016. “How Many Biological Replicates Are Needed in an RNA-seq Experiment and Which Differential Expression Tool Should You Use?” RNA 22 (6): 839–51.

Soneson, Charlotte, Michael I. Love, and Mark Robinson. 2015. “Differential analyses for RNA-seq: transcript-level estimates improve gene-level inferences.” F1000Research 4 (1521). doi:10.12688/f1000research.7563.1.

Tonner, Peter D, Cynthia L Darnell, Barbara E Engelhardt, and Amy K Schmid. 2017. “Detecting Differential Growth of Microbial Populations with Gaussian Process Regression.” Genome Res. 27 (2): 320–33.

Trapnell, Cole, David G Hendrickson, Martin Sauvageau, Loyal Goff, John L Rinn, and Lior Pachter. 2013. “Differential analysis of gene regulation at transcript resolution with RNA-seq.” Nature Biotechnology. doi:10.1038/nbt.2450.

Wickham, Hadley. 2009. ggplot2. New York, NY: Springer New York. doi:10.1007/978-0-387-98141-3.

Witten, Daniela M. 2011. “Classification and clustering of sequencing data using a Poisson model.” The Annals of Applied Statistics 5 (4): 2493–2518. doi:10.1214/11-AOAS493.

Wu, Hao, Chi Wang, and Zhijin Wu. 2013. “A new shrinkage estimator for dispersion improves differential expression detection in RNA-seq data.” Biostatistics 14 (2). Oxford University Press: 232–43. doi:10.1093/biostatistics/kxs033.