-
Notifications
You must be signed in to change notification settings - Fork 74
/
importSeqc.R
248 lines (218 loc) · 8.87 KB
/
importSeqc.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
.constructSCEFromSeqcOutputs <- function(
sampleName,
matrix,
features,
barcodes) {
coln <- paste(sampleName, barcodes[[1]], sep = "_")
rownames(matrix) <- features[[1]]
sce <- SingleCellExperiment::SingleCellExperiment(
assays = list(counts = matrix))
SummarizedExperiment::rowData(sce) <- features
SummarizedExperiment::colData(sce) <- S4Vectors::DataFrame(
cell_barcode = as.character(barcodes[[1]]),
column_name = coln,
sample = sampleName,
row.names = coln)
return(sce)
}
.unionGeneMatrix <- function(geneUnion, matrix){
missGene <- geneUnion[!geneUnion %in% rownames(matrix)]
missMat <- Matrix::Matrix(0, nrow = length(missGene), ncol = ncol(matrix),
dimnames = list(missGene, NULL))
matb <- methods::as(matrix, "CsparseMatrix")
rownames(matb) <- rownames(matrix)
mat <- rbind(matb, missMat)
if (anyDuplicated(rownames(mat))) {
mat <- mat[!duplicated(rownames(mat)), ]
warning("Duplicated genes exist in count matrix. Filtered",
" duplicated genes.")
}
return(mat)
}
.getGeneUnion <- function(geneList){
gene <- geneList
for (i in seq_along(geneList)){
gene[[i]] <- geneList[[i]][[1]]
}
geneUnion <- base::Reduce(union, gene)
return(geneUnion)
}
.readBarcodesSEQC <- function(path) {
res <- data.table::fread(path, header = FALSE, sep=",", colClasses = "character")
res <- res[,-1,drop = FALSE]
colnames(res) <- "cell_barcode"
return(res)
}
.readFeaturesSEQC <- function(path) {
res <- data.table::fread(path, header = FALSE, sep=",", colClasses = "character")
res <- res[,-1,drop = FALSE]
colnames(res) <- "feature_name"
return(res)
}
.importSEQC <- function(
seqcDirs,
samples,
prefix,
gzipped,
class,
delayedArray,
cbNotFirstCol,
feNotFirstCol,
combinedSample,
rowNamesDedup) {
if (length(seqcDirs) != length(samples)) {
stop("'seqcDirs' and 'samples' have unequal lengths!")
}
if (length(seqcDirs) != length(prefix)) {
stop("'seqcDirs' and 'prefix' have unequal lengths!")
}
res <- vector("list", length = length(seqcDirs))
cb <- vector("list", length = length(seqcDirs))
fe <- vector("list", length = length(seqcDirs))
mat <- vector("list", length = length(seqcDirs))
for (i in seq_along(seqcDirs)) {
dir <- seqcDirs[i]
matrixFile <- paste(prefix[i], 'sparse_molecule_counts.mtx', sep = "_")
featuresFile <- paste(prefix[i], 'sparse_counts_genes.csv', sep = "_")
barcodesFile <- paste(prefix[i], 'sparse_counts_barcodes.csv',
sep = "_")
cb[[i]] <- .readBarcodesSEQC(file.path(dir, barcodesFile))
fe[[i]] <- .readFeaturesSEQC(file.path(dir, featuresFile))
mat[[i]] <- .readMatrixMM(file.path(dir, matrixFile),
gzipped = gzipped, class = class, delayedArray = delayedArray)
mat[[i]] <- t(mat[[i]])
rownames(mat[[i]]) <- fe[[i]][[1]]
}
if (isTRUE(combinedSample) & length(seqcDirs) > 1) {
geneUnion <- .getGeneUnion(fe)
for (i in seq_along(seqcDirs)) {
matrix <- .unionGeneMatrix(geneUnion = geneUnion, matrix = mat[[i]])
matrix <- matrix[geneUnion, ]
feature <- S4Vectors::DataFrame('feature_name' = rownames(matrix))
scei <- .constructSCEFromSeqcOutputs(
sampleName = samples[i],
matrix = matrix,
features = feature,
barcodes = cb[[i]])
res[[i]] <- scei
}
sce <- do.call(SingleCellExperiment::cbind, res)
if (isTRUE(rowNamesDedup)) {
if (any(duplicated(rownames(sce)))) {
message("Duplicated gene names found, adding '-1', '-2', ",
"... suffix to them.")
}
sce <- dedupRowNames(sce)
}
return(sce)
} else {
for (i in seq_along(seqcDirs)) {
scei <- .constructSCEFromSeqcOutputs(
sampleName = samples[i],
matrix = mat[[i]],
features = fe[[i]],
barcodes = cb[[i]])
res[[i]] <- scei
}
if (length(seqcDirs) == 1) {
sce <- res[[1]]
if (isTRUE(rowNamesDedup)) {
if (any(duplicated(rownames(sce)))) {
message("Duplicated gene names found, adding '-1', '-2', ",
"... suffix to them.")
}
sce <- dedupRowNames(sce)
}
return(sce)
} else {
return(res)
}
}
}
#' @name importSEQC
#' @rdname importSEQC
#' @title Construct SCE object from seqc output
#' @description Read the filtered barcodes, features, and matrices for all
#' samples from (preferably a single run of) seqc output. Import and combine
#' them as one big \link[SingleCellExperiment]{SingleCellExperiment} object.
#' @param seqcDirs A vector of paths to seqc output files. Each sample
#' should have its own path. For example: \code{"./pbmc_1k_50x50"}. Must have
#' the same length as \code{samples}.
#' @param samples A vector of user-defined sample names for the samples to be
#' imported. Must have the same length as \code{seqcDirs}.
#' @param prefix A vector containing the prefix of file names within each
#' sample directory. It cannot be null and the vector should have the same
#' length as \emph{samples}.
#' @param gzipped Boolean. \code{TRUE} if the seqc output files
#' (sparse_counts_barcode.csv, sparse_counts_genes.csv, and
#' sparse_molecule_counts.mtx) were gzip compressed. \code{FALSE} otherwise.
#' Default seqc outputs are not gzipped. Default \code{FALSE}.
#' @param class Character. The class of the expression matrix stored in the SCE
#' object. Can be one of \code{"Matrix"} (as returned by \link{readMM}
#' function), or \code{"matrix"} (as returned by \link[base]{matrix} function).
#' Default \code{"Matrix"}.
#' @param delayedArray Boolean. Whether to read the expression matrix as
#' \link{DelayedArray} object or not. Default \code{FALSE}.
#' @param feNotFirstCol Boolean. \code{TRUE} if first column of
#' sparse_counts_genes.csv is row index and it will be removed. \code{FALSE}
#' the first column will be kept.
#' @param cbNotFirstCol Boolean. \code{TRUE} if first column of
#' sparse_counts_barcode.csv is row index and it will be removed. \code{FALSE}
#' the first column will be kept.
#' @param combinedSample Boolean. If \code{TRUE}, \code{importSEQC} returns a
#' \code{SingleCellExperiment} object containing the combined count matrix,
#' feature annotations and the cell annotations. If \code{FALSE},
#' \code{importSEQC} returns a list containing multiple
#' \code{SingleCellExperiment} objects. Each \code{SingleCellExperiment}
#' contains count matrix, feature annotations and cell annotations for
#' each sample.
#' @param rowNamesDedup Boolean. Whether to deduplicate rownames. Only applied
#' if \code{combinedSample} is \code{TRUE} or only one \code{seqcDirs}
#' specified. Default \code{TRUE}.
#' @details
#' \code{importSEQC} imports output from seqc. The default
#' sparse_counts_barcode.csv or sparse_counts_genes.csv from seqc output
#' contains two columns. The first column is row index and the second column
#' is cell-barcode or gene symbol. \code{importSEQC} will remove first column.
#' Alternatively, user can call
#' \code{cbNotFirstCol} or \code{feNotFirstCol} as FALSE to keep the first
#' column of these files. When \code{combinedSample} is TRUE, \code{importSEQC}
#' will combined count matrix with genes detected in at least one sample.
#' @return A \code{SingleCellExperiment} object containing the combined count
#' matrix, the feature annotations, and the cell annotation.
#' @examples
#' # Example #1
#' # The following filtered feature, cell, and matrix files were downloaded from
#' # https://support.10xgenomics.com/single-cell-gene-expression/datasets/
#' # 3.0.0/pbmc_1k_v3
#' # The top 50 hg38 genes are included in this example.
#' # Only the top 50 cells are included.
#' sce <- importSEQC(
#' seqcDirs = system.file("extdata/pbmc_1k_50x50", package = "singleCellTK"),
#' samples = "pbmc_1k_50x50",
#' prefix = "pbmc_1k",
#' combinedSample = FALSE)
#' @export
importSEQC <- function(
seqcDirs = NULL,
samples = NULL,
prefix = NULL,
gzipped = FALSE,
class = c("Matrix", "matrix"),
delayedArray = FALSE,
cbNotFirstCol = TRUE,
feNotFirstCol = TRUE,
combinedSample = TRUE,
rowNamesDedup = TRUE) {
class <- match.arg(class)
.importSEQC(seqcDirs = seqcDirs,
samples = samples,
prefix = prefix,
gzipped = gzipped,
class = class,
delayedArray = delayedArray,
cbNotFirstCol = cbNotFirstCol,
feNotFirstCol = feNotFirstCol,
combinedSample = combinedSample,
rowNamesDedup = rowNamesDedup)
}