-
Notifications
You must be signed in to change notification settings - Fork 1
/
AllClasses.R
310 lines (301 loc) · 9.07 KB
/
AllClasses.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
## DESeqAnalysis ===============================================================
#' DESeq2 differential expression analysis
#'
#' Class containing all elements generated during differential expression
#' analysis with DESeq2. This class is essentially a `list` with validity checks
#' to ensure `DESeqTransform` and `DESeqResults` correspond to the
#' `DESeqDataSet`.
#'
#' @section DESeqDataSet:
#'
#' We recommend generating the `DESeqDataSet` by coercion from `bcbioRNASeq`
#' object using `as(dds, "bcbioRNASeq")`. Don't use the `DESeq2::DESeqDataSet()`
#' or `DESeq2::DESeqDataSetFromMatrix()` constructors to generate the
#' `DESeqDataSet` object.
#'
#' @section DESeqTransform:
#'
#' Object containing variance-stabilized counts. We recommend slotting the
#' return from either `DESeq2::varianceStabilizingTransformation()` or
#' `DESeq2::rlog()`.
#'
#' @section DESeqResults:
#'
#' Don't modify any of the `DESeqResults` objects manually. This includes
#' rearranging the rows or dropping genes without adjusted P values. We'll take
#' care of this automatically in supported methods.
#'
#' @author Michael Steinbaugh
#' @export
#' @note Updated 2023-12-18.
#'
#' @slot data `DESeqDataSet`.
#'
#' @slot transform `DESeqTransform`.
#'
#' @slot results `list`.
#' One or more unshrunken `DESeqResults`.
#'
#' @slot lfcShrink `list`.
#' *Optional*. One or more shrunken `DESeqResults`. If set, must correspond to
#' those defined in `results`. Otherwise, can set as empty list (`list()`).
#'
#' @seealso
#' - `help(topic = "Annotated-class", package = "S4Vectors")`.
#' - `showClass("Annotated")`.
#'
#' @return `DESeqAnalysis`.
setClass(
Class = "DESeqAnalysis",
contains = "Annotated",
slots = c(
data = "DESeqDataSet",
transform = "DESeqTransform",
results = "list",
lfcShrink = "list"
),
prototype = prototype(
lfcShrink = list()
)
)
setValidity(
Class = "DESeqAnalysis",
method = function(object) {
data <- slot(object, "data")
transform <- slot(object, "transform")
results <- slot(object, "results")
lfcShrink <- slot(object, "lfcShrink")
ok <- validate(
identical(dimnames(data), dimnames(transform)),
msg = "DESeqDataSet and DESeqTransform must correspond."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- validate(
is.list(results),
is.list(lfcShrink),
msg = "results and lfcShrink must be list."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- validate(
all(bapply(
X = results,
FUN = function(x) {
identical(rownames(x), rownames(data))
}
)),
msg = "DESeqDataSet and DESeqResults must correspond."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- validate(
hasValidNames(results),
msg = "DESeqResults list must be named."
)
if (!isTRUE(ok)) {
return(ok)
}
## Alpha levels in the slotted results must be identical.
alphas <- vapply(
X = results,
FUN = function(x) {
metadata(x)[["alpha"]]
},
FUN.VALUE = numeric(1L)
)
ok <- validate(length(unique(alphas)) == 1L)
if (!isTRUE(ok)) {
return(ok)
}
## Note that `lfcShrink` slot is currently optional.
if (hasLength(lfcShrink)) {
ok <- validate(
identical(names(results), names(lfcShrink)),
all(as.logical(Map(
unshrunken = results,
shrunken = lfcShrink,
f = function(unshrunken, shrunken) {
identical(rownames(unshrunken), rownames(shrunken))
}
))),
msg = "Unshrunken and shrunken DESeqResults must correspond."
)
if (!isTRUE(ok)) {
return(ok)
}
## Ensure that DESeqResults slotted into `lfcShrink` is actually
## shrunken using the `lfcShrink()` function. This also checks to
## ensure that the same method was used for all contrasts.
shrinkTypes <- vapply(
X = lfcShrink,
FUN = lfcShrinkType,
FUN.VALUE = character(1L)
)
ok <- validate(
length(unique(shrinkTypes)) == 1L,
msg = "Invalid shrink type."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- validate(
identical(
vapply(
X = results,
FUN = function(x) {
metadata(x)[["alpha"]]
},
FUN.VALUE = numeric(1L)
),
vapply(
X = lfcShrink,
FUN = function(x) {
metadata(x)[["alpha"]]
},
FUN.VALUE = numeric(1L)
)
),
msg = "lfcShrink alpha must match the results alpha."
)
if (!isTRUE(ok)) {
return(ok)
}
}
TRUE
}
)
## DESeqAnalysisList ===========================================================
#' List containing related DESeq2 analyses
#'
#' @author Michael Steinbaugh
#' @export
#' @note Updated 2021-03-09.
#' @return `DESeqAnalysisList`.
setClass(
Class = "DESeqAnalysisList",
contains = "SimpleList"
)
setValidity(
Class = "DESeqAnalysisList",
method = function(object) {
## Currently allowing an empty list.
if (!hasLength(object)) {
return(TRUE)
}
## Require that all objects in list are named, without duplicates.
ok <- validate(hasValidNames(object))
if (!isTRUE(ok)) {
return(ok)
}
## Check that all of the elements in the list are DESeqAnalysis.
ok <- validate(
all(bapply(
X = object,
FUN = function(object) {
is(object, "DESeqAnalysis")
}
)),
msg = "Not a list of DESeqAnalysis objects."
)
if (!isTRUE(ok)) {
return(ok)
}
## Ensure that all slotted DESeqAnalysis objects are valid.
## This step can be slow for large objects.
ok <- validate(
all(bapply(object, validObject)),
msg = "Not all DESeqAnalysis objects are valid."
)
if (!isTRUE(ok)) {
return(ok)
}
## Check that all rownames in slotted DESeqResults are identical.
rn <- rownames(object[[1L]])
ok <- validate(
isCharacter(rn),
msg = "Row names in first DESeqAnalysis are invalid."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- bapply(
X = object,
rn = rn,
FUN = function(x, rn) {
identical(rownames(x), rn)
}
)
if (!all(ok)) {
return(sprintf(
"Row names mismatch in: %s.",
toInlineString(names(ok)[!ok], n = 5L)
))
}
TRUE
}
)
## DESeqResultsList ============================================================
#' List containing related DESeqResults objects
#'
#' @author Michael Steinbaugh
#' @export
#' @note Updated 2023-12-18.
#' @return `DESeqResultsList`.
setClass(
Class = "DESeqResultsList",
contains = "SimpleList"
)
setValidity(
Class = "DESeqResultsList",
method = function(object) {
## Currently allowing an empty list.
if (!hasLength(object)) {
return(TRUE)
}
## Require that all objects in list are named.
ok <- validate(hasValidNames(object))
if (!isTRUE(ok)) {
return(ok)
}
## Check that all of the elements in the list are DESeqAnalysis.
ok <- validate(
all(bapply(
X = object,
FUN = function(object) {
is(object, "DESeqResults")
}
)),
msg = "Not a list of DESeqResults objects."
)
if (!isTRUE(ok)) {
return(ok)
}
## Check that all rownames in slotted DESeqResults are identical.
rn <- rownames(object[[1L]])
ok <- validate(
isCharacter(rn),
msg = "Row names in first DESeqResults are invalid."
)
if (!isTRUE(ok)) {
return(ok)
}
ok <- bapply(
X = object,
rn = rn,
FUN = function(x, rn) {
identical(rownames(x), rn)
}
)
if (!all(ok)) {
return(sprintf(
"Row names mismatch in: %s.",
toInlineString(names(ok)[!ok], n = 5L)
))
}
TRUE
}
)