graph_test | R Documentation |
We are often interested in finding genes that are differentially expressed across a single-cell trajectory. Monocle3 introduces a new approach for finding such genes that draws on a powerful technique in spatial correlation analysis, the Moran’s I test. Moran’s I is a measure of multi-directional and multi-dimensional spatial autocorrelation. The statistic tells you whether cells at nearby positions on a trajectory will have similar (or dissimilar) expression levels for the gene being tested. Although both Pearson correlation and Moran’s I ranges from -1 to 1, the interpretation of Moran’s I is slightly different: +1 means that nearby cells will have perfectly similar expression; 0 represents no correlation, and -1 means that neighboring cells will be anti-correlated.
graph_test(
cds,
neighbor_graph = c("knn", "principal_graph"),
reduction_method = "UMAP",
k = 25,
method = c("Moran_I"),
alternative = "greater",
expression_family = "quasipoisson",
cores = 1,
verbose = FALSE,
nn_control = list()
)
cds |
a cell_data_set object upon which to perform this operation |
neighbor_graph |
String indicating what neighbor graph to use. "principal_graph" and "knn" are supported. Default is "knn", but "principal_graph" is recommended for trajectory analysis. |
reduction_method |
character, the method used to reduce dimension. Currently only supported for "UMAP". |
k |
Number of nearest neighbors used for building the kNN graph which is passed to knn2nb function during the Moran's I (Geary's C) test procedure. |
method |
a character string specifying the method (currently only 'Moran_I' is supported) for detecting significant genes showing correlation along the principal graph embedded in the low dimensional space. |
alternative |
a character string specifying the alternative hypothesis, must be one of greater (default), less or two.sided. |
expression_family |
a character string specifying the expression family function used for the test. |
cores |
the number of cores to be used while testing each gene for differential expression. |
verbose |
Whether to show spatial test (Moran's I) errors and warnings. Only valid for cores = 1. |
nn_control |
An optional list of parameters used to make the nearest neighbor index. See the set_nn_control help for detailed information. |
a data frame containing the p values and q-values from the Moran's I test on the parallel arrays of models.
moran.test
geary.test
expression_matrix <- readRDS(system.file('extdata',
'worm_l2/worm_l2_expression_matrix.rds',
package='monocle3'))
cell_metadata <- readRDS(system.file('extdata',
'worm_l2/worm_l2_coldata.rds',
package='monocle3'))
gene_metadata <- readRDS(system.file('extdata',
'worm_l2/worm_l2_rowdata.rds',
package='monocle3'))
cds <- new_cell_data_set(expression_data=expression_matrix,
cell_metadata=cell_metadata,
gene_metadata=gene_metadata)
cds <- preprocess_cds(cds, num_dim = 100)
cds <- reduce_dimension(cds)
cds <- cluster_cells(cds, resolution=1e-5)
colData(cds)$assigned_cell_type <- as.character(partitions(cds))
colData(cds)$assigned_cell_type <- dplyr::recode(colData(cds)$assigned_cell_type,
"1"="Germline",
"2"="Body wall muscle",
"3"="Unclassified neurons",
"4"="Vulval precursors",
"5"="Failed QC",
"6"="Seam cells",
"7"="Pharyngeal epithelia",
"8"="Coelomocytes",
"9"="Am/PH sheath cells",
"10"="Failed QC",
"11"="Touch receptor neurons",
"12"="Intestinal/rectal muscle",
"13"="Pharyngeal neurons",
"14"="NA",
"15"="flp-1(+) interneurons",
"16"="Canal associated neurons",
"17"="Ciliated sensory neurons",
"18"="Other interneurons",
"19"="Pharyngeal gland",
"20"="Failed QC",
"21"="Ciliated sensory neurons",
"22"="Oxygen sensory neurons",
"23"="Ciliated sensory neurons",
"24"="Ciliated sensory neurons",
"25"="Ciliated sensory neurons",
"26"="Ciliated sensory neurons",
"27"="Oxygen sensory neurons",
"28"="Ciliated sensory neurons",
"29"="Unclassified neurons",
"30"="Socket cells",
"31"="Failed QC",
"32"="Pharyngeal gland",
"33"="Ciliated sensory neurons",
"34"="Ciliated sensory neurons",
"35"="Ciliated sensory neurons",
"36"="Failed QC",
"37"="Ciliated sensory neurons",
"38"="Pharyngeal muscle")
neurons_cds <- cds[,grepl("neurons", colData(cds)$assigned_cell_type, ignore.case=TRUE)]
pr_graph_test_res <- graph_test(cds, neighbor_graph="knn")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.