author: Michael I. Love, Simon Anders, Wolfgang Huber date: August 14, 2015 title: “RNA-Seq workflow: gene-level exploratory analysis and differential expression” output: html_document
Michael Love [1], Simon Anders [2], Wolfgang Huber [2]
[1] Department of Biostatistics, Dana-Farber Cancer Institute and Harvard School of Public Health, Boston, US;
[2] European Molecular Biology Laboratory (EMBL), Heidelberg, Germany.
This lab will walk you through an end-to-end gene-level RNA-Seq differential expression workflow, using DESeq2 along with other Bioconductor packages. We will start from the FASTQ files, show how these were aligned to the reference genome, prepare a count matrix which tallies the number of RNA-seq reads/fragments within each gene for each sample, perform exploratory data analysis (EDA), perform differential gene expression analysis, and visually explore the results.
We note that a number of other Bioconductor packages can also be used for statistical inference of differential expression at the gene level including edgeR, limma, DSS, EBSeq, and BaySeq.
If you have questions about this workflow, please post them 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 deseq2
. Note the
posting guide
for crafting an optimal question for the support site.
The data used in this workflow is stored in the airway package, which summarizes an RNA-Seq experiment wherein airway smooth muscle cells were treated with dexamethasone, a synthetic glucocorticoid steroid with anti-inflammatory effects. Glucocorticoids are used, for example, in asthma patients 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. The reference for the experiment is:
Himes BE, Jiang X, Wagner P, Hu R, Wang Q, Klanderman B, Whitaker RM, Duan Q, Lasky-Su J, Nikolos C, Jester W, Johnson M, Panettieri R Jr, Tantisira KG, Weiss ST, Lu Q. “RNA-Seq Transcriptome Profiling Identifies CRISPLD2 as a Glucocorticoid Responsive Gene that Modulates Cytokine Function in Airway Smooth Muscle Cells.” PLoS One. 2014 Jun 13;9(6):e99625. PMID: 24926665. GEO: GSE52778.
As input, the DESeq2 statistical method expects count data as obtained, e.g., from RNA-Seq or another high-throughput sequencing experiment, in the form of a matrix of integer values. The value in the i-th row and the j-th column of the matrix tells how many reads (or fragments, for paired-end RNA-seq) can be uniquely assigend 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 must be raw counts of sequencing reads/fragments. This is important for DESeq2's statistical model to hold, as only the raw counts allow assessing the measurement precision correctly. It is important to not provide counts which were pre-normalized for sequencing depth or library size, as the statistical model is most powerful when applied to raw counts, and is designed to account for library size differences internally.
The computational analysis of an RNA-Seq experiment begins earlier: we first obtain a set of 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. 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 the development is too rapid for this document to provide an up-to-date list. We recommend consulting benchmarking papers that discuss the advantages and disadvantages of each software, which include accuracy, sensitivity in aligning reads over splice junctions, speed, memory footprint, and many other features.
The reads for this experiment were aligned to the Ensembl release 75 human reference genome using the STAR read aligner:
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 was used to generate BAM files.
cat files | parallel -j 7 samtools view -bS aligned/{}.Aligned.out.sam -o aligned/{}.bam
The BAM files for a number of sequencing runs can then be used to generate count matrices, as described in the following section.
Besides the count matrix, which 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. 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, which is part of the airway package.
dir <- system.file("extdata", package="airway", mustWork=TRUE)
In this directory, we find the eight BAM files (and some other files):
list.files(dir)
## [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, and which 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(dir,"sample_table.csv")
(sampleTable <- read.csv(csvfile,row.names=1))
## SampleName cell dex albut Run avgLength Experiment Sample BioSample
## SRR1039508 GSM1275862 N61311 untrt untrt SRR1039508 126 SRX384345 SRS508568 SAMN02422669
## SRR1039509 GSM1275863 N61311 trt untrt SRR1039509 126 SRX384346 SRS508567 SAMN02422675
## SRR1039512 GSM1275866 N052611 untrt untrt SRR1039512 126 SRX384349 SRS508571 SAMN02422678
## SRR1039513 GSM1275867 N052611 trt untrt SRR1039513 87 SRX384350 SRS508572 SAMN02422670
## SRR1039516 GSM1275870 N080611 untrt untrt SRR1039516 120 SRX384353 SRS508575 SAMN02422682
## SRR1039517 GSM1275871 N080611 trt untrt SRR1039517 126 SRX384354 SRS508576 SAMN02422673
## SRR1039520 GSM1275874 N061011 untrt untrt SRR1039520 101 SRX384357 SRS508579 SAMN02422683
## SRR1039521 GSM1275875 N061011 trt untrt SRR1039521 98 SRX384358 SRS508580 SAMN02422677
Note: here and elsewhere in the workflow, the parentheses ()
around the
entire code of the last line above is an R trick to print the output of the function
in addition to saving it to sampleTable
. This is equivalent to assigning
and then showing the object in two steps:
sampleTable <- read.csv(csvfile,row.names=1)
sampleTable
Once the reads have been aligned, there are a number of tools which can be used to count the number of reads/fragments which can be uniquely 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.
The following tools can be used generate count matrices:
function | package | framework | output | DESeq2 input function |
---|---|---|---|---|
summarizeOverlaps | GenomicAlignments | R/Bioconductor | SummarizedExperiment | DESeqDataSet |
featureCounts | Rsubread | R/Bioconductor | matrix | DESeqDataSetFromMatrix |
htseq-count | HTSeq | Python | files | DESeqDataSetFromHTSeq |
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(dir, 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. 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. 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>
Next, we need to read in the gene model which will be used for counting reads/fragments. We will read the gene model from a GTF file, 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 which 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, 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. Finally, the makeTxDbFromBiomart function can be used to automatically pull a gene model from Biomart.
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(dir,"Homo_sapiens.GRCh37.75_subset.gtf")
(txdb <- makeTxDbFromGFF(gtffile, format="gtf", circ_seqs=character()))
## TxDb object:
## # Db type: TxDb
## # Supporting package: GenomicFeatures
## # Data source: /var/lib/jenkins/R/x86_64-unknown-linux-gnu-library/3.2/airway/extdata/Homo_sapiens.GRCh37.75_subset.gtf
## # Organism: 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: 2015-08-16 06:57:08 -0700 (Sun, 16 Aug 2015)
## # GenomicFeatures version at creation time: 1.20.2
## # RSQLite version at creation time: 1.0.0
## # DBSCHEMAVERSION: 1.1
The following line produces a GRangesList of all the exons grouped by gene.
(ebg <- exonsBy(txdb, by="gene"))
## 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] - | 103 ENSE00003520086
## ... ... ... ... ... ... ...
## [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
After these preparations, the actual counting is easy. The function summarizeOverlaps from the GenomicAlignments package will do this. This produces a SummarizedExperiment object, which 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 20 million 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 it would also run in serial.
register(SerialParam())
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 r
Biocpkg("GenomicAlignments")
package. 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.
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
. The fragments
argument can be used when singleEnd=FALSE
to specify if unpaired
reads should be counted (yes if fragments=TRUE
).
Here we show the component parts of a SummarizedExperiment object,
and the classes derived from it, such as the DESeqDataSet which is
explained in the next section. The assay(s)
(pink block) contains
the matrix (or matrices) of data. In our case we have created a single
matrix which contains the fragment counts for each gene and sample.
The rowRanges
(blue block) contains information about the genomic
ranges, and the colData
(green block) contains information about the
samples or experiments. 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
).
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: SummarizedExperiment
## dim: 20 8
## exptData(0):
## assays(1): counts
## rownames(20): ENSG00000009724 ENSG00000116649 ... ENSG00000271794 ENSG00000271895
## rowRanges metadata column 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))
## SRR1039508_subset.bam SRR1039509_subset.bam SRR1039512_subset.bam
## ENSG00000009724 38 28 66
## ENSG00000116649 1004 1255 1122
## ENSG00000120942 218 256 233
## ENSG00000120948 2751 2080 3353
## ENSG00000171819 4 50 19
## ENSG00000171824 869 1075 1115
## SRR1039513_subset.bam SRR1039516_subset.bam SRR1039517_subset.bam
## ENSG00000009724 24 42 41
## ENSG00000116649 1313 1100 1879
## ENSG00000120942 252 269 465
## ENSG00000120948 1614 3519 3716
## ENSG00000171819 543 1 10
## ENSG00000171824 1051 944 1405
## SRR1039520_subset.bam SRR1039521_subset.bam
## ENSG00000009724 47 36
## ENSG00000116649 745 1536
## ENSG00000120942 207 400
## ENSG00000120948 2220 1990
## ENSG00000171819 14 1067
## ENSG00000171824 748 1590
colSums(assay(se))
## SRR1039508_subset.bam SRR1039509_subset.bam SRR1039512_subset.bam SRR1039513_subset.bam
## 6478 6501 7699 6801
## SRR1039516_subset.bam SRR1039517_subset.bam SRR1039520_subset.bam SRR1039521_subset.bam
## 8009 10849 5254 9168
The rowRanges
have their own accessor function:
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] - | 103 ENSE00003520086
## ... ... ... ... ... ... ...
## [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
Note that the rowRanges
slot is a GRangesList, which
contains all the information about the exons for each gene, i.e., for each row
of the count matrix. It 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 14
## ..$ Db type : chr "TxDb"
## ..$ Supporting package : chr "GenomicFeatures"
## ..$ Data source : chr "/var/lib/jenkins/R/x86_64-unknown-linux-gnu-library/3.2/airway/extdata/Homo_sapiens.GRCh37.75_subset.gtf"
## ..$ Organism : 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 "2015-08-16 06:57:08 -0700 (Sun, 16 Aug 2015)"
## ..$ GenomicFeatures version at creation time: chr "1.20.2"
## ..$ RSQLite version at creation time : chr "1.0.0"
## ..$ DBSCHEMAVERSION : chr "1.1"
The colData
also has its own accessor function:
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))
## DataFrame with 8 rows and 9 columns
## SampleName cell dex albut Run avgLength Experiment Sample
## <factor> <factor> <factor> <factor> <factor> <integer> <factor> <factor>
## SRR1039508 GSM1275862 N61311 untrt untrt SRR1039508 126 SRX384345 SRS508568
## SRR1039509 GSM1275863 N61311 trt untrt SRR1039509 126 SRX384346 SRS508567
## SRR1039512 GSM1275866 N052611 untrt untrt SRR1039512 126 SRX384349 SRS508571
## SRR1039513 GSM1275867 N052611 trt untrt SRR1039513 87 SRX384350 SRS508572
## SRR1039516 GSM1275870 N080611 untrt untrt SRR1039516 120 SRX384353 SRS508575
## SRR1039517 GSM1275871 N080611 trt untrt SRR1039517 126 SRX384354 SRS508576
## SRR1039520 GSM1275874 N061011 untrt untrt SRR1039520 101 SRX384357 SRS508579
## SRR1039521 GSM1275875 N061011 trt untrt SRR1039521 98 SRX384358 SRS508580
## BioSample
## <factor>
## SRR1039508 SAMN02422669
## SRR1039509 SAMN02422675
## SRR1039512 SAMN02422678
## SRR1039513 SAMN02422670
## SRR1039516 SAMN02422682
## SRR1039517 SAMN02422673
## SRR1039520 SAMN02422683
## SRR1039521 SAMN02422677
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, limma, DSS, EBSeq, BaySeq, and more. We will continue, using DESeq2. 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.
Bioconductor software packages often define and use a custom class for storing data, which 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 accidentally 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)
which
specifies which of two (or more groups) the samples belong to. For
the airway experiment, we will specify ~ cell + dex
, which
means 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:
se$dex <- relevel(se$dex, "untrt")
se$dex
## [1] untrt trt untrt trt untrt trt untrt trt
## Levels: untrt trt
For running DESeq2 models, you can use R's formula notation to
express any fixed-effects experimental design.
Note that DESeq2 uses 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 a time in a
time series 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.
We now use R's data command to load a prepared
SummarizedExperiment that was generated from the publicly available
sequencing data files associated with the Himes et al. paper,
described above. The steps we used to produce this object were
equivalent to those you worked through in the previous sections,
except that we used all the reads and all the genes. For more details
on the exact steps used to create this object, type
vignette("airway")
into your R session.
data("airway")
se <- airway
We can quickly check the millions of fragments which 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 SRR1039521
## 20.6 18.8 25.3 15.2 24.4 30.8 19.1 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 Experiment Sample
## <factor> <factor> <factor> <factor> <factor> <integer> <factor> <factor>
## SRR1039508 GSM1275862 N61311 untrt untrt SRR1039508 126 SRX384345 SRS508568
## SRR1039509 GSM1275863 N61311 trt untrt SRR1039509 126 SRX384346 SRS508567
## SRR1039512 GSM1275866 N052611 untrt untrt SRR1039512 126 SRX384349 SRS508571
## SRR1039513 GSM1275867 N052611 trt untrt SRR1039513 87 SRX384350 SRS508572
## SRR1039516 GSM1275870 N080611 untrt untrt SRR1039516 120 SRX384353 SRS508575
## SRR1039517 GSM1275871 N080611 trt untrt SRR1039517 126 SRX384354 SRS508576
## SRR1039520 GSM1275874 N061011 untrt untrt SRR1039520 101 SRX384357 SRS508579
## SRR1039521 GSM1275875 N061011 trt untrt SRR1039521 98 SRX384358 SRS508580
## BioSample
## <factor>
## SRR1039508 SAMN02422669
## SRR1039509 SAMN02422675
## SRR1039512 SAMN02422678
## SRR1039513 SAMN02422670
## SRR1039516 SAMN02422682
## SRR1039517 SAMN02422673
## SRR1039520 SAMN02422683
## SRR1039521 SAMN02422677
Here we see that this object already contains an informative
colData
slot – because we have already prepared it for you, as
described in the airway vignette.
However, when you work with your own data, you will have to add the
pertinent sample / phenotypic information for the experiment at this stage.
We highly recommend keeping this information in a comma-separated
value (CSV) or tab-separated value (TSV) file, which can be exported
from an Excel spreadsheet, and the assign this to the colData
slot,
making sure that the rows correspond to the columns of the
SummarizedExperiment. We made sure of this correspondence 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, which 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)
If we only wanted to perform transformations and exploratory data analysis
we could use a ~ 1
for the design, but we would need to remember to
substitute a real design, e.g. ~ condition
, before we run DESeq
for differential testing or else we would only be testing the intercept.
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 contruct 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.
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)
## SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520
## ENSG00000000003 679 448 873 408 1138 1047 770
## ENSG00000000005 0 0 0 0 0 0 0
## ENSG00000000419 467 515 621 365 587 799 417
## ENSG00000000457 260 211 263 164 245 331 233
## ENSG00000000460 60 55 40 35 78 63 76
## ENSG00000000938 0 0 2 0 1 0 0
## SRR1039521
## ENSG00000000003 572
## ENSG00000000005 0
## ENSG00000000419 508
## ENSG00000000457 229
## ENSG00000000460 60
## ENSG00000000938 0
In this count matrix, each row represents an Ensembl gene, each column a sequenced RNA library, and the values give the raw numbers of 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 importnat to 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 countscoldata
: a table with information about the samplesTo 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))
## class: DESeqDataSet
## dim: 64102 8
## exptData(0):
## assays(1): counts
## rownames(64102): ENSG00000000003 ENSG00000000005 ... LRG_98 LRG_99
## rowRanges metadata column names(0):
## colnames(8): SRR1039508 SRR1039509 ... SRR1039520 SRR1039521
## colData names(9): SampleName cell ... Sample BioSample
We will continue with the object generated from the SummarizedExperiment section.
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 which have no or nearly no information about the amount of gene expression. Here we remove rows of the DESeqDataSet which have no counts, or only a single count across all samples:
dds <- dds[ rowSums(counts(dds)) > 1, ]
Many common statistical methods for exploratory analysis of multidimensional data, especially methods for clustering and ordination (e.g., principal-component analysis and the like), work best for (at least approximately) homoskedastic data; this means that the variance of an observed quantity (here, the expression strength of a gene) does not depend on the mean. In RNA-Seq data, however, variance grows with the mean. For example, if one performs PCA (principal components analysis) directly on a matrix of normalized read counts, the result typically depends only on the few most strongly expressed genes 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 small pseudocount; however, now the genes with low counts tend to dominate the results because, due to the strong Poisson noise inherent to small count values, they show the strongest relative differences between samples.
As a solution, DESeq2 offers the regularized-logarithm transformation,
or rlog for short. For genes with high counts, the rlog
transformation differs not much from an ordinary log2
transformation. For genes with lower counts, however, the values are
shrunken towards the genes' averages across all samples. Using
an empirical Bayesian prior on inter-sample differences in the form of
a ridge penalty, this is done such that the rlog-transformed data
are approximately homoskedastic. See the help for ?rlog
for more
information and options. Another transformation, the variance
stabilizing transformation, is discussed alongside the rlog in the
DESeq2 vignette.
Note: the rlog transformation is 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 a SummarizedExperiment object which contains the rlog-transformed values in its assay slot:
rld <- rlog(dds)
head(assay(rld))
## SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520
## ENSG00000000003 9.399171 9.142606 9.501676 9.320847 9.757103 9.512160 9.617316
## ENSG00000000419 8.901330 9.113940 9.032562 9.063908 8.981945 9.108498 8.894880
## ENSG00000000457 7.949872 7.882373 7.834296 7.916446 7.773868 7.886645 7.946387
## ENSG00000000460 5.849514 5.882338 5.487152 5.770373 5.940353 5.663945 6.107597
## ENSG00000000938 -1.638030 -1.637430 -1.558279 -1.636020 -1.597596 -1.639307 -1.637555
## ENSG00000000971 11.927571 12.159625 12.387352 12.571706 12.469397 12.784626 12.454894
## SRR1039521
## ENSG00000000003 9.315363
## ENSG00000000419 9.052290
## ENSG00000000457 7.908329
## ENSG00000000460 5.907786
## ENSG00000000938 -1.637671
## ENSG00000000971 12.855556
To show the effect of the transformation, 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-transformed values. For the log2 method, we need estimate size factors to account for sequencing depth (this is done automatically for the rlog method).
par( mfrow = c( 1, 2 ) )
dds <- estimateSizeFactors(dds)
plot(log2( 1 + counts(dds, normalized=TRUE)[ , 1:2] ),
pch=16, cex=0.3)
plot(assay(rld)[ , 1:2],
pch=16, cex=0.3)
Note that, in order to make it easier to see where several points are
plotted on top of each other, we set the plotting color to a
semi-transparent black and changed the points to solid circles
(pch=16
) with reduced size (cex=0.3
).
We can see how genes with low counts seem to be excessively variable on the ordinary logarithmic scale, while the rlog transform compresses differences for genes for which the data cannot provide good information anyway.
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 avoid that the distance measure is dominated by a few highly variable genes, and have a roughly equal contribution from all genes, we use it on the rlog-transformed data:
sampleDists <- dist( t( assay(rld) ) )
sampleDists
## SRR1039508 SRR1039509 SRR1039512 SRR1039513 SRR1039516 SRR1039517 SRR1039520
## SRR1039509 40.83890
## SRR1039512 37.30507 50.02512
## SRR1039513 55.69711 41.43769 43.56257
## SRR1039516 41.94152 53.53777 40.94830 57.05538
## SRR1039517 57.64886 47.54369 53.47818 46.09035 42.06359
## SRR1039520 37.01611 51.75673 34.81494 52.49906 43.16856 57.08938
## SRR1039521 55.99655 41.41332 51.85482 34.77774 58.35742 47.85851 44.73485
Note the use of the function t to transpose the data matrix. We need this because dist calculates distances between data rows and our samples constitute the columns.
We visualize the distances in a heatmap, 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 those distances in the matrix,
we manually provide the 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.
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)
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, implemented in the CRAN package
PoiClaClu. Similar to the transformations offered in
DESeq2, this measure of dissimilarity also takes the 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 can plot the heatmap as before:
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)
Another way to visualize sample-to-sample distances is a principal-components analysis (PCA). In this ordination method, the data points (i.e., here, the samples) are projected onto the 2D plane such that they spread out in the two directions which explain most of the differences in the data. The x-axis is the direction (or principal component) which separates the data points the most. The amount of the total variance which is contained in the direction is printed in the axis label.
plotPCA(rld, intgroup = c("dex", "cell"))
Here, we have used the function plotPCA which 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
ggplot2. 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.
(data <- plotPCA(rld, intgroup = c( "dex", "cell"), returnData=TRUE))
## PC1 PC2 group dex cell name
## SRR1039508 -14.324282 -4.206237 untrt : N61311 untrt N61311 SRR1039508
## SRR1039509 6.750634 -2.244202 trt : N61311 trt N61311 SRR1039509
## SRR1039512 -8.127015 -3.950826 untrt : N052611 untrt N052611 SRR1039512
## SRR1039513 14.498811 -2.940964 trt : N052611 trt N052611 SRR1039513
## SRR1039516 -11.885082 13.727739 untrt : N080611 untrt N080611 SRR1039516
## SRR1039517 8.370297 17.814698 trt : N080611 trt N080611 SRR1039517
## SRR1039520 -9.961634 -10.008929 untrt : N061011 untrt N061011 SRR1039520
## SRR1039521 14.678270 -8.191278 trt : N061011 trt N061011 SRR1039521
percentVar <- round(100 * attr(data, "percentVar"))
We can then use this data to build up the plot, specifying that the color of the points should reflect dexamethasone treatment and the shape should reflect the cell line.
library("ggplot2")
ggplot(data, aes(PC1, PC2, color=dex, shape=cell)) + geom_point(size=3) +
xlab(paste0("PC1: ",percentVar[1],"% variance")) +
ylab(paste0("PC2: ",percentVar[2],"% variance"))
From both visualizations, we see that the differences between cells are
considerable, though not stronger than the differences due to
treatment with dexamethasone. 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 by using the design formula ~ cell + dex
when setting up the
data object in the beginning.
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 the original data, but only a matrix of distances. Here we have the MDS plot for the distances calculated from the rlog transformed counts:
mds <- data.frame(cmdscale(sampleDistMatrix))
mds <- cbind(mds, as.data.frame(colData(rld)))
ggplot(mds, aes(X1,X2,color=dex,shape=cell)) + geom_point(size=3)
And here from the PoissonDistance:
mdspois <- data.frame(cmdscale(samplePoisDistMatrix))
mdspois <- cbind(mds, as.data.frame(colData(dds)))
ggplot(mdspois, aes(X1,X2,color=dex,shape=cell)) + geom_point(size=3)
Finally, we are ready to run the differential expression pipeline. With the data object prepared, the DESeq2 analysis can now be run with a single call to the function DESeq:
dds <- DESeq(dds)
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 (which control for differences in the
library size of the sequencing experiments), the estimation of
dispersion for each gene, and fitting a generalized linear model.
A DESeqDataSet is returned which contains all the fitted information within it, and the following section describes how to extract out results tables of interest from this object.
Calling results without any arguments will extract the estimated log2 fold changes and p values for the last variable in the design formula. If there are more than 2 levels for this variable, results will extract the results table for a comparison of the last level over the first level.
(res <- results(dds))
## log2 fold change (MAP): dex untrt vs trt
## Wald test p-value: dex untrt vs trt
## DataFrame with 29391 rows and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 708.6021697 0.37423027 0.09872592 3.7905979 0.0001502851 0.001224147
## ENSG00000000419 520.2979006 -0.20214241 0.10929202 -1.8495625 0.0643766347 0.189223547
## ENSG00000000457 237.1630368 -0.03624420 0.13682871 -0.2648874 0.7910961823 0.907794194
## ENSG00000000460 57.9326331 0.08520813 0.24645453 0.3457357 0.7295413492 0.875201476
## ENSG00000000938 0.3180984 0.11522626 0.14589379 0.7897955 0.4296472220 NA
## ... ... ... ... ... ... ...
## ENSG00000273485 1.2864477 -0.03490687 0.2986167 -0.1168952 0.9069431 NA
## ENSG00000273486 15.4525365 0.09662405 0.3385221 0.2854290 0.7753155 0.8990371
## ENSG00000273487 8.1632350 -0.56255486 0.3731295 -1.5076666 0.1316399 0.3177048
## ENSG00000273488 8.5844790 -0.10794133 0.3680474 -0.2932811 0.7693073 0.8960855
## ENSG00000273489 0.2758994 -0.11249629 0.1420249 -0.7920882 0.4283092 NA
As res
is a DataFrame object, it carries metadata
with information on the meaning of the columns:
mcols(res, use.names=TRUE)
## DataFrame with 6 rows and 2 columns
## type description
## <character> <character>
## baseMean intermediate mean of normalized counts for all samples
## log2FoldChange results log2 fold change (MAP): dex untrt vs trt
## lfcSE results standard error: dex untrt vs trt
## stat results Wald statistic: dex untrt vs trt
## pvalue results Wald test p-value: dex untrt vs trt
## padj results BH adjusted p-values
The first column, baseMean
, is a just the average of the normalized
count values, dividing by size factors, taken over all samples. The
remaining four columns refer to a specific contrast, namely the
comparison of the trt
level over the untrt
level for the factor
variable dex
. See the help page for results (by typing ?results
)
for information on how to obtain other contrasts.
The column log2FoldChange
is the effect size estimate. It tells us
how much the gene's expression seems to have changed due to treatment
with dexamethasone in comparison to untreated samples. This value is
reported on a logarithmic scale to base 2: for example, a log2 fold
change of 1.5 means that the gene's expression is increased by a
multiplicative factor of \( 2^{1.5} \approx 2.82 \).
Of course, this estimate has an uncertainty associated with it, which
is available in the column lfcSE
, the standard error estimate for
the log2 fold change estimate. We can also express the uncertainty of
a particular effect size estimate as the result of a statistical
test. The purpose of a test for differential expression is to test
whether the data provides sufficient evidence to conclude that this
value is really different from zero. DESeq2 performs for each gene a
hypothesis test to see whether evidence is sufficient to decide
against the null hypothesis that there is no effect of the treatment
on the gene and that the observed difference between treatment and
control was merely caused by experimental variability (i.e., the type
of variability that you can just as well expect between different
samples in the same treatment group). As usual in statistics, the
result of this test is reported as a p value, and it is found in the
column pvalue
. (Remember that a p value indicates the probability
that a fold change as strong as the observed one, or even stronger,
would be seen under the situation described by the null hypothesis.)
We can also summarize the results with the following line of code, which reports some additional information, which will be covered in later sections.
summary(res)
##
## out of 29391 with nonzero total read count
## adjusted p-value < 0.1
## LFC > 0 (up) : 2250, 7.7%
## LFC < 0 (down) : 2647, 9%
## outliers [1] : 0, 0%
## low counts [2] : 11756, 40%
## (mean count < 5.2)
## [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:
padj
in
the results table)lfcThreshold
argument of results. See the DESeq2 vignette for a demonstration
of the use of this argument.If we lower the false discovery rate threshold, we should also
tell this value to results()
, so that the function will use an
alternative threshold for the optimal independent filtering step:
res.05 <- results(dds, alpha=.05)
table(res.05$padj < .05)
##
## FALSE TRUE
## 12095 4070
Sometimes a subset of the p values in res
will be NA
(“not
available”). This is DESeq's way of reporting that all counts for
this gene were zero, and hence not test was applied. In addition, p
values can be assigned NA
if the gene was excluded from analysis
because it contained an extreme count outlier. For more information,
see the outlier detection section of the vignette.
If you use DESeq2 results in published research, please cite our
methods paper which you can find by typing citation("DESeq2")
(and the same holds for other Bioconductor packages).
This helps to support the individuals who put time into open source
software for genomic analysis.
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 in the numerator, and the name of the
level in 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 (MAP): cell N061011 vs N61311
## Wald test p-value: cell N061011 vs N61311
## DataFrame with 29391 rows and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 708.6021697 0.29054171 0.13600021 2.13633278 0.0326523 0.1981039
## ENSG00000000419 520.2979006 -0.05069310 0.14916364 -0.33984894 0.7339703 0.9238903
## ENSG00000000457 237.1630368 0.01474318 0.18161982 0.08117606 0.9353019 0.9862379
## ENSG00000000460 57.9326331 0.20241839 0.28064506 0.72126120 0.4707488 0.8108444
## ENSG00000000938 0.3180984 0.00000000 0.07169692 0.00000000 1.0000000 NA
## ... ... ... ... ... ... ...
## ENSG00000273485 1.2864477 -0.180248107 0.16456445 -1.095304057 0.2733835 NA
## ENSG00000273486 15.4525365 -0.029979349 0.30827915 -0.097247410 0.9225299 NA
## ENSG00000273487 8.1632350 -0.001914496 0.28117903 -0.006808817 0.9945674 NA
## ENSG00000273488 8.5844790 0.380608538 0.29209485 1.303030644 0.1925643 NA
## ENSG00000273489 0.2758994 0.000000000 0.06955643 0.000000000 1.0000000 NA
If results for an interaction term are desired, the name
argument of results should be used. Please see the
help for the results function for more details.
Novices in high-throughput biology often assume that thresholding these p values at a low value, say 0.05, as is often done in other settings, would be appropriate – but it is not. We briefly explain why:
There are 5722 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] 5722
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 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 / 5722 = 26% false positives.
DESeq2 uses the Benjamini-Hochberg (BH) adjustment as described in
the base R p.adjust function; in brief, this method calculates for
each gene an adjusted p value which answers the following question:
if one called significant all genes with a p value less than or
equal to this gene's p value threshold, what would be the fraction
of false positives (the false discovery rate, FDR) among them (in
the sense of the calculation outlined above)? These values, called the
BH-adjusted p values, are given in the column padj
of the res
object.
Hence, if we consider a fraction of 10% false positives acceptable, we can consider all genes with an adjusted p value below 10% = 0.1 as significant. How many such genes are there?
sum(res$padj < 0.1, na.rm=TRUE)
## [1] 4897
We subset the results table to these genes and then sort it by the log2 fold change estimate to get the significant genes with the strongest down-regulation.
resSig <- subset(res, padj < 0.1)
head(resSig[ order( resSig$log2FoldChange ), ])
## log2 fold change (MAP): dex untrt vs trt
## Wald test p-value: dex untrt vs trt
## DataFrame with 6 rows and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000179593 67.24305 -4.880507 0.3308118 -14.75312 2.937579e-49 9.418948e-47
## ENSG00000109906 385.07103 -4.860888 0.3321621 -14.63408 1.702536e-48 5.176591e-46
## ENSG00000152583 997.43977 -4.315374 0.1723804 -25.03402 2.606482e-138 4.596531e-134
## ENSG00000250978 56.31819 -4.090157 0.3288245 -12.43872 1.610577e-35 2.679483e-33
## ENSG00000163884 561.10717 -4.078073 0.2103212 -19.38974 9.420848e-84 1.038354e-80
## ENSG00000168309 159.52692 -3.991152 0.2547745 -15.66543 2.606752e-55 1.178720e-52
…and with the strongest upregulation. The order function gives
the indices in increasing order, so a simple way to ask for decreasing
order is to add a -
sign. Alternatively, you can use the argument
decreasing=TRUE
.
head(resSig[ order( -resSig$log2FoldChange ), ])
## log2 fold change (MAP): dex untrt vs trt
## Wald test p-value: dex untrt vs trt
## DataFrame with 6 rows and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000162692 508.17023 3.452453 0.1763753 19.574472 2.552673e-85 3.462799e-82
## ENSG00000146006 46.80760 2.856273 0.3366876 8.483453 2.186078e-17 1.073858e-15
## ENSG00000105989 333.21469 2.850961 0.1754634 16.248175 2.301133e-59 1.193544e-56
## ENSG00000214814 243.27698 2.759538 0.2224910 12.402921 2.519693e-35 4.114332e-33
## ENSG00000267339 26.23357 2.743928 0.3511985 7.813039 5.582532e-15 2.182881e-13
## ENSG00000013293 244.49733 2.646115 0.1981218 13.356003 1.092766e-40 2.240806e-38
A quick way to visualize the counts for a particular gene is to use the plotCounts function, which takes as arguments the DESeqDataSet, a gene name, and the group over which to plot the counts.
topGene <- rownames(res)[which.min(res$padj)]
plotCounts(dds, gene=topGene, intgroup=c("dex"))
We can also make more customizable plots using the ggplot function from the ggplot2 package:
data <- plotCounts(dds, gene=topGene, intgroup=c("dex","cell"), returnData=TRUE)
ggplot(data, aes(x=dex, y=count, color=cell)) +
scale_y_log10() +
geom_point(position=position_jitter(width=.1,height=0), size=3)
Here we use a more structural arrangement instead of random jitter, and color by the treatment.
ggplot(data, aes(x=dex, y=count, fill=dex)) +
scale_y_log10() +
geom_dotplot(binaxis="y", stackdir="center")
Note that the DESeq test actually takes into account the cell line effect, so a more detailed plot would also show the cell lines.
ggplot(data, aes(x=dex, y=count, color=cell, group=cell)) +
scale_y_log10() + geom_point(size=3) + geom_line()
An “MA-plot” provides a useful overview for an experiment with a two-group comparison. The log2 fold change for a particular comparison is plotted on the y-axis and the average of the counts normalized by size factor is shown on the x-axis (“M” for minus, because a log ratio is equal to log minus log, and “A” for average).
plotMA(res, ylim=c(-5,5))
Each gene is represented with a dot. Genes with an adjusted \(p\) value below a threshold (here 0.1, the default) are shown in red. The DESeq2 package incorporates a prior on log2 fold changes, resulting in moderated log2 fold changes from genes with low counts and highly variable counts, as can be seen by the narrowing of spread of points on the left side of the plot. This plot demonstrates that only genes with a large average normalized count contain sufficient information to yield a significant call.
We can label individual points on the MA plot as well. Here we use the
with R function to plot a circle and text for a selected row of the
results object. Within the with function, only the baseMean
and
log2FoldChange
values for the selected rows of res
are used.
plotMA(res, ylim=c(-5,5))
with(res[topGene, ], {
points(baseMean, log2FoldChange, col="dodgerblue", cex=2, lwd=2)
text(baseMean, log2FoldChange, topGene, pos=2, col="dodgerblue")
})
Whether a gene is called significant depends not only on its LFC but also on its within-group variability, which DESeq2 quantifies as the dispersion. For strongly expressed genes, the dispersion can be understood as a squared coefficient of variation: a dispersion value of 0.01 means that the gene's expression tends to differ by typically \(\sqrt{0.01} = 10\%\) between samples of the same treatment group. For weak genes, the Poisson noise is an additional source of noise.
Another useful diagnostic plot is the histogram of the p values.
hist(res$pvalue, breaks=20, col="grey50", border="white")
This plot becomes a bit smoother by excluding genes with very small counts:
hist(res$pvalue[res$baseMean > 1], breaks=20, col="grey50", border="white")
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 signal, one usually carries it out only for a subset of most highly variable genes. Here, for demonstration, let us select the 35 genes with the highest variance across samples. We will work with the rlog transformed counts:
library("genefilter")
topVarGenes <- head(order(-rowVars(assay(rld))),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. We provide the column side colors to help identify the treated samples (in blue) from the untreated samples (in grey).
mat <- assay(rld)[ topVarGenes, ]
mat <- mat - rowMeans(mat)
df <- as.data.frame(colData(rld)[,c("cell","dex")])
pheatmap(mat, annotation_col=df)
We can now see blocks of genes which covary across patients. Note that a set of genes at the top of the heatmap are separating the N061011 cell line from the others. At the bottom of the heatmap, we see a set of genes for which the treated samples have higher gene expression.
The MA plot highlights an important property of RNA-Seq data. For weakly expressed genes, we have no chance of seeing differential expression, because the low read counts suffer from so high Poisson noise that any biological effect is drowned in the uncertainties from the read counting. We can also show this by examining the ratio of small p values (say, less than, 0.01) for genes binned by mean normalized count:
# create bins using the quantile function
qs <- c(0, quantile(res$baseMean[res$baseMean > 0], 0:7/7))
# cut the genes into the bins
bins <- cut(res$baseMean, qs)
# rename the levels of the bins using the middle point
levels(bins) <- paste0("~",round(.5*qs[-1] + .5*qs[-length(qs)]))
# calculate the ratio of $p$ values less than .01 for each bin
ratios <- tapply(res$pvalue, bins, function(p) mean(p < .01, na.rm=TRUE))
# plot these ratios
barplot(ratios, xlab="mean normalized count", ylab="ratio of small p values")
At first sight, there may seem to be little benefit in filtering out these genes. After all, the test found them to be non-significant anyway. However, these genes have an influence on the multiple testing adjustment, whose performance improves if such genes are removed. By removing the weakly-expressed genes from the input to the FDR procedure, we can find more genes to be significant among those which we keep, and so improved the power of our test. This approach is known as independent filtering.
The DESeq2 software automatically performs independent filtering
which maximizes the number of genes which will have 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. We can observe how the number of
rejections changes for various cutoffs based on mean normalized
count. The following optimal threshold and table of possible values is
stored as an attribute of the results object.
attr(res,"filterThreshold")
## 40%
## 5.203248
plot(attr(res,"filterNumRej"), type="b",
xlab="quantiles of 'baseMean'",
ylab="number of rejections")
The term independent highlights an important caveat. Such filtering is permissible only if the filter criterion is independent of the actual test statistic. Otherwise, the filtering would invalidate the test and consequently the assumptions of the BH procedure. This is why we filtered on the average over all samples: this filter is blind to the assignment of samples to the treatment and control group and hence independent. The independent filtering software used inside DESeq2 comes from the genefilter package, which contains a reference to a paper describing the statistical foundation for independent filtering.
Our result table only uses Ensembl gene IDs, but gene names may be more informative. Bioconductor's annotation packages help with mapping various ID schemes to each other.
We load the AnnotationDbi package and the annotation package org.Hs.eg.db:
library("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] "ENTREZID" "PFAM" "IPI" "PROSITE" "ACCNUM" "ALIAS"
## [7] "CHR" "CHRLOC" "CHRLOCEND" "ENZYME" "MAP" "PATH"
## [13] "PMID" "REFSEQ" "SYMBOL" "UNIGENE" "ENSEMBL" "ENSEMBLPROT"
## [19] "ENSEMBLTRANS" "GENENAME" "UNIPROT" "GO" "EVIDENCE" "ONTOLOGY"
## [25] "GOALL" "EVIDENCEALL" "ONTOLOGYALL" "OMIM" "UCSCKG"
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] "ENTREZID" "PFAM" "IPI" "PROSITE" "ACCNUM" "ALIAS"
## [7] "CHR" "CHRLOC" "CHRLOCEND" "ENZYME" "MAP" "PATH"
## [13] "PMID" "REFSEQ" "SYMBOL" "UNIGENE" "ENSEMBL" "ENSEMBLPROT"
## [19] "ENSEMBLTRANS" "GENENAME" "UNIPROT" "GO" "EVIDENCE" "ONTOLOGY"
## [25] "GOALL" "EVIDENCEALL" "ONTOLOGYALL" "OMIM" "UCSCKG"
We can use the mapIds function to add invidual columns to our results table. 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 untrt vs trt
## Wald test p-value: dex untrt vs trt
## DataFrame with 6 rows and 8 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000152583 997.4398 -4.315374 0.1723804 -25.03402 2.606482e-138 4.596531e-134
## ENSG00000165995 495.0929 -3.188413 0.1277306 -24.96201 1.581652e-137 1.394622e-133
## ENSG00000101347 12703.3871 -3.617791 0.1499256 -24.13057 1.194487e-128 6.471353e-125
## ENSG00000120129 3409.0294 -2.871106 0.1190242 -24.12204 1.467843e-128 6.471353e-125
## ENSG00000189221 2341.7673 -3.230291 0.1373498 -23.51871 2.625339e-122 9.259569e-119
## ENSG00000211445 12285.6151 -3.552498 0.1589749 -22.34628 1.312364e-110 3.857258e-107
## symbol entrez
## <character> <character>
## ENSG00000152583 SPARCL1 8404
## ENSG00000165995 CACNB2 783
## ENSG00000101347 SAMHD1 25939
## ENSG00000120129 DUSP1 1843
## ENSG00000189221 MAOA 4128
## ENSG00000211445 GPX3 2878
You can easily save the results table in a CSV file, which you can then load with a spreadsheet program such as Excel. The call to as.data.frame is necessary to convert the DataFrame object (IRanges package) to a data.frame object which can be processed by write.csv.
write.csv(as.data.frame(resOrdered), file="results.csv")
Another more sophisticated package for exporting results from various Bioconductor analysis packages is the ReportingTools package. 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 more details.
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 location of the genes. We
can therefore easily plot our differential expression results in
genomic space. While the results function by default outputs a
DataFrame, using the format
argument, we can ask for GRanges or
GRangesList output.
(resGR <- results(dds, format="GRanges"))
## GRanges object with 29391 ranges and 6 metadata columns:
## seqnames ranges strand | baseMean log2FoldChange
## <Rle> <IRanges> <Rle> | <numeric> <numeric>
## ENSG00000000003 X [ 99883667, 99894988] - | 708.602169691234 0.374230273808482
## ENSG00000000419 20 [ 49551404, 49575092] - | 520.297900552084 -0.202142413648455
## ENSG00000000457 1 [169818772, 169863408] - | 237.163036796015 -0.0362442040444828
## ENSG00000000460 1 [169631245, 169823221] + | 57.9326331250967 0.0852081315366819
## ENSG00000000938 1 [ 27938575, 27961788] - | 0.318098378392895 0.11522625736902
## ... ... ... ... ... ... ...
## ENSG00000273485 10 [105209953, 105210609] + | 1.28644765243289 -0.0349068685144913
## ENSG00000273486 3 [136556180, 136557863] - | 15.4525365439045 0.0966240513958813
## ENSG00000273487 1 [ 92654794, 92656264] + | 8.1632349843654 -0.562554859380289
## ENSG00000273488 3 [100080031, 100080481] + | 8.58447903624707 -0.107941329357511
## ENSG00000273489 7 [131178723, 131182453] - | 0.275899382507797 -0.112496288893145
## lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 0.0987259223705911 3.7905978979231 0.000150285083105386 0.00122414662381685
## ENSG00000000419 0.10929201676454 -1.84956248070665 0.0643766346624258 0.189223547475715
## ENSG00000000457 0.136828705328857 -0.264887429559264 0.791096182308222 0.907794194104991
## ENSG00000000460 0.246454533269651 0.345735703889261 0.729541349240061 0.875201475772005
## ENSG00000000938 0.14589379275318 0.789795475150598 0.429647222034592 <NA>
## ... ... ... ... ...
## ENSG00000273485 0.298616693212256 -0.116895234954864 0.906943074221519 <NA>
## ENSG00000273486 0.338522146535403 0.285429040270416 0.77531546155328 0.89903706080337
## ENSG00000273487 0.373129485902887 -1.50766658930488 0.131639880655068 0.317704844033409
## ENSG00000273488 0.368047392946519 -0.293281059521583 0.769307329767728 0.896085519184537
## ENSG00000273489 0.142024947940281 -0.792088224812784 0.428309235519245 <NA>
## -------
## seqinfo: 722 sequences (1 circular) from an unspecified genome
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 which have a fold change (exclude genes with no counts). We add the gene symbol as a name, if the symbol exists or is not duplicated in our subset.
window <- resGR[topGene] + 1e6
strand(window) <- "*"
hasLFC <- !is.na(resGR$log2FoldChange)
resGRsub <- resGR[resGR %over% window & hasLFC]
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 false discovery rate.
sig <- factor(ifelse(is.na(resGRsub$padj) | resGRsub$padj > .1,"notsig","sig"))
We can then plot the results using Gviz functions. We create an axis track specifying our location in the genome, a track which will show the genes and their names, colored by significance, and a data track which 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=sig)
d <- DataTrack(resGRsub, data="log2FoldChange", baseline=0,
type="h", name="log2 fold change", strand="+")
plotTracks(list(g,d,a), groupAnnotation="group", notsig="lightblue", sig="pink")
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 which might be affecting many or all of the genes in the dataset. We can use statistical methods from the sva package to detect such groupings of the samples, and then we can add these to the DESeqDataSet's design, in order to account for them. The SVA software uses the term surrogate variables for the estimated variables which we want to account for in our analysis.
library("sva")
Below we get 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
\Robject{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 help at ?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:
par(mfrow=c(2,1),mar=c(5,5,1,1))
stripchart(svseq$sv[,1] ~ dds$cell,vertical=TRUE)
abline(h=0)
stripchart(svseq$sv[,2] ~ dds$cell,vertical=TRUE)
abline(h=0)
Finally, in order to use SVA to remove any effect on the counts from our surrogate variables, we simply add these two as columns to the DESeqDataSet and add them to the design.
ddssva <- dds
ddssva$SV1 <- svseq$sv[,1]
ddssva$SV2 <- svseq$sv[,2]
design(ddssva) <- ~ SV1 + SV2 + dex
ddssva <- DESeq(ddssva)
head(results(ddssva), 4)
## log2 fold change (MAP): dex untrt vs trt
## Wald test p-value: dex untrt vs trt
## DataFrame with 4 rows and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## ENSG00000000003 708.60217 0.37499573 0.1062538 3.5292447 0.0004167476 0.003254457
## ENSG00000000419 520.29790 -0.19957419 0.1163350 -1.7155125 0.0862513038 0.239356623
## ENSG00000000457 237.16304 -0.03449066 0.1415061 -0.2437397 0.8074324328 0.915682986
## ENSG00000000460 57.93263 0.10130192 0.2718768 0.3726023 0.7094444358 0.866426237
DESeq2 can be used to analyze time series experiments, for example
to find those genes which react in a condition specific manner over
time. Here we demonstrate a basic time series analysis with the
fission data package,
which contains gene counts for an RNA-Seq time course of fission
yeast. The yeast were exposed to oxidative stress, and half of the
samples contain a deletion of the gene atf21.
We use a design which 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 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 which 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$pvalue),],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 padj symbol
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric> <character>
## SPBC2F12.09c 174.6712 -2.65763737 0.7498270 99.23199 7.671942e-20 5.186233e-16 atf21
## SPAC1002.18 444.5050 -0.05118463 0.2030554 57.72116 3.590886e-11 1.213719e-07 urg3
## SPAC1002.19 336.3732 -0.39267927 0.5749887 43.26296 3.268243e-08 7.364441e-05 urg1
## SPAC1002.17c 261.7731 -1.13882844 0.6072772 39.13718 2.228530e-07 3.766216e-04 urg2
This is just one of the tests which 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.
We can plot the counts for the groups over time using ggplot2, for the gene with the smallest p value, testing for condition-dependent time profile and accounting for differences at time 0. 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.
data <- plotCounts(ddsTC, which.min(resTC$pvalue),
intgroup=c("minute","strain"), returnData=TRUE)
ggplot(data, aes(x=minute, y=count, color=strain, group=strain)) +
geom_point() + stat_smooth(se=FALSE,method="loess") + scale_y_log10()
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" "minute_30_vs_0"
## [5] "minute_60_vs_0" "minute_120_vs_0" "minute_180_vs_0" "strainmut.minute15"
## [9] "strainmut.minute30" "strainmut.minute60" "strainmut.minute120" "strainmut.minute180"
res30 <- results(ddsTC, name="strainmut.minute30", test="Wald")
res30[which.min(resTC$pvalue),]
## log2 fold change (MLE): strainmut.minute30
## Wald test p-value: strainmut.minute30
## DataFrame with 1 row and 6 columns
## baseMean log2FoldChange lfcSE stat pvalue padj
## <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
## SPBC2F12.09c 174.6712 -2.601034 0.6314737 -4.11899 3.805364e-05 0.2572426
We can further more cluster significant genes by their profiles.
betas <- coef(ddsTC)
colnames(betas)
## [1] "Intercept" "strain_mut_vs_wt" "minute_15_vs_0" "minute_30_vs_0"
## [5] "minute_60_vs_0" "minute_120_vs_0" "minute_180_vs_0" "strainmut.minute15"
## [9] "strainmut.minute30" "strainmut.minute60" "strainmut.minute120" "strainmut.minute180"
library("pheatmap")
topGenes <- head(order(resTC$pvalue),20)
mat <- betas[topGenes, -c(1,2)]
thr <- 3 # threshold for plotting
mat[mat < -thr] <- -thr
mat[mat > thr] <- thr
pheatmap(mat, breaks=seq(from=-thr, to=thr, length=101),
cluster_col=FALSE)
As last part of this document, we call the function sessionInfo, which reports the version numbers of R and all the packages used in this session. It is good practice to always keep such a record as it will help to trace down what has happened in case that an R script ceases to work because the functions have been changed in a newer version of a package. The session information should also always be included in any emails to the Bioconductor support site along with all code used in the analysis.
sessionInfo()
## R version 3.2.1 (2015-06-18)
## Platform: x86_64-unknown-linux-gnu (64-bit)
## Running under: Ubuntu precise (12.04.4 LTS)
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C LC_TIME=en_US.UTF-8
## [4] LC_COLLATE=C LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C LC_ADDRESS=C
## [10] LC_TELEPHONE=C LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] grid stats4 parallel stats graphics grDevices utils datasets methods
## [10] base
##
## other attached packages:
## [1] fission_0.102.0 sva_3.14.0 mgcv_1.8-6 nlme_3.1-120
## [5] Gviz_1.12.1 org.Hs.eg.db_3.1.2 RSQLite_1.0.0 DBI_0.3.1
## [9] genefilter_1.50.0 ggplot2_1.0.1 PoiClaClu_1.0.2 RColorBrewer_1.1-2
## [13] pheatmap_1.0.7 DESeq2_1.8.1 RcppArmadillo_0.5.300.4 Rcpp_0.12.0
## [17] BiocParallel_1.2.20 GenomicAlignments_1.4.1 GenomicFeatures_1.20.2 AnnotationDbi_1.30.1
## [21] Biobase_2.28.0 Rsamtools_1.20.4 Biostrings_2.36.3 XVector_0.8.0
## [25] airway_0.102.0 GenomicRanges_1.20.5 GenomeInfoDb_1.4.2 IRanges_2.2.7
## [29] S4Vectors_0.6.3 BiocGenerics_0.14.0 rmarkdown_0.7 knitr_1.10.5
## [33] BiocStyle_1.6.0
##
## loaded via a namespace (and not attached):
## [1] splines_3.2.1 Formula_1.2-1 latticeExtra_0.6-26
## [4] BSgenome_1.36.3 lattice_0.20-31 biovizBase_1.16.0
## [7] digest_0.6.8 colorspace_1.2-6 htmltools_0.2.6
## [10] Matrix_1.2-1 plyr_1.8.3 XML_3.98-1.3
## [13] biomaRt_2.24.0 zlibbioc_1.14.0 xtable_1.7-4
## [16] scales_0.2.5 annotate_1.46.1 nnet_7.3-9
## [19] proto_0.3-10 survival_2.38-1 magrittr_1.5
## [22] evaluate_0.7.2 MASS_7.3-40 foreign_0.8-63
## [25] tools_3.2.1 formatR_1.2 matrixStats_0.14.2
## [28] stringr_1.0.0 munsell_0.4.2 locfit_1.5-9.1
## [31] cluster_2.0.1 lambda.r_1.1.7 futile.logger_1.4.1
## [34] RCurl_1.95-4.7 dichromat_2.0-0 VariantAnnotation_1.14.11
## [37] bitops_1.0-6 labeling_0.3 gtable_0.1.2
## [40] reshape2_1.4.1 gridExtra_2.0.0 rtracklayer_1.28.8
## [43] Hmisc_3.16-0 futile.options_1.0.0 stringi_0.5-5
## [46] geneplotter_1.46.0 rpart_4.1-9 acepack_1.3-3.3