-
-
Notifications
You must be signed in to change notification settings - Fork 12
/
Copy pathSingleCellExperiment.R
767 lines (680 loc) · 26.2 KB
/
SingleCellExperiment.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
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
# nolint start: line_length_linter
# Sources used to understand how to convert between SingleCellExperiment and AnnData
# https://bioconductor.org/packages/release/bioc/vignettes/SingleCellExperiment/inst/doc/intro.html#2_Creating_SingleCellExperiment_instances
# https://bioconductor.org/packages/3.20/bioc/vignettes/SummarizedExperiment/inst/doc/SummarizedExperiment.html#column-sample-data
# https://github.com/ivirshup/sc-interchange/issues
# https://github.com/ivirshup/sc-interchange/issues/2
# https://www.bioconductor.org/packages/devel/bioc/vignettes/SingleCellExperiment/inst/doc/intro.html#3_Adding_low-dimensional_representations
# nolint end: line_length_linter
#' Convert an AnnData object to a SingleCellExperiment object
#'
#' `to_SingleCellExperiment()` converts an AnnData object
#' to a SingleCellExperiment object.
#'
#' @param adata an AnnData object, e.g., InMemoryAnnData
#'
#' @param assays_mapping A named vector mapping `layers` in `adata` to `assay`
#' names in the created SingleCellExperiment object. The names of the vector
#' should be the names of the `assays` in the resulting SingleCellExperiment
#' object, and the values should be the names of the `layers` in `adata`, and
#' can include the `X` matrix as well. If `X` is not in the list, it will be
#' added as `counts` or `data`.
#' @param colData_mapping A named vector mapping `obs` in `adata` to `colData`
#' in the created SingleCellExperiment object. The names of the vector should
#' be the names of the `colData` columns in the resulting SingleCellExperiment
#' object. The values should be the names of the `obs` columns in `adata`.
#' @param rowData_mapping A named vector mapping `var` names in `adata` to
#' `rowData` in the created SingleCellExperiment object. The names of the
#' vector should be the names of the `rowData` columns in the resulting
#' SingleCellExperiment object. The values should be the names of the `var`
#' columns in `adata`.
#' @param reduction_mapping A named vector or list mapping reduction names in
#' `adata` to reduction names in the created SingleCellExperiment object. For
#' the simpler named vector format, the names should be the names of the
#' `reducedDims` in the resulting SingleCellExperiment object, and the values
#' should be the names of the `obsm` in `adata`.
#'
#' For more advanced mapping, use the list format where each item has the
#' following keys:
#'
#' - `obsm`: the name of the `obsm` slot in `adata`
#' - `varm`: the name of the `varm` slot in `adata` (optional)
#' - `uns`: the name of the `uns` slot in `adata` (optional)
#'
#' If `'varm'` or `'uns'` is given, a
#' [SingleCellExperiment::LinearEmbeddingMatrix] will be created for that item
#' with `adata$varm[[varm]]` passed to the `featureLoadings` argument and
#' `adata$uns[[uns]]` passed as `metadata`
#' @param colPairs_mapping A named vector mapping obsp names in `adata` to
#' colPairs names in the created SingleCellExperiment object. The names of the
#' vector should be the names of the `colPairs` in the resulting
#' SingleCellExperiment object. The values should be the names of the `obsp`
#' in `adata`.
#' @param rowPairs_mapping A named vector mapping varp names in `adata` to
#' rowPairs names in the created SingleCellExperiment object. The names of the
#' vector should be the names of the `rowPairs` in the resulting
#' SingleCellExperiment object. The values should be the names of the `varp`
#' in `adata`.
#' @param metadata_mapping A named vector mapping uns names in `adata` to
#' metadata names in the created SingleCellExperiment object. The names of the
#' vector should be the names of the `metadata` in the resulting
#' SingleCellExperiment object. The values should be the names of the `uns` in
#' `adata`.
#'
#' @return `to_SingleCellExperiment()` returns a SingleCellExperiment
#' representing the content of `adata`.
#'
#' @details
#' If an unnamed vector is provided to a mapping argument the values will be
#' used as names
#'
#' @examples
#' if (interactive()) {
#' ## useful when interacting with the SingleCellExperiment !
#' library(SingleCellExperiment)
#' }
#' ad <- AnnData(
#' X = matrix(1:5, 3L, 5L),
#' layers = list(
#' A = matrix(5:1, 3L, 5L),
#' B = matrix(letters[1:5], 3L, 5L)
#' ),
#' obs = data.frame(row.names = LETTERS[1:3], cell = 1:3),
#' var = data.frame(row.names = letters[1:5], gene = 1:5)
#' )
#'
#' ## construct a SingleCellExperiment from an AnnData object
#' sce <- to_SingleCellExperiment(ad)
#' sce
#' @export
# nolint start: cyclocomp_linter
to_SingleCellExperiment <- function(
# nolint end: cyclocomp_linter
adata,
assays_mapping = NULL,
colData_mapping = NULL, # nolint
rowData_mapping = NULL, # nolint
reduction_mapping = NULL,
colPairs_mapping = NULL, # nolint
rowPairs_mapping = NULL, # nolint
metadata_mapping = NULL
) {
check_requires(
"Converting AnnData to SingleCellExperiment",
"SingleCellExperiment",
"Bioc"
)
if (!(inherits(adata, "AbstractAnnData"))) {
cli_abort(
"{.arg adata} must be a {.cls AbstractAnnData} but has class {.cls {class(adata)}}"
)
}
# guess mappings if not provided
# nolint start object_name_linter
assays_mapping <- self_name(assays_mapping) %||% .to_SCE_guess_assays(adata)
colData_mapping <- self_name(colData_mapping) %||%
.to_SCE_guess_all(adata, "obs")
rowData_mapping <- self_name(rowData_mapping) %||%
.to_SCE_guess_all(adata, "var")
reduction_mapping <- self_name(reduction_mapping) %||%
.to_SCE_guess_reduction(adata)
colPairs_mapping <- self_name(colPairs_mapping) %||%
.to_SCE_guess_all(adata, "obsp")
rowPairs_mapping <- self_name(rowPairs_mapping) %||%
.to_SCE_guess_all(adata, "varp")
metadata_mapping <- self_name(metadata_mapping) %||%
.to_SCE_guess_all(adata, "uns")
# nolint end object_name_linter
# trackstatus: class=SingleCellExperiment, feature=get_X, status=done
# trackstatus: class=SingleCellExperiment, feature=get_layers, status=done
sce_assays <- vector("list", length(assays_mapping))
names(sce_assays) <- names(assays_mapping)
for (i in seq_along(assays_mapping)) {
from <- assays_mapping[[i]]
to <- names(assays_mapping)[[i]]
if (from != "X") {
sce_assays[[to]] <- to_R_matrix(adata$layers[[from]])
} else {
sce_assays[[to]] <- to_R_matrix(adata$X)
}
}
# construct colData
# FIXME: probably better way to make a dataframe from a list of vectors
# trackstatus: class=SingleCellExperiment, feature=get_obs, status=done
# trackstatus: class=SingleCellExperiment, feature=get_obs_names, status=done
col_data <- .to_SCE_process_simple_mapping(adata, colData_mapping, "obs")
# construct rowData
# trackstatus: class=SingleCellExperiment, feature=get_var, status=done
# trackstatus: class=SingleCellExperiment, feature=get_var_names, status=done
row_data <- .to_SCE_process_simple_mapping(adata, rowData_mapping, "var")
# construct colPairs
# trackstatus: class=SingleCellExperiment, feature=get_obsp, status=done
col_pairs <- .to_SCE_process_simple_mapping(adata, colPairs_mapping, "obsp")
# construct rowPairs
# trackstatus: class=SingleCellExperiment, feature=get_varp, status=done
row_pairs <- .to_SCE_process_simple_mapping(adata, rowPairs_mapping, "varp")
# construct metadata
# trackstatus: class=SingleCellExperiment, feature=get_uns, status=done
metadata <- .to_SCE_process_simple_mapping(adata, metadata_mapping, "uns")
# construct reducedDims
reduceddims <- .to_SCE_process_reduction_mapping(adata, reduction_mapping)
arguments <- list(
assays = sce_assays,
colPairs = col_pairs,
rowPairs = row_pairs,
metadata = metadata,
checkDimnames = TRUE
)
# add col_data if not empty list
if (length(col_data) > 0) {
arguments$colData <- as(col_data, "DataFrame")
}
# add row_data if not empty list
if (length(row_data) > 0) {
arguments$rowData <- as(row_data, "DataFrame")
}
# construct output object
sce <- do.call(SingleCellExperiment::SingleCellExperiment, arguments)
rownames(sce) <- rownames(adata$var)
colnames(sce) <- rownames(adata$obs)
SingleCellExperiment::reducedDims(sce) <- reduceddims # only here to ensure that the dimensions are right
sce
}
# nolint start: object_length_linter object_name_linter
.to_SCE_guess_assays <- function(adata) {
# nolint end: object_length_linter object_name_linter
if (!(inherits(adata, "AbstractAnnData"))) {
cli_abort(
"{.arg adata} must be a {.cls AbstractAnnData} but has class {.cls {class(adata)}}"
)
}
layers <- list()
if (!is.null(adata$X)) {
layer_name_for_x <-
if (!"counts" %in% names(adata$layers)) {
# could expand checks, to check if integers
"counts"
} else {
"data"
}
layers[[layer_name_for_x]] <- "X"
}
for (layer_name in names(adata$layers)) {
layers[[layer_name]] <- layer_name
}
layers
}
# nolint start: object_length_linter object_name_linter
.to_SCE_guess_all <- function(adata, slot) {
# nolint end: object_length_linter object_name_linter
if (!(inherits(adata, "AbstractAnnData"))) {
cli_abort(
"{.arg adata} must be a {.cls AbstractAnnData} but has class {.cls {class(adata)}}"
)
}
self_name(names(adata[[slot]]))
}
# nolint start: object_length_linter object_name_linter
.to_SCE_guess_reduction <- function(adata) {
# nolint end: object_length_linter object_name_linter
.to_Seurat_guess_reductions(adata) # nolint object_usage_linter
}
# nolint start: object_length_linter object_name_linter
.to_SCE_process_simple_mapping <- function(adata, mapping, slot) {
# nolint end: object_length_linter object_name_linter
if (rlang::is_empty(mapping)) {
return(list())
}
purrr::map(mapping, function(.item) {
adata[[slot]][[.item]]
})
}
# nolint start: object_length_linter object_name_linter
.to_SCE_process_reduction <- function(
# nolint end: object_length_linter object_name_linter
adata,
obsm_key,
varm_key,
uns_key
) {
if (!(obsm_key %in% adata$obsm_keys())) {
cli_abort(
c(
"{.val {obsm_key}} is not an item in {.code adata$obsm}",
"i" = "{.code adata$obsm_keys()}: {.val {adata$obsm_keys()}}"
)
)
}
embedding <- adata$obsm[[obsm_key]]
rownames(embedding) <- adata$obs_names
# If there is no extra info to add, return the embedding matrix
if (is.null(varm_key) && is.null(uns_key)) {
return(embedding)
}
if (!is.null(varm_key)) {
if (!(varm_key %in% adata$varm_keys())) {
cli_abort(
c(
"{.val {varm_key}} is not an item in {.code adata$varm}",
"i" = "{.code adata$varm_keys()}: {.val {adata$varm_keys()}}"
)
)
}
loadings <- adata$varm[[varm_key]]
rownames(loadings) <- colnames(embedding)
} else {
loadings <- matrix(nrow = 0, ncol = ncol(embedding))
}
if (!is.null(uns_key)) {
if (!(uns_key %in% adata$uns_keys())) {
cli_abort(
c(
"{.val {uns_key}} is not an item in {.code adata$uns}",
"i" = "{.code adata$uns_keys()}: {.val {adata$uns_keys()}}"
)
)
}
metadata <- adata$uns[[uns_key]]
} else {
metadata <- list()
}
lem <- SingleCellExperiment::LinearEmbeddingMatrix(
sampleFactors = embedding,
featureLoadings = loadings,
metadata = metadata
)
rownames(lem) <- rownames(embedding)
colnames(lem) <- colnames(loadings)
lem
}
# nolint start: object_length_linter object_name_linter
.to_SCE_process_reduction_mapping <- function(adata, reduction_mapping) {
# nolint end: object_length_linter object_name_linter
# If the mapping is a vector expand it to the list format
if (is.atomic(reduction_mapping)) {
reduction_mapping <- purrr::map(reduction_mapping, function(.obsm) {
c(obsm = .obsm)
})
} else {
purrr::walk(names(reduction_mapping), function(.name) {
mapping <- reduction_mapping[[.name]]
if (!(is.character(mapping) && ("obsm" %in% names(mapping)))) {
cli_abort(c(
paste(
"If it is a list, each item in {.arg reduction_mapping} must be a
{.cls character} vector with the name {.val obsm} and optional",
"names {.val varm} and {.val uns}"
),
"i" = paste(
"{.code reduction_mapping[[{.val { .name }}]]}:",
"{style_vec(mapping, wrap = TRUE)}"
)
))
}
})
}
purrr::map(reduction_mapping, function(.map) {
varm <- if ("varm" %in% names(.map)) .map[["varm"]] else NULL
uns <- if ("uns" %in% names(.map)) .map[["uns"]] else NULL
.to_SCE_process_reduction(
adata,
obsm_key = .map[["obsm"]],
varm_key = varm,
uns_key = uns
)
})
}
#' Convert a SingleCellExperiment object to an AnnData object
#'
#' `from_SingleCellExperiment()` converts a
#' SingleCellExperiment to an AnnData object.
#'
#' @param sce An object inheriting from SingleCellExperiment.
#'
#' @param output_class Name of the AnnData class. Must be one of `"HDF5AnnData"`
#' or `"InMemoryAnnData"`.
#'
#' @param x_mapping Name of the assay in `sce` to use as the `X` matrix in the AnnData object.
#' @param layers_mapping A named vector or list mapping layer names in AnnData to assay names in the
#' SingleCellExperiment object. Each name corresponds to the layer name in AnnData, and each value to the
#' assay name in SCE.
#' @param obs_mapping A named vector or list mapping column names in AnnData's obs to column names in SCE's colData.
#' Each name corresponds to a column name in AnnData's obs, and each value to a column name in SCE's colData.
#' @param var_mapping A named vector or list mapping column names in AnnData's var to column names in SCE's rowData.
#' Each name corresponds to a column name in AnnData's var, and each value to a column name in SCE's rowData.
#' @param obsm_mapping A named vector mapping obsm keys in AnnData to reducedDims names in SCE.
#' Each name corresponds to a key in AnnData's obsm, and each value to the name of a reducedDim in SCE.
#' Example: `obsm_mapping = c(X_pca = "pca", X_umap = "umap")`.
#' @param varm_mapping A named vector mapping varm keys in AnnData to reducedDims names in SCE.
#' Each name corresponds to a key in AnnData's varm, and each value to the name of a reducedDim in SCE.
#' Example: `varm_mapping = c(PCs = "pca")`.
#' @param obsp_mapping A named vector or list mapping obsp keys in AnnData to colPairs in the
#' SingleCellExperiment object. Each name corresponds to a key in AnnData's obsp, and each value to a name in
#' SCE's colPairs.
#' @param varp_mapping A named vector or list mapping varp keys in AnnData to rowPairs in the
#' SingleCellExperiment object. Each name corresponds to a key in AnnData's varp, and each value to a
#' name in SCE's rowPairs.
#' @param uns_mapping A named vector or list mapping uns keys in AnnData to metadata in the
#' SingleCellExperiment object. Each name corresponds to a key in AnnData's uns, and each value to a
#' name in SCE's metadata.
#' @param ... Additional arguments to pass to the generator function.
#'
#' @return `from_SingleCellExperiment()` returns an AnnData object
#' (e.g., InMemoryAnnData) representing the content of `sce`.
#'
#' @examples
#' ## construct an AnnData object from a SingleCellExperiment
#' library(SingleCellExperiment)
#' sce <- SingleCellExperiment(
#' assays = list(counts = matrix(1:5, 5L, 3L)),
#' colData = DataFrame(cell = 1:3, row.names = paste0("Cell", 1:3)),
#' rowData = DataFrame(gene = 1:5, row.names = paste0("Gene", 1:5))
#' )
#' from_SingleCellExperiment(sce, "InMemory")
#'
#' @export
# nolint start: object_name_linter
from_SingleCellExperiment <- function(
# nolint end: object_name_linter
sce,
output_class = c("InMemory", "HDF5AnnData"),
x_mapping = NULL,
layers_mapping = NULL,
obs_mapping = NULL,
var_mapping = NULL,
obsm_mapping = NULL,
varm_mapping = NULL,
obsp_mapping = NULL,
varp_mapping = NULL,
uns_mapping = NULL,
...
) {
check_requires(
"Converting SingleCellExperiment to AnnData",
"SingleCellExperiment",
"Bioc"
)
output_class <- match.arg(output_class)
if (!(inherits(sce, "SingleCellExperiment"))) {
cli_abort(
"{.arg sce} must be a {.cls SingleCellExperiment} but has class {.cls {sce}}"
)
}
# For any mappings that are not set, using the guessing function
layers_mapping <- self_name(layers_mapping) %||%
.from_SCE_guess_layers(sce, x_mapping)
obs_mapping <- self_name(obs_mapping) %||%
.from_SCE_guess_all(
sce,
SingleCellExperiment::colData
)
var_mapping <- self_name(var_mapping) %||%
.from_SCE_guess_all(
sce,
SingleCellExperiment::rowData
)
obsm_mapping <- self_name(obsm_mapping) %||% .from_SCE_guess_obsm(sce)
varm_mapping <- self_name(varm_mapping) %||% .from_SCE_guess_varm(sce)
obsp_mapping <- self_name(obsp_mapping) %||%
.from_SCE_guess_obspvarp(
sce,
SingleCellExperiment::colPairs
)
varp_mapping <- self_name(varp_mapping) %||%
.from_SCE_guess_obspvarp(
sce,
SingleCellExperiment::rowPairs
)
uns_mapping <- self_name(uns_mapping) %||%
.from_SCE_guess_all(sce, S4Vectors::metadata)
generator <- get_anndata_constructor(output_class)
adata <- generator$new(shape = rev(dim(sce)), ...)
# Fill in slots in the object
.from_SCE_process_obs(adata, sce, obs_mapping)
.from_SCE_process_var(adata, sce, var_mapping)
# trackstatus: class=SingleCellExperiment, feature=set_X, status=wip
if (!is.null(x_mapping)) {
adata$X <- .from_SCE_convert(
SummarizedExperiment::assay(sce, x_mapping, withDimnames = FALSE)
)
}
.from_SCE_process_layers(adata, sce, layers_mapping)
.from_SCE_process_obsm(adata, sce, obsm_mapping)
.from_SCE_process_varm(adata, sce, varm_mapping)
.from_SCE_process_obsp(adata, sce, obsp_mapping)
.from_SCE_process_varp(adata, sce, varp_mapping)
.from_SCE_process_uns(adata, sce, uns_mapping)
adata
}
# nolint start: object_length_linter object_name_linter
.from_SCE_guess_all <- function(sce, slot) {
# nolint end: object_length_linter object_name_linter
self_name(names(slot(sce)))
}
# nolint start: object_length_linter object_name_linter
.from_SCE_guess_layers <- function(sce, x_mapping) {
# nolint end: object_length_linter object_name_linter
layers_mapping <- self_name(SummarizedExperiment::assayNames(sce))
if (!is.null(x_mapping)) {
layers_mapping <- layers_mapping[layers_mapping != x_mapping]
}
if (rlang::is_empty(layers_mapping)) {
c()
} else {
layers_mapping
}
}
# nolint start: object_length_linter object_name_linter
.from_SCE_guess_obsm <- function(sce) {
# nolint end: object_length_linter object_name_linter
obsm_mapping <- self_name(SingleCellExperiment::reducedDimNames(sce))
if (rlang::is_empty(obsm_mapping)) {
c()
} else {
obsm_mapping
}
}
# nolint start: object_length_linter object_name_linter
.from_SCE_guess_varm <- function(sce) {
# nolint end: object_length_linter object_name_linter
varm_mapping <- c()
for (reduction_name in names(SingleCellExperiment::reducedDims(sce))) {
reduction <- SingleCellExperiment::reducedDim(sce, reduction_name)
if (inherits(reduction, "LinearEmbeddingMatrix")) {
varm_mapping[reduction_name] <- reduction_name
}
}
varm_mapping
}
# If sce is a SummarizedExperiment, return an empty mapping
# nolint start: object_length_linter object_name_linter
.from_SCE_guess_obspvarp <- function(sce, slot) {
# nolint end: object_length_linter object_name_linter
.from_SCE_guess_all(sce, slot)
}
# Convert BioConductor-specific objects to base R objects
# Convert matrices
# Convert dgCMatrix to dgRMatrix
# nolint start: object_length_linter object_name_linter
.from_SCE_convert <- function(object, transpose = TRUE) {
# nolint end: object_length_linter object_name_linter
if (inherits(object, "DataFrame")) {
as.data.frame(object)
} else if (inherits(object, "SimpleList")) {
as.list(object)
} else if (inherits(object, "matrix") || inherits(object, "Matrix")) {
if (inherits(object, "denseMatrix")) {
object <- as.matrix(object)
}
if (transpose) {
to_py_matrix(object)
} else {
object
}
} else {
object
}
}
# trackstatus: class=SingleCellExperiment, feature=set_obs_names, status=wip
# trackstatus: class=SingleCellExperiment, feature=set_obs, status=wip
# nolint start: object_name_linter
.from_SCE_process_obs <- function(adata, sce, obs_mapping) {
# nolint end: object_name_linter
if (!rlang::is_empty(obs_mapping)) {
adata$obs <- SummarizedExperiment::colData(sce) |>
as.data.frame() |>
setNames(names(obs_mapping))
} else {
# Store an empty data.frame to keep the obs names
if (is.null(colnames(sce))) {
adata$obs <- data.frame(matrix(nrow = ncol(sce), ncol = 0))
} else {
adata$obs <- data.frame(row.names = colnames(sce))
}
}
}
# trackstatus: class=SingleCellExperiment, feature=set_var_names, status=wip
# trackstatus: class=SingleCellExperiment, feature=set_var, status=wip
# nolint start: object_name_linter
.from_SCE_process_var <- function(adata, sce, var_mapping) {
# nolint end: object_name_linter
if (!rlang::is_empty(var_mapping)) {
adata$var <- SummarizedExperiment::rowData(sce) |>
as.data.frame() |>
setNames(names(var_mapping))
} else {
# Store an empty data.frame to keep the var names
if (is.null(rownames(sce))) {
adata$var <- data.frame(matrix(nrow = nrow(sce), ncol = 0))
} else {
adata$var <- data.frame(row.names = rownames(sce))
}
}
}
# trackstatus: class=SingleCellExperiment, feature=set_layers, status=wip
# nolint start: object_length_linter object_name_linter
.from_SCE_process_layers <- function(adata, sce, layers_mapping) {
# nolint end: object_length_linter object_name_linter
if (rlang::is_empty(layers_mapping)) {
return(invisible())
}
adata$layers <- purrr::map(layers_mapping, function(.layer) {
.from_SCE_convert(SummarizedExperiment::assay(sce, .layer))
})
}
# trackstatus: class=SingleCellExperiment, feature=set_obsm, status=wip
# nolint start: object_name_linter
.from_SCE_process_obsm <- function(adata, sce, obsm_mapping) {
# nolint end: object_name_linter
if (rlang::is_empty(obsm_mapping)) {
return(invisible())
}
adata$obsm <- purrr::imap(obsm_mapping, function(reduction_name, obsm_key) {
# Check if the reduction exists
if (!reduction_name %in% names(SingleCellExperiment::reducedDims(sce))) {
cli_abort(c(
"Reduction {.val {reduction_name}} not found in SCE object.",
"i" = "Available reductions: {.val {names(SingleCellExperiment::reducedDims(sce))}}"
))
}
reduction <- SingleCellExperiment::reducedDim(sce, reduction_name)
if (inherits(reduction, "LinearEmbeddingMatrix")) {
SingleCellExperiment::sampleFactors(reduction)
} else {
reduction
}
})
}
# trackstatus: class=SingleCellExperiment, feature=set_varm, status=wip
# nolint start: object_name_linter
.from_SCE_process_varm <- function(adata, sce, varm_mapping) {
# nolint end: object_name_linter
if (rlang::is_empty(varm_mapping)) {
return(invisible())
}
adata$varm <- purrr::map(varm_mapping, function(reduction_name) {
# Check if the reduction exists
if (!reduction_name %in% names(SingleCellExperiment::reducedDims(sce))) {
cli_abort(c(
"Reduction {.val {reduction_name}} not found in SCE object.",
"i" = "Available reductions: {.val {names(SingleCellExperiment::reducedDims(sce))}}"
))
}
reduction <- SingleCellExperiment::reducedDim(sce, reduction_name)
if (!inherits(reduction, "LinearEmbeddingMatrix")) {
cli_abort(paste(
"{.code reducedDim(sce, {.val {reduction_name}}} must be a {.cls LinearEmbeddingMatrix}",
"but has class {.cls {class(reduction)[1]}}"
))
}
loadings <- SingleCellExperiment::featureLoadings(reduction)
if (nrow(loadings) != adata$n_vars()) {
cli_abort(paste(
"The number of rows ({.val {nrow(loadings)}}) in the loadings for ",
"reduction {.val {reduction_name}} does not match {.code adata$n_vars()}",
"({.val {adata$n_vars()}})"
))
}
loadings
})
}
# trackstatus: class=SingleCellExperiment, feature=set_obsp, status=wip
# nolint start: object_length_linter object_name_linter
.from_SCE_process_obsp <- function(adata, sce, obsp_mapping) {
# nolint end: object_length_linter object_name_linter
if (rlang::is_empty(obsp_mapping)) {
return(invisible())
}
adata$obsp <- purrr::map(obsp_mapping, function(colpair_name) {
# Check if the colPair exists
if (!(colpair_name %in% SingleCellExperiment::colPairNames(sce))) {
cli_abort(c(
"colPair {.val {colpair_name}} not found in SCE object.",
"i" = "Available colPairs: {.val {names(SingleCellExperiment::colPairs(sce))}}"
))
}
.from_SCE_convert(
SingleCellExperiment::colPair(sce, colpair_name, asSparse = TRUE),
transpose = FALSE
)
})
}
# trackstatus: class=SingleCellExperiment, feature=set_varp, status=wip
# nolint start: object_length_linter object_name_linter
.from_SCE_process_varp <- function(adata, sce, varp_mapping) {
# nolint end: object_length_linter object_name_linter
if (rlang::is_empty(varp_mapping)) {
return(invisible())
}
adata$varp <- purrr::map(varp_mapping, function(rowpair_name) {
# Check if the rowPair exists
if (!(rowpair_name %in% SingleCellExperiment::rowPairNames(sce))) {
cli_abort(c(
"rowPair {.val {rowpair_name}} not found in SCE object.",
"i" = "Available rowPairs: {.val {names(SingleCellExperiment::rowPairs(sce))}}"
))
}
.from_SCE_convert(
SingleCellExperiment::rowPair(sce, rowpair_name, asSparse = TRUE),
transpose = FALSE
)
})
}
# trackstatus: class=SingleCellExperiment, feature=set_uns, status=wip
# nolint start: object_length_linter object_name_linter
.from_SCE_process_uns <- function(adata, sce, uns_mapping) {
# nolint end: object_length_linter object_name_linter
if (rlang::is_empty(uns_mapping)) {
return(invisible())
}
adata$uns <- purrr::map(uns_mapping, function(meta_name) {
# Check if the metadata exists
if (!(meta_name %in% names(S4Vectors::metadata(sce)))) {
cli_abort(c(
"Metadata {.val {meta_name}} not found in SCE object.",
"i" = "Available metadata: {.val {names(S4Vectors::metadata(sce))}}"
))
}
S4Vectors::metadata(sce)[[meta_name]]
})
}