Skip to content

Commit

Permalink
fix v5 import
Browse files Browse the repository at this point in the history
  • Loading branch information
alexvpickering committed May 22, 2024
1 parent 34aa1dd commit 2b1fe59
Show file tree
Hide file tree
Showing 2 changed files with 139 additions and 11 deletions.
13 changes: 9 additions & 4 deletions R/seurat_to_sce.R
Original file line number Diff line number Diff line change
Expand Up @@ -5,8 +5,13 @@ seurat_to_sce <- function(sdata, dataset_name) {
Seurat::DefaultAssay(sdata) <- 'RNA'

# join layers if v5
if (methods::is(sdata[['RNA']], 'Assay5')) {
# use meta.features (v3) or meta.data (v5)
meta_name <- 'meta.features'
is.v5 <- methods::is(sdata[['RNA']], 'Assay5')

if (is.v5) {
sdata[['RNA']] <- SeuratObject::JoinLayers(sdata[['RNA']])
meta_name <- 'meta.data'
}

# normalize if need to
Expand All @@ -20,7 +25,7 @@ seurat_to_sce <- function(sdata, dataset_name) {
# meta.features need to have same row.names as assay
for (assay_name in Seurat::Assays(sdata)) {
assay <- sdata[[assay_name]]
assay@meta.features <- assay@meta.features[row.names(assay), ]
slot(assay, meta_name) <- slot(assay, meta_name)[row.names(assay), ]
sdata[[assay_name]] <- assay
}

Expand Down Expand Up @@ -51,14 +56,14 @@ seurat_to_sce <- function(sdata, dataset_name) {
sce$cluster <- unname(Seurat::Idents(sdata))

# default HVGs unless integrated assay
hvgs <- sdata[['RNA']]@var.features
hvgs <- SeuratObject::VariableFeatures(sdata, assay = 'RNA')

is.integrated <- have.samples | have.integrated
if (is.integrated) {

# get HVGs and transfer reductions from integrated assay if present
if (have.integrated) {
hvgs <- sdata[['integrated']]@var.features
hvgs <- SeuratObject::VariableFeatures(sdata, assay = 'integrated')
SingleCellExperiment::reducedDims(sce) <-
SingleCellExperiment::reducedDims(SingleCellExperiment::altExp(sce, 'integrated'))
}
Expand Down
137 changes: 130 additions & 7 deletions tests/testthat/test_import_scseq.R
Original file line number Diff line number Diff line change
Expand Up @@ -89,13 +89,6 @@ test_that("cellranger .h5 files with multiple assays can be imported", {
tx2gene_dir <- file.path(tempdir(), 'tx2gene')
dir.create(tx2gene_dir)

mock_uploaded_data(dataset_name,
sc_dir,
uploaded_data_dir,
type = 'h5',
version = '3',
gene.type = c(rep('Gene Expression', 199), 'Antibody Capture'))

expect_error(
suppressWarnings(import_scseq(dataset_name, uploaded_data_dir, sc_dir, tx2gene_dir)),
NA
Expand All @@ -104,3 +97,133 @@ test_that("cellranger .h5 files with multiple assays can be imported", {
# cleanup
unlink(c(sc_dir, uploaded_data_dir, tx2gene_dir), recursive = TRUE)
})


test_that("Seurat files with RNA Assay5 can be imported", {
# setup
dataset_name <- 'test'
sc_dir <- file.path(tempdir(), 'single-cell')
dir.create(sc_dir)
uploaded_data_dir <- file.path(tempdir(), 'uploads')
dir.create(uploaded_data_dir)

tx2gene_dir <- file.path(tempdir(), 'tx2gene')
dir.create(tx2gene_dir)

# create Seurat object
pbmc_raw <- read.table(
file = system.file("extdata", "pbmc_raw.txt", package = "Seurat"),
as.is = TRUE
)

# can upload non-split Assay5
pbmc_raw <- as(as.matrix(pbmc_raw), 'dgCMatrix')
scdata <- Seurat::CreateSeuratObject(counts = pbmc_raw)
expect_true(methods::is(scdata[['RNA']], 'Assay5'))


qs::qsave(scdata, file.path(uploaded_data_dir, 'scdata.qs'))

expect_error(
suppressWarnings(import_scseq(
dataset_name,
species = 'Homo sapiens',
uploaded_data_dir,
sc_dir,
tx2gene_dir,
metrics = NULL)),
NA
)

# cleanup
unlink(c(sc_dir, uploaded_data_dir, tx2gene_dir), recursive = TRUE)
})

test_that("Seurat files with split RNA Assay5 can be imported", {
# setup
dataset_name <- 'test'
sc_dir <- file.path(tempdir(), 'single-cell')
dir.create(sc_dir)
uploaded_data_dir <- file.path(tempdir(), 'uploads')
dir.create(uploaded_data_dir)

tx2gene_dir <- file.path(tempdir(), 'tx2gene')
dir.create(tx2gene_dir)

# create Seurat object
pbmc_raw <- read.table(
file = system.file("extdata", "pbmc_raw.txt", package = "Seurat"),
as.is = TRUE
)

# can upload non-split Assay5
pbmc_raw <- as(as.matrix(pbmc_raw), 'dgCMatrix')
scdata <- Seurat::CreateSeuratObject(counts = pbmc_raw)
expect_true(methods::is(scdata[['RNA']], 'Assay5'))

scdata$batch <- c('a', 'b', 'c')
scdata[["RNA"]] <- split(scdata[["RNA"]], f = scdata$batch)

expect_equal(SeuratObject::Layers(scdata), c("counts.a", "counts.b", "counts.c"))


qs::qsave(scdata, file.path(uploaded_data_dir, 'scdata.qs'))

expect_error(
suppressWarnings(import_scseq(
dataset_name,
species = 'Homo sapiens',
uploaded_data_dir,
sc_dir,
tx2gene_dir,
metrics = NULL)),
NA
)

# cleanup
unlink(c(sc_dir, uploaded_data_dir, tx2gene_dir), recursive = TRUE)
})


test_that("Seurat object with RNA Assay v3 can be imported", {

options(Seurat.object.assay.version = "v3")
# setup
dataset_name <- 'test'
sc_dir <- file.path(tempdir(), 'single-cell')
dir.create(sc_dir)
uploaded_data_dir <- file.path(tempdir(), 'uploads')
dir.create(uploaded_data_dir)

tx2gene_dir <- file.path(tempdir(), 'tx2gene')
dir.create(tx2gene_dir)

# create Seurat object
pbmc_raw <- read.table(
file = system.file("extdata", "pbmc_raw.txt", package = "Seurat"),
as.is = TRUE
)

# can upload Assay
pbmc_raw <- as(as.matrix(pbmc_raw), 'dgCMatrix')
scdata <- Seurat::CreateSeuratObject(counts = pbmc_raw)
expect_true(methods::is(scdata[['RNA']], 'Assay'))


qs::qsave(scdata, file.path(uploaded_data_dir, 'scdata.qs'))

expect_error(
suppressWarnings(import_scseq(
dataset_name,
species = 'Homo sapiens',
uploaded_data_dir,
sc_dir,
tx2gene_dir,
metrics = NULL)),
NA
)

# cleanup
options(Seurat.object.assay.version = "v5")
unlink(c(sc_dir, uploaded_data_dir, tx2gene_dir), recursive = TRUE)
})

0 comments on commit 2b1fe59

Please sign in to comment.