/
plotCountsPerBiotype-methods.R
210 lines (189 loc) · 6.1 KB
/
plotCountsPerBiotype-methods.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
#' @name plotCountsPerBiotype
#' @author Michael Steinbaugh, Rory Kirchner
#' @inherit AcidGenerics::plotCountsPerBiotype
#' @note Updated 2023-12-04.
#'
#' @inheritParams AcidRoxygen::params
#' @param ... Additional arguments.
#'
#' @param biotypeCol `character(1)`.
#' Biotype column name defined in `colData()`.
#'
#' @examples
#' data(
#' RangedSummarizedExperiment,
#' SingleCellExperiment_splatter,
#' package = "AcidTest"
#' )
#'
#' ## SummarizedExperiment ====
#' object <- RangedSummarizedExperiment
#' plotCountsPerBiotype(object)
#'
#' ## SingleCellExperiment ====
#' object <- SingleCellExperiment_splatter
#' plotCountsPerBiotype(object)
NULL
## Updated 2021-09-10.
`plotCountsPerBiotype,SE` <- # nolint
function(object,
assay = 1L,
biotypeCol = "geneBiotype",
n = 9L,
interestingGroups = NULL,
geom = c("violin", "boxplot"),
trans = c("identity", "log2", "log10"),
labels = list(
"title" = "Counts per biotype",
"subtitle" = NULL,
"sampleAxis" = NULL,
"countAxis" = "counts"
)) {
assert(
validObject(object),
isScalar(assay),
isString(biotypeCol),
isInt(n)
)
geom <- match.arg(geom)
trans <- match.arg(trans)
labels <- matchLabels(labels)
interestingGroups(object) <-
matchInterestingGroups(object, interestingGroups)
interestingGroups <- interestingGroups(object)
rowData <- rowData(object)
if (!isSubset(biotypeCol, colnames(rowData))) {
alertWarning(sprintf(
"{.fun %s} does not contain biotypes defined in {.val %s}.",
"rowData", biotypeCol
))
return(invisible(NULL))
}
rowData <- decode(rowData)
rowData[["rowname"]] <- as.factor(rownames(object))
## Get the top biotypes from the row data.
biotypes <- table(rowData[[biotypeCol]])
## Requiring at least 10 genes per biotype.
biotypes <- biotypes[which(biotypes > 10L)]
biotypes <- sort(biotypes, decreasing = TRUE)
biotypes <- head(biotypes, n = n)
biotypes <- names(biotypes)
## Melt the count matrix into long format.
data <- melt(
object = object,
assay = assay,
min = 1L,
minMethod = "perRow",
trans = trans
)
data <- decode(data)
## Prepare the minimal data frame required for plotting.
data <- leftJoin(x = data, y = rowData, by = "rowname")
keep <- which(data[[biotypeCol]] %in% biotypes)
data <- data[keep, , drop = FALSE]
data <- data[, c("value", "interestingGroups", biotypeCol)]
## Sanitize the biotype column to appear nicer in plots.
data[[biotypeCol]] <- gsub("_", " ", data[[biotypeCol]])
## Plot.
p <- ggplot(
data = as.data.frame(data),
mapping = aes(
x = .data[["interestingGroups"]],
y = .data[["value"]]
)
)
switch(
EXPR = geom,
"violin" = {
p <- p +
geom_violin(
mapping = aes(
fill = .data[["interestingGroups"]]
),
color = NA,
scale = "area",
trim = TRUE
)
},
"boxplot" = {
p <- p +
geom_boxplot(
mapping = aes(
color = .data[["interestingGroups"]]
),
fill = NA
)
}
)
p <- p +
scale_y_continuous(
breaks = pretty_breaks(),
labels = prettyNum
) +
facet_wrap(facets = vars(.data[[biotypeCol]]), scales = "free_y")
## Labels.
if (!identical(trans, "identity")) {
labels[["countAxis"]] <- paste(trans, labels[["countAxis"]])
}
labels[["color"]] <- paste(interestingGroups, collapse = ":\n")
labels[["fill"]] <- labels[["color"]]
names(labels)[names(labels) == "sampleAxis"] <- "x"
names(labels)[names(labels) == "countAxis"] <- "y"
p <- p + do.call(what = labs, args = labels)
## Color palette.
p <- p + acid_scale_color_discrete()
p <- p + acid_scale_fill_discrete()
## Return.
p
}
## Updated 2021-09-08.
`plotCountsPerBiotype,SCE` <- # nolint
`plotCountsPerBiotype,SE`
## Updated 2021-09-09.
`plotCountsPerBroadClass,SE` <- # nolint
function(object,
...,
labels = list(
"title" = "Counts per broad class biotype",
"subtitle" = NULL,
"sampleAxis" = NULL,
"countAxis" = "counts"
)) {
plotCountsPerBiotype(
object = object,
biotypeCol = "broadClass",
labels = matchLabels(labels),
...
)
}
`plotCountsPerBroadClass,SCE` <- # nolint
`plotCountsPerBroadClass,SE`
#' @rdname plotCountsPerBiotype
#' @export
setMethod(
f = "plotCountsPerBiotype",
signature = signature(object = "SingleCellExperiment"),
definition = `plotCountsPerBiotype,SCE`
)
#' @rdname plotCountsPerBiotype
#' @export
setMethod(
f = "plotCountsPerBiotype",
signature = signature(object = "SummarizedExperiment"),
definition = `plotCountsPerBiotype,SE`
)
#' @rdname plotCountsPerBiotype
#' @export
setMethod(
f = "plotCountsPerBroadClass",
signature = signature(object = "SingleCellExperiment"),
definition = `plotCountsPerBroadClass,SCE`
)
## NOTE This is currently used in bcbioRNASeq R Markdown template.
#' @rdname plotCountsPerBiotype
#' @export
setMethod(
f = "plotCountsPerBroadClass",
signature = signature(object = "SummarizedExperiment"),
definition = `plotCountsPerBroadClass,SE`
)