/
Deseq2_edgeR.R
439 lines (415 loc) · 13.7 KB
/
Deseq2_edgeR.R
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
################################################################################
#' Plot edgeR results for a phyloseq or a edgeR object.
#'
#' `r lifecycle::badge("maturing")`
#'
#' @inheritParams clean_pq
#' @param contrast (required):This argument specifies what comparison
#' to extract from the object to build a results table.
#' See \code{\link[DESeq2]{results}} man page for more details.
#' @param pval (default: 0.05): the significance cutoff used for optimizing
#' the independent filtering. If the adjusted p-value cutoff (FDR) will be a
#' value other than 0.05, pval should be set to that value.
#' @param taxolev taxonomic level of interest
#' @param color_tax taxonomic level used for
#' color assignation
#' @param verbose (logical): whether the function print some
#' information during the computation
#' @param ... Additional arguments passed on to \code{\link[edgeR]{exactTest}}
#' or \code{\link[ggplot2]{ggplot}}
#'
#' @export
#' @examplesIf tolower(Sys.info()[["sysname"]]) != "windows"
#' data("GlobalPatterns", package = "phyloseq")
#' GP_archae <- subset_taxa(GlobalPatterns, GlobalPatterns@tax_table[, 1] == "Archaea")
#' \donttest{
#' if (requireNamespace("edgeR")) {
#' plot_edgeR_pq(GP_archae, c("SampleType", "Soil", "Feces"),
#' color_tax = "Kingdom"
#' )
#'
#' plot_edgeR_pq(GP_archae, c("SampleType", "Soil", "Feces"),
#' taxolev = "Class", color_tax = "Kingdom"
#' )
#' }
#' }
#' @author Adrien Taudière
#'
#' @return A \code{\link{ggplot}}2 plot representing edgeR results
#'
#' @seealso \code{\link[edgeR]{exactTest}}
#' @seealso \code{\link{plot_deseq2_pq}}
plot_edgeR_pq <-
function(physeq,
contrast = NULL,
pval = 0.05,
taxolev = "Genus",
color_tax = "Phylum",
verbose = TRUE,
...) {
if (!inherits(physeq, "phyloseq")) {
stop("data must be an object of class 'phyloseq'")
}
if (verbose) {
message("Conversion to edgeR format")
}
data_edger <- phyloseq_to_edgeR(physeq, group = contrast[1])
if (verbose) {
message("Perform edgeR binary test")
}
et <-
edgeR::exactTest(data_edger, pair = c(contrast[2], contrast[3]), ...)
tt <-
edgeR::topTags(
et,
n = nrow(et$table),
adjust.method = "BH",
sort.by = "PValue"
)
res <- tt@.Data[[1]]
sigtab <- res[(res$FDR < pval), ]
sigtab <- cbind(methods::as(sigtab, "data.frame"))
sigtabgen <- subset(sigtab, !is.na(taxolev))
d <-
tapply(sigtabgen$logFC, sigtabgen[, color_tax], function(x) {
max(x)
})
d <- sort(d, TRUE)
sigtabgen$col_tax <-
factor(as.character(sigtabgen[, color_tax]), levels = names(d))
d <-
tapply(sigtabgen$logFC, sigtabgen[, taxolev], function(x) {
max(x)
})
d <- sort(d, TRUE)
sigtabgen$tax <-
factor(as.character(sigtabgen[, taxolev]), levels = names(d))
p <-
ggplot(sigtabgen, aes(x = tax, y = logFC, color = col_tax), ...) +
geom_point(size = 6) +
theme(axis.text.x = element_text(
angle = -90,
hjust = 0,
vjust = 0.5
)) +
labs(
title = paste(
"Change in abundance for ",
contrast[1],
" (",
contrast[2],
" vs ",
contrast[3],
")",
sep = ""
)
)
return(p)
}
################################################################################
################################################################################
# Plot the result of a DESeq2 test
################################################################################
#' Plot DESeq2 results for a phyloseq or a DESeq2 object.
#' @description
#' `r lifecycle::badge("experimental")`
#'
#' @param data (required) a \code{\link{phyloseq-class}} or a
#' \code{\link[DESeq2]{DESeqDataSet-class}} object.
#' @param tax_table Required if data is a
#' \code{\link[DESeq2]{DESeqDataSet-class}} object.
#' The taxonomic table used to find the \code{taxa} and \code{color_taxa}
#' arguments. If data is a \code{\link{phyloseq-class}} object, data@tax_table
#' is used.
#' @param contrast (required) contrast specifies what comparison to extract
#' from the object to build a results table. See \code{\link[DESeq2]{results}}
#' man page for more details.
#' @param pval (default: 0.05) the significance cutoff used for optimizing
#' the independent filtering. If the adjusted p-value cutoff (FDR) will be a
#' value other than 0.05, pval should be set to that value.
#' @param taxolev taxonomic level of interest
#' @param select_taxa Either the name of the taxa (in the form of [DESeq2::results()])
#' or a logical vector (length of the results from [DESeq2::results()]) to select taxa
#' to plot.
#' @param color_tax taxonomic level used for color or a
#' color vector.
#' @param tax_depth Taxonomic depth to test for differential
#' distribution among contrast. If Null the analysis is done at the OTU
#' (i.e. Species) level. If not Null, data need to be a column name in
#' the `tax_table` slot of the \code{\link{phyloseq-class}} object.
#' @param verbose whether the function print some information during
#' the computation
#' @param jitter_width width for the jitter positioning
#' @param ... Additional arguments passed on to \code{\link[DESeq2]{DESeq}}
#' or \code{\link[ggplot2]{ggplot}}
#'
#' @export
#' @details
#' Please cite `DESeq2` package if you use chis function.
#' @examples
#' \donttest{
#'
#' data("GlobalPatterns", package = "phyloseq")
#' GP <- subset_taxa(GlobalPatterns, GlobalPatterns@tax_table[, 1] == "Archaea")
#' GP <- subset_samples(GP, SampleType %in% c("Soil", "Skin"))
#' if (requireNamespace("DESeq2")) {
#' res <- DESeq2::DESeq(phyloseq_to_deseq2(GP, ~SampleType),
#' test = "Wald", fitType = "local"
#' )
#' plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
#' tax_table = GP@tax_table, color_tax = "Kingdom"
#' )
#' plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
#' tax_table = GP@tax_table, color_tax = "Kingdom",
#' pval = 0.7
#' )
#' plot_deseq2_pq(res, c("SampleType", "Soil", "Skin"),
#' tax_table = GP@tax_table, color_tax = "Class",
#' select_taxa = c("522457", "271582")
#' )
#' }
#' }
#' @author Adrien Taudière
#'
#' @return A \code{\link{ggplot}}2 plot representing DESeq2 results
#'
#' @seealso \code{\link[DESeq2]{DESeq}}
#' @seealso \code{\link[DESeq2]{results}}
#' @seealso \code{\link{plot_edgeR_pq}}
plot_deseq2_pq <-
function(data,
contrast = NULL,
tax_table = NULL,
pval = 0.05,
taxolev = "Genus",
select_taxa = NULL,
color_tax = "Phylum",
tax_depth = NULL,
verbose = TRUE,
jitter_width = 0.1,
...) {
if (!inherits(data, "phyloseq")) {
if (!inherits(data, "DESeqDataSet")) {
stop("data must be an object of class 'phyloseq' or 'DESeqDataSet'")
}
} else {
# Calculate new dataset given the Taxa depth if tax_depth is not null
if (!is.null(tax_depth)) {
data_tax <- data
data_tax@otu_table <-
otu_table(
apply(data@otu_table, 2, function(x) {
tapply(
x,
data@tax_table[, tax_depth], sum
)
}),
taxa_are_rows = TRUE
)
data_tax@tax_table <-
tax_table(apply(
data@tax_table[, 1:match(
tax_depth,
colnames(data@tax_table)
)],
2, function(x) {
tapply(
x, data@tax_table[, tax_depth],
function(xx) {
xx[1]
}
)
}
))
data_tax@refseq <- NULL
data <- data_tax
if (is.na(match(taxolev, colnames(data@tax_table)))) {
taxolev <- tax_depth
}
}
if (is.null(tax_table) && inherits(data, "phyloseq")) {
tax_table <- data@tax_table
}
if (verbose) {
message("Conversion to Deseq2 format.")
}
data_deseq2 <-
phyloseq_to_deseq2(data, stats::as.formula(paste("~", contrast[1])))
if (verbose) {
message("Calculation of Deseq2 results.")
}
data <-
DESeq2::DESeq(
data_deseq2,
test = "Wald",
fitType = "parametric",
quiet = !verbose,
...
)
}
# Calcul deseq2 results
res <- DESeq2::results(data, contrast = contrast)
if (!is.null(select_taxa[1])) {
res <- res[select_taxa, ]
}
d <- res[which(res$padj < pval), ]
if (dim(d)[1] == 0) {
message("None taxa present significant distribution pattern through
contrast.")
return("None taxa present significant distribution pattern through
contrast.")
}
d <-
cbind(methods::as(d, "data.frame"), methods::as(tax_table[rownames(d), ], "matrix"))
# Compute colors
are_colors <- function(x) {
sapply(x, function(xx) {
tryCatch(
is.matrix(grDevices::col2rgb(xx)),
error = function(e) {
FALSE
}
)
})
}
if (!sum(are_colors(color_tax)) > 0) {
x <- tapply(d$log2FoldChange, d[, color_tax], function(x) {
max(x)
})
x <- sort(x, TRUE)
d$col_tax <-
factor(as.character(d[, color_tax]), levels = names(x))
} else {
d$col_tax <- rep(color_tax, length = dim(d)[1])
}
# Compute log2FoldChange values
x <- tapply(d$log2FoldChange, d[, taxolev], function(x) {
max(x)
})
x <- sort(x, TRUE)
d$tax <- factor(as.character(d[, taxolev]), levels = names(x))
if (!sum(are_colors(color_tax)) > 0) {
p <-
ggplot(d, aes(x = tax, y = log2FoldChange, color = col_tax), ...) +
geom_point(
size = 6,
position = position_jitter(width = jitter_width, height = 0)
) +
theme(axis.text.x = element_text(
angle = -90,
hjust = 0,
vjust = 0.5
)) +
labs(
title = paste(
"Change in abundance for ",
contrast[1],
" (top:",
contrast[2],
" vs down:",
contrast[3],
")",
sep = ""
)
)
} else {
p <-
ggplot(d, aes(x = tax, y = log2FoldChange), ...) +
geom_point(
size = 6, color = d$col_tax,
position = position_jitter(width = jitter_width, height = 0)
) +
theme(axis.text.x = element_text(
angle = -90,
hjust = 0,
vjust = 0.5
)) +
labs(
title = paste(
"Change in abundance for ",
contrast[1],
" (top:",
contrast[2],
" vs down:",
contrast[3],
")",
sep = ""
)
)
}
return(p)
}
################################################################################
################################################################################
#' Convert phyloseq OTU count data into DGEList for edgeR package
#'
#'
#' @inheritParams clean_pq
#'
#' @param group (required) A character vector or factor giving the experimental
#' group/condition for each sample/library. Alternatively, you may provide
#' the name of a sample variable. This name should be among the output of
#' \code{sample_variables(physeq)}, in which case
#' \code{get_variable(physeq, group)} would return either a character vector or
#' factor.
#' This is passed on to \code{\link[edgeR]{DGEList}},
#' and you may find further details or examples in its documentation.
#'
#' @param method The label of the edgeR-implemented normalization
#' to use.
#' See \code{\link[edgeR]{calcNormFactors}} for supported options and details.
#' The default option is \code{"RLE"}, which is a scaling factor method
#' proposed by Anders and Huber (2010).
#' At time of writing, the \link[edgeR]{edgeR} package supported
#' the following options to the \code{method} argument:
#'
#' \code{c("TMM", "RLE", "upperquartile", "none")}.
#'
#' @return A DGEList object. See [edgeR::estimateTagwiseDisp()] for more details.
#'
#' @param ... Additional arguments passed on to \code{\link[edgeR]{DGEList}}
#' @export
#'
phyloseq_to_edgeR <- function(physeq, group, method = "RLE", ...) {
verify_pq(physeq)
# Enforce orientation.
if (!taxa_are_rows(physeq)) {
otu_table(physeq) <- otu_table(t(as.matrix(unclass(physeq@otu_table))),
taxa_are_rows = TRUE
)
}
x <- methods::as(otu_table(physeq), "matrix")
# Add one to protect against overflow, log(0) issues.
x <- x + 1
# Check `group` argument
if (identical(all.equal(length(group), 1), TRUE) &&
nsamples(physeq) > 1) {
# Assume that group was a sample variable name (must be categorical)
group <- get_variable(physeq, group)
}
# Define gene annotations (`genes`) as tax_table
taxonomy <- tax_table(physeq, errorIfNULL = FALSE)
if (!is.null(taxonomy)) {
taxonomy <- data.frame(methods::as(taxonomy, "matrix"))
}
# Now turn into a DGEList
y <- edgeR::DGEList(
counts = x,
group = group,
genes = taxonomy,
remove.zeros = TRUE,
...
)
# Calculate the normalization factors
z <- edgeR::calcNormFactors(y, method = method)
# Check for division by zero inside `calcNormFactors`
if (!all(is.finite(z$samples$norm.factors))) {
stop(
"Something wrong with edgeR::calcNormFactors on this data,
non-finite $norm.factors, consider changing `method` argument"
)
}
# Estimate dispersions
return(edgeR::estimateTagwiseDisp(edgeR::estimateCommonDisp(z)))
}
################################################################################