A step-by-step guide to learning the core functionalities of the ACTIONet framework

Here we provide an in-depth, step-by-step guide to familiarize you with ACTIONet building blocks. ACTIONet framework, at its core, consists of 8 major components/steps, which are depicted below. In what follows, We will illustrate and test the use case of each functionality. For the sake of compatibility, we will use the PBMC 3k dataset used in the Seurat and Scanpy for this tutorial. All test datasets are available for download from Datasets

The overall ACTIONet workflow

Setup

Loading ACTIONet package

To start, you need to load the ACTIONet package. For the ACTIONet installation instruction, see ACTIONet repository.

require(ACTIONet)

input_path = "./input/PBMC3k/"

Importing data

We have preprocessed and stored the PBMC 3k dataset as an SCE object, which is one of the main datatype used as input to the ACTIONet framework. For more information on how to convert different datatypes to SingleCellExperiment format, please consult our tutorial on PortingData: Import/export options in the ACTIONet framework.

ace = readRDS(file.path(input_path, "input_sce.RDS"))

ace
class: SingleCellExperiment 
dim: 1838 2638 
metadata(0):
assays(1): logcounts
rownames(1838): TNFRSF4 CPSF3L ... S100B PRMT2
rowData names(1): n_cells
colnames(2638): AAACATACAACCAC-1 AAACATTGAGCTAC-1 ... TTTGCATGAGAGGC-1
  TTTGCATGCCTCAC-1
colData names(4): n_genes percent_mito n_counts louvain
reducedDimNames(4): draw_graph_fr pca tsne umap
spikeNames(0):

Reducing the ACE object

Next, we use the reduce.ace() function to both normalize the counts’ matrix (lib-size normalization + log-transformation, by default) and compute a factorized (reduced) form of the kernel matrix (50 by default). By default, no batch correction is performed, and the simplest depth-normalization method is used. We then save the output for future use.

dim = 50

ace = reduce.ace(ace = ace, reduced_dim = dim)

S = logcounts(ace)  # Normalized logcount matrix
S_r = colMaps(ace)[["ACTION"]]  # Reduced kernel matrix

PS: ACTIONet supports additional normalization methods, including linnorm and scater However, the default method is the fastest and performs reasonably well. It is possible to renormalize the counts matrix afterwards.

PS: ACTIONet also seamlessly integrates with Harmony and MNN batch-correction method. For more information on how to use batch correction please consult To batch correct or not to batch correct, that is the question!.

Creating an ACTIONetExperiment (ACE) object to hold results

We have developed a new class, ACTIONetExperiment, that extends the SingleCellExperiment class to incorporate additional slots to store multidimensional and structural meta-data for rows and columns. The motivation is to merge benefits of both SingleCellExperiment and AnnData into the same class and provide interoperability between R and Python. For more information on the ACTIONetExperiment object, please visit Introducing ACTIONetExperiment(ACE) Object – Ins And Outs.

As it can be seen, the new slots are rowNets, colNets, rowMaps, and colMaps. Each slot has its on getter/setter functions, which we will use through out the rest of the tutorial to store results.


Running ACTIONet

Running multi-level ACTION decomposition

This will run ACTIONet with increasing number of archetypes up to a max of depth. Depth of 30 works well in most cases, but you can decrease this number in favor of speed, or increase it in favor of identifying rare cell types.

depth = 30

ACTION.out = run_ACTION(S_r, k_max = depth)

Prune nonspecific and/or unreliable archetypes

To remove unreliable archetypes as a first pass filtering, we use prune_archetypes() function. min_cell parameter indicates the minimum number of “influential cells” for an archetype to be retained. If an archetype is defined solely by a single cell, it is likely that is it unreliable. Setting min_cell=2 removes such archetypes. This value is appropriate for the majority of cases and is provided here for illustration purposes only.

min_cells = 5

pruning.out = prune_archetypes(ACTION.out$C, ACTION.out$H, min_cells = min_cells)


C_stacked = pruning.out$C_stacked
H_stacked = pruning.out$H_stacked


colMaps(ace)[["C_stacked"]] = Matrix::t(C_stacked)
colMaps(ace)[["H_stacked"]] = H_stacked

colMapTypes(ace)[["C_stacked"]] = "internal"
colMapTypes(ace)[["H_stacked"]] = "internal"

C_stacked and H_stacked are the concatenated C and H matrices across different levels, respectively, after pruning noisy archetypes.

Building ACTIONet graph

To build ACTIONet, we use the approximate nearest neighbor search HNSW algorithm, as implemented in hnswlib. We modified the library to incorporate ideas from k*-NN to construct a density-dependent nearest-neighbor graph. Density parameter here is inversely proportional to the LC parameter in the k*-NN algorithm.

density = 1

set.seed(0)
G = build_ACTIONet(H_stacked = H_stacked, density = density)

To store computed graph as a network associated with cells (columns) in the ACE object, we can use the colNets of the ACE object. Column-networks, or colNets, are networks associated with columns of ace object (which are cells here) and directly map to obsp field of AnnData objects. We are adding one network, calling ACTIONet, and save the adjacent matrix ‘G’ inside ACTIONet (which is inside ace). (You can add networks as many as you want, and here we are adding one network.)

colNets(ace)$ACTIONet = G

Layout ACTIONet

We have adopted and modified the stochastic-gradient descent (SGD)-based layout algorithm utilized in the UMAP for our visualization, but modified it to use ACTIONet network as input instead of constructing a kNN graph. We also use ACTION reduction (S_r) for our initialization. Finally, we siplified all parameters into a compactness value can be between 1 and 100 and indicates the compactness of the layout.

compactness = 50
# initial.coordinates = t(scale(t(S_r)))

vis.out = layout_ACTIONet(G, S_r = S_r, compactness_level = compactness, thread_no = 1)

We then store the coordinates and colors as colMaps() of the ACE object

X = Matrix::t(vis.out$coordinates)
rownames(X) = c("x", "y")
colnames(X) = colnames(ace)
colMaps(ace)$ACTIONet2D = X
colMapTypes(ace)[["ACTIONet2D"]] = "embedding"

X = Matrix::t(vis.out$coordinates_3D)
rownames(X) = c("x", "y", "z")
colnames(X) = colnames(ace)
colMaps(ace)$ACTIONet3D = X
colMapTypes(ace)[["ACTIONet3D"]] = "embedding"

X = Matrix::t(vis.out$colors)
rownames(X) = c("r", "g", "b")
colnames(X) = colnames(ace)
colMaps(ace)$denovo_color = X
colMapTypes(ace)[["denovo_color"]] = "embedding"

Identify equivalent classes of archetypes and group them together

There is a large redundancy between the set of archetypes due to repeated decomposition using increasing values of archetypes. unify_archetypes() tries to coalesce these archetypes and partition them into equivalent classes. Resolution adjust the final number of archetypes retained (larger resolution, more archetypes):

resolution = 1

unification.out = unify_archetypes(S_r, pruning.out$C_stacked, pruning.out$H_stacked, 
    resolution = resolution)

And store the results back in the ACE object:

colMaps(ace)[["H_unified"]] = as(unification.out$H_unified, "sparseMatrix")
colMaps(ace)[["C_unified"]] = as(Matrix::t(unification.out$C_unified), "sparseMatrix")

colMapTypes(ace)[["H_unified"]] = "internal"
colMapTypes(ace)[["C_unified"]] = "internal"

ace$assigned_archetype = unification.out$assigned_archetype

Use graph core of induced subgraphs by archetypes to infer centrality/quality of each cell

We have developed a novel technique to estimate the quality of cells based on the strategic placement in the ACTIONet, w.r.t to each archetype. For each archetype, we identify the subgraph induced by its most closely associated cells. Then, within each subgraph, we compute coreness of vertices in the network and smooth it out using a localized network diffusion:

ace$node_centrality = compute_archetype_core_centrality(G, ace$assigned_archetype)

Compute gene specificity for each archetype

This computes the “marker-ness” of each genes for different archetypes. The key idea of it is to compare the mean observed expression of a gene under distribuion of an archetype across all cells, compared to the uniform distribution, while accounting to the sparsity of genes/cells:

S = logcounts(ace)

# We use compute_archetype_feature_specificity() and
# compute_archetype_feature_specificity_full() for sparse and dense S matrices,
# respectively.
specificity.out = compute_archetype_feature_specificity_full(S, unification.out$H_unified)
specificity.out = lapply(specificity.out, function(specificity.scores) {
    rownames(specificity.scores) = rownames(ace)
    colnames(specificity.scores) = paste("A", 1:ncol(specificity.scores))
    return(specificity.scores)
})
rowMaps(ace)[["H_unified_profile"]] = specificity.out[["archetypes"]]
rowMaps(ace)[["H_unified_upper_significance"]] = specificity.out[["upper_significance"]]
rowMaps(ace)[["H_unified_lower_significance"]] = specificity.out[["lower_significance"]]

rowMapTypes(ace)[["H_unified_profile"]] = "internal"
rowMapTypes(ace)[["H_unified_upper_significance"]] = "internal"
rowMapTypes(ace)[["H_unified_lower_significance"]] = "internal"

Creating a trace of all that we have done

trace = list(ACTION.out = ACTION.out, pruning.out = pruning.out, vis.out = vis.out, 
    unification.out = unification.out)
trace$log = list(genes = rownames(ace), cells = colnames(ace), time = Sys.time())

Visualization

The main function to visualize results is plot.ACTIONet(). It takes the resulting ACE object. By default, it uses ACTIONet2D coordinate_slot, but we can provide any arbitrary embedding to it. Also, in addition to colors, we use transparency to emphasize different characteristics of cells, including estimated quality scores.

ACTIONet embedding

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "ACTIONet2D", title = "ACTIONet")

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "ACTIONet2D", title = "ACTIONet (with node focusing)", 
    transparency.attr = ace$node_centrality)

Comparison with other embeddings

And now we plot other embeddings for comparison purposes:

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "pca", title = "pca")

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "draw_graph_fr", title = "ForceDirected")

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "tsne", title = "tSNE")

plot.ACTIONet(ace, ace$louvain, coordinate_slot = "umap", title = "UMAP")

Shahin Mohammadi

Jose Davila-Velderrain

Sebastian Pineda

2020-06-23