Main function for MAGIC imputation of single-cell data.
Performs MAGIC (Markov Affinity-based Graph Imputation of Cells) for denoising and imputation of single-cell data.
Run MAGIC on a Seurat object and store results in a new assay.
Run MAGIC on a SingleCellExperiment object.
Usage
magic(data, ...)
# Default S3 method
magic(
data,
genes = NULL,
knn = 5,
knn_max = NULL,
decay = 1,
t = 3,
npca = 100,
solver = c("exact", "approximate"),
init = NULL,
t_max = 20,
knn_dist = c("euclidean", "cosine", "manhattan", "correlation"),
n_jobs = 1,
seed = NULL,
verbose = TRUE,
...
)
# S3 method for class 'matrix'
magic(
data,
genes = NULL,
knn = 5,
knn_max = NULL,
decay = 1,
t = 3,
npca = 100,
solver = c("exact", "approximate"),
init = NULL,
t_max = 20,
knn_dist = c("euclidean", "cosine", "manhattan", "correlation"),
n_jobs = 1,
seed = NULL,
verbose = TRUE,
...
)
# S3 method for class 'dgCMatrix'
magic(data, ...)
# S3 method for class 'data.frame'
magic(data, ...)
# S3 method for class 'Seurat'
magic(
data,
assay = NULL,
slot = "data",
genes = "all_genes",
new_assay_name = NULL,
knn = 5,
knn_max = NULL,
decay = 1,
t = 3,
npca = 100,
solver = c("exact", "approximate"),
init = NULL,
t_max = 20,
knn_dist = c("euclidean", "cosine", "manhattan", "correlation"),
n_jobs = 1,
seed = NULL,
verbose = TRUE,
...
)
# S3 method for class 'SingleCellExperiment'
magic(
data,
assay_name = "logcounts",
genes = "all_genes",
new_assay_name = "MAGIC",
knn = 5,
knn_max = NULL,
decay = 1,
t = 3,
npca = 100,
solver = c("exact", "approximate"),
init = NULL,
t_max = 20,
knn_dist = c("euclidean", "cosine", "manhattan", "correlation"),
n_jobs = 1,
seed = NULL,
verbose = TRUE,
...
)Arguments
- data
A SingleCellExperiment object.
- ...
Additional arguments.
- genes
Gene selection. Default is "all_genes".
- knn
Number of neighbors. Default is 5.
- knn_max
Maximum neighbors. Default is NULL.
- decay
Alpha decay. Default is 1.
- t
Diffusion steps. Default is 3.
- npca
PCA components. Default is 100.
- solver
Solver type. Default is "exact".
- init
Previous result. Default is NULL.
- t_max
Maximum t. Default is 20.
- knn_dist
Distance metric. Default is "euclidean".
- n_jobs
Parallel workers. Default is 1.
- seed
Random seed. Default is NULL.
- verbose
Verbosity. Default is TRUE.
- assay
Character. Name of the assay to use as input. Default is the default assay of the object.
- slot
Character. Slot to use for input data. Default is "data" (normalized data). For Seurat v5, this corresponds to the layer name.
- new_assay_name
Name for the new assay. Default is "MAGIC".
- assay_name
Character. Name of the assay to use. Default is "logcounts".
Value
A "magic" object (for matrices) or modified input object (for Seurat/SCE).
The Seurat object with a new assay containing MAGIC-imputed values. The new assay is named "MAGIC_<original_assay>" by default.
SingleCellExperiment object with MAGIC results in a new assay.
Details
MAGIC denoises single-cell data using diffusion on a cell similarity graph:
Optional PCA dimension reduction
Build k-nearest neighbor graph
Construct alpha-decaying kernel
Create row-stochastic diffusion operator
Apply diffusion operator t times
The imputation formula is: \(X_{imputed} = P^t X\)
This function:
Extracts normalized data from the specified assay
Transposes to cells x genes format
Runs MAGIC imputation
Creates a new assay with the imputed data
Stores MAGIC parameters in the object's misc slot
The function is compatible with both Seurat v4 and v5. For v4, it uses GetAssayData/SetAssayData. For v5, it uses LayerData where available.
References
van Dijk, D., et al. (2018). Recovering gene interactions from single-cell data using data diffusion. Cell, 174(3), 716-729. doi:10.1016/j.cell.2018.05.061
Examples
if (FALSE) { # \dontrun{
# Basic usage
set.seed(42)
data <- matrix(rpois(5000, 2), nrow = 500, ncol = 10)
result <- magic(data, t = 3)
imputed <- as.matrix(result)
# Automatic t selection
result <- magic(data, t = "auto")
# Specific genes
result <- magic(data, genes = c(1, 2, 3))
} # }
if (FALSE) { # \dontrun{
library(Seurat)
# Create example Seurat object
counts <- matrix(rpois(3000, lambda = 5), nrow = 30, ncol = 100)
rownames(counts) <- paste0("Gene", 1:30)
colnames(counts) <- paste0("Cell", 1:100)
seurat_obj <- CreateSeuratObject(counts = counts)
seurat_obj <- NormalizeData(seurat_obj)
# Run MAGIC
seurat_obj <- magic(seurat_obj, t = 3)
# Access MAGIC results
GetAssayData(seurat_obj, assay = "MAGIC_RNA", slot = "data")
# Make MAGIC assay the default
DefaultAssay(seurat_obj) <- "MAGIC_RNA"
} # }