## ----include = FALSE---------------------------------------------------------- knitr::opts_chunk$set( collapse = TRUE, comment = "#>", out.width = "100%" ) ## ----installation, eval=FALSE------------------------------------------------- # if (!require("BiocManager", quietly = TRUE)) # install.packages("BiocManager") # BiocManager::install("spacexr") ## ----setup, message = FALSE--------------------------------------------------- library(ggplot2) library(SpatialExperiment) library(SummarizedExperiment) library(spacexr) ## ----reference---------------------------------------------------------------- # Load simulated data. data("rctdSim") # Create SummarizedExperiment. reference_se <- SummarizedExperiment( assays = list(counts = rctdSim$reference_counts), colData = rctdSim$reference_cell_types ) # Examine reference. (optional) print(dim(assay(reference_se))) # Gene expression matrix dimensions table(colData(reference_se)$cell_type) # Number of occurrences of each cell type ## ----spatial, fig.align='center', fig.width=5, fig.height=5, out.width='50%'---- # Create SpatialExperiment. spatial_spe <- SpatialExperiment( assay = rctdSim$spatial_rna_counts, spatialCoords = rctdSim$spatial_rna_coords ) # Examine pixels. (optional) print(dim(assay(spatial_spe))) # Gene expression matrix dimensions ggplot(spatialCoords(spatial_spe), aes(x, y)) + geom_point(size = 5) + theme(axis.text = element_blank(), axis.ticks = element_blank() ) ## ----trueProportions---------------------------------------------------------- unique(t(assay(rctdSim$proportions_spe))) ## ----trueProportionsPlot, fig.height = 5, fig.width = 7, out.width='75%'------ plotAllWeights( rctdSim$proportions_spe, r = 0.05, lwd = 0, title = "Ground Truth" ) ## ----createRctd, message = FALSE---------------------------------------------- # Preprocess data rctd_data <- createRctd(spatial_spe, reference_se) ## ----runRctd, message = FALSE------------------------------------------------- results_spe <- runRctd(rctd_data, rctd_mode = "doublet", max_cores = 1) ## ----doubletWeights----------------------------------------------------------- # Cell type mixture 1 assay(results_spe, "weights")[, 1:6] # Cell type mixture 2 assay(results_spe, "weights")[, 7:12] ## ----classification----------------------------------------------------------- classification_df <- data.frame( pixel = colnames(assay(results_spe)), spot_class = colData(results_spe)$spot_class, first_type = colData(results_spe)$first_type, second_type = colData(results_spe)$second_type ) classification_df ## ----doubletVisualization, fig.height = 5, fig.width = 7, out.width='75%'----- plotAllWeights(results_spe, r = 0.05, lwd = 0, title = "Doublet Results") ## ----fullVisualization, fig.height = 5, fig.width = 7, out.width='75%'-------- plotAllWeights( results_spe, "weights_full", r = 0.05, lwd = 0, title = "Full Results" ) ## ----densityVisualization, fig.height = 5, fig.width = 7, out.width='75%'----- plotCellTypeWeight( results_spe, "ct1", "weights_full", size = 5, stroke = 0.5, title = "\nDensity of Cell Type 1\n" ) ## ----sessionInfo-------------------------------------------------------------- sessionInfo()