How to Test and Store Multiple Parameters or Analyses?¶
The default Giotto workflow is similar to other scRNA-seq workflows and does not require you to provide a custom name for each analysis (e.g. PCA, UMAP, …), but running an analysis twice will overwrite the previous results with a warning.
However, there are situations where being able to run and store multiple analyses can be advantageous:
Test multiple parameters for a single analysis
Test multiple combinations across functions (See Example: hvg->pca->umap)
Use different output results as input for downstream analyses (See Example: spatial genes)
data:image/s3,"s3://crabby-images/f5c45/f5c45fcbb8208576ae63839b8db6d6765aa4e325" alt="Multiple Analysis"
We will use the seqFish+ somatosensory cortex as an example dataset after creating and processing a Giotto object.
1. Calculate Highly Variable Genes (2 Methods)¶
# using the loess method
VC_test <- calculateHVG(gobject = VC_test,
method = 'cov_loess', difference_in_cov = 0.1,
HVGname = 'loess_hvg')
data:image/s3,"s3://crabby-images/87683/876835c42c4b5eebd51dd843c20677b6c5757414" alt="Loess"
# using the expression groups method
VC_test <- calculateHVG(gobject = VC_test
, method = 'cov_group', zscore_threshold = 1,
HVGname = 'group_hvg')
data:image/s3,"s3://crabby-images/78513/78513f04f6d0d6f5d1908f0c3baf63b926f49b42" alt="Group"
# compare the highly variable genes between two methods
gene_metadata = fDataDT(VC_test)
mytable = table(loess = gene_metadata$loess_hvg, group = gene_metadata$group_hvg)
data:image/s3,"s3://crabby-images/96254/9625459297fe9f28435ffa07a6220d4c2181ed42" alt="Group"
2. Perform Multiple PCAs¶
Using the 2 different HVG sets (loess_genes and group_genes)
Store PCA results using custom names (‘pca_loess’ and ‘pca_group’)
Plot PCA results
## 2. PCA ##
# pca with genes from loess
loess_genes = gene_metadata[loess_hvg == 'yes']$gene_ID
VC_test <- runPCA(gobject = VC_test, genes_to_use = loess_genes, name = 'pca_loess', scale_unit = F)
plotPCA(gobject = VC_test, dim_reduction_name = 'pca_loess')
data:image/s3,"s3://crabby-images/090ab/090ab413651b95f83c98a293940b1a040c3db772" alt="Group"
# pca with genes from group
group_genes = gene_metadata[group_hvg == 'yes']$gene_ID
VC_test <- runPCA(gobject = VC_test, genes_to_use = group_genes, name = 'pca_group', scale_unit = F)
plotPCA(gobject = VC_test, dim_reduction_name = 'pca_group')
data:image/s3,"s3://crabby-images/83522/83522104829b511e1d3ecaffe0077937d5ed102e" alt="Group"
3. Create Multiple UMAPs¶
Using the 2 different PCA results (‘pca_loess’ and ‘pca_group’)
Store UMAP results using custom names (‘umap_loess’ and ‘umap_group’)
Plot UMAP results
## 3. UMAP ##
VC_test <- runUMAP(VC_test, dim_reduction_to_use = 'pca', dim_reduction_name = 'pca_loess',
name = 'umap_loess', dimensions_to_use = 1:30)
plotUMAP(gobject = VC_test, dim_reduction_name = 'umap_loess')
data:image/s3,"s3://crabby-images/a6b92/a6b925b68cd1b27937a02652d09c6712d6dd98f8" alt="Group"
VC_test <- runUMAP(VC_test, dim_reduction_to_use = 'pca', dim_reduction_name = 'pca_group',
name = 'umap_group', dimensions_to_use = 1:30)
plotUMAP(gobject = VC_test, dim_reduction_name = 'umap_group')
data:image/s3,"s3://crabby-images/cb20a/cb20a5c32f664124d627c88a7f976014e29802ef" alt="Group"
4. Create Multiple Spatial Networks¶
Create spatial with multiple k’s and other parameters (k=5, k=10, k=100 & maximum_distance=200)
Subset field 1
Visualize network on field 1 (‘spatial_network’, ‘large_network’, ‘distance_work’)
Spatial Network¶
## 4. spatial network
VC_test <- createSpatialNetwork(gobject = VC_test, method = 'kNN', k = 5) # standard name: 'spatial_network'
VC_test <- createSpatialNetwork(gobject = VC_test, method = 'kNN', k = 10, name = 'large_network') VC_test <- createSpatialNetwork(gobject = VC_test, method = 'kNN', k = 100, maximum_distance_knn = 200, minimum_k = 2, name = 'distance_network')
## visualize different spatial networks on first field (~ layer 1)
cell_metadata = pDataDT(VC_test)
field1_ids = cell_metadata[Field_of_View == 0]$cell_ID
subVC_test = subsetGiotto(VC_test, cell_ids = field1_ids)
spatPlot(gobject = subVC_test, show_network = T,
network_color = 'blue', spatial_network_name = 'spatial_network')
data:image/s3,"s3://crabby-images/7ca0b/7ca0beee70649a062da2e012296cbffac48530c6" alt="Group"
Large Network¶
spatPlot(gobject = subVC_test, show_network = T,
network_color = 'blue', spatial_network_name = 'large_network')
data:image/s3,"s3://crabby-images/62855/62855fa1aea324f39dd6de1e20119c3d0f47dd3a" alt="Group"
Distance Network¶
spatPlot(gobject = subVC_test, show_network = T,
network_color = 'blue', spatial_network_name = 'distance_network')
data:image/s3,"s3://crabby-images/e0472/e04728b6de4346c000f0f9510e6f02fc65f0ba51" alt="Group"
5. Find Spatial Genes (Multiple Methods)¶
Use the different spatial networks as input to identify spatial genes with the rank method
Visualize top spatial genes for 2 methods
Large Network Spatial Genes¶
## 5. spatial genes
# the provided spatial_network_name can be given to downstream analyses
# spatial genes based on large network
ranktest_large = binSpect(VC_test,
subset_genes = loess_genes,
bin_method = 'rank',
spatial_network_name = 'large_network')
spatGenePlot(VC_test,
expression_values = 'scaled',
genes = ranktest_large$genes[1:6], cow_n_col = 2, point_size = 1,
genes_high_color = 'red', genes_mid_color = 'white', genes_low_color = 'darkblue', midpoint = 0)
data:image/s3,"s3://crabby-images/21602/21602d749749a5de054f70033737388616299157" alt="Group"
Distance Network Spatial Genes¶
# spatial genes based on distance network
ranktest_dist = binSpect(VC_test,
subset_genes = loess_genes,
bin_method = 'rank',
spatial_network_name = 'distance_network')
spatGenePlot(VC_test,
expression_values = 'scaled',
genes = ranktest_dist$genes[1:6], cow_n_col = 2, point_size = 1,
genes_high_color = 'red', genes_mid_color = 'white', genes_low_color = 'darkblue', midpoint = 0)
data:image/s3,"s3://crabby-images/cbb6e/cbb6ebf08239a7096be319203e901e5f4caf05fd" alt="Group"