--- title: "PANORAMIC tutorial" output: BiocStyle::html_document vignette: > %\VignetteIndexEntry{PANORAMIC tutorial} %\VignetteEncoding{UTF-8} %\VignetteEngine{knitr::rmarkdown} editor_options: markdown: wrap: 72 --- ```{r setup, include=FALSE} knitr::opts_chunk$set(collapse = TRUE, comment = "#>") ``` # Introduction PANORAMIC is designed for *multi-sample* spatial colocalization analysis. It estimates sample-level spatial effects, then pools them with multilevel meta-analysis to test group-level differences. Compared with single-sample analyses, PANORAMIC gives you: - uncertainty-aware pooling across samples - explicit between-sample / between-patient heterogeneity modeling - direct two-group contrasts (`beta_diff`, `p_diff`, `fdr_diff`) - consistent table and plotting utilities for downstream interpretation This vignette focuses on practical usage: 1. Build a toy multi-sample dataset. 2. Run PANORAMIC stepwise (prepare -\> spatial stats -\> pooling). 3. Run the one-call workflow helper (`panoramic_analyze()`). 4. Interpret outputs and generate publication-style plots. # Installation ```{r installation, eval=FALSE} if (!requireNamespace("BiocManager", quietly = TRUE)) { install.packages("BiocManager") } BiocManager::install("panoramic") ``` # Load packages ```{r libraries, message=FALSE, warning=FALSE} library(panoramic) library(SummarizedExperiment) library(BiocParallel) library(dplyr) library(ggplot2) ``` # Workflow At A Glance | Step | Function | Primary output | Why it matters | |:------------|:---------------------------|:------------------------------------------------------------------|:-----------------------------------------| | 1 | `panoramic_prepare()` | list of prepared `SpatialExperiment`s with cached spatial objects | standardizes per-sample inputs | | 2 | `panoramic_spatialstats()` | `SummarizedExperiment` with sample-level `yi` / `vi` | quantifies within-sample spatial effects | | 3 | `panoramic_meta_mv()` | pooled multilevel results + differential contrasts | tests group-level biology | | Convenience | `panoramic_analyze()` | list with `prep`, `stats`, `pooled`, `tables` | one-call reproducible pipeline | # Simulate Example Data We simulate two groups of samples with different spatial structure: - `group1`: mostly random cell mixing - `group2`: stronger local colocalization for early cell types ```{r simulate_data, message=FALSE, warning=FALSE} set.seed(123) toy <- panoramic_simulate_dataset( n_group1 = 3, n_group2 = 3, n_cells_group1 = 200, n_cells_group2 = 350, group_labels = c("group1", "group2"), scenario_group1 = "random", scenario_group2 = "colocalized", seed = 123 ) spe_list <- toy$spe_list design <- toy$design design ``` Quick geometric sanity check: ```{r plot_samples, message=FALSE, warning=FALSE} plot_df <- do.call(rbind, lapply(spe_list, function(spe) { coords <- SpatialExperiment::spatialCoords(spe) data.frame( x = coords[, 1], y = coords[, 2], cell_type = SummarizedExperiment::colData(spe)$cell_type, sample_id = SummarizedExperiment::colData(spe)$sample_id, stringsAsFactors = FALSE ) })) ggplot(plot_df, aes(x = x, y = y, color = cell_type)) + geom_point(size = 0.7, alpha = 0.7) + facet_wrap(~ sample_id, nrow = 2) + coord_equal() + theme_bw() + labs(x = "x", y = "y", color = "Cell type") ``` # Step 1: Prepare Samples ```{r prepare, message=FALSE, warning=FALSE} prep <- panoramic_prepare( spe_list = spe_list, design = design, cell_type = "cell_type", min_cells = 5, window = "concave" ) names(S4Vectors::metadata(prep[[1]])$panoramic) ``` # Step 2: Compute Sample-Level Spatial Effects Default PANORAMIC statistic is `stat = "local_comp_enrichment"` (percentage-point enrichment vs local null expectation). ```{r spatial_stats, message=FALSE, warning=FALSE} radii_um <- c(25) se <- panoramic_spatialstats( prep = prep, radii_um = radii_um, stat = "local_comp_enrichment", nsim = 30, seed = 123, BPPARAM = BiocParallel::SerialParam() ) dim(se) head(as.data.frame(rowData(se))[, c("ct1", "ct2", "radius_um", "stat")]) ``` # Step 3: Pool Across Samples And Test Group Differences `panoramic_meta_mv()` pools sample-level effects and computes differential contrasts between groups. ```{r meta_mv, message=FALSE, warning=FALSE} # Here patient_col and sample_col are both "sample" in toy data. # In real data, patient_col can differ from sample_col (e.g., multiple cores per patient). se_meta <- panoramic_meta_mv( se, patient_col = "sample", group_col = "group", sample_col = "sample", tau_structure = "patient", method = "REML", group_tau2 = "separate" ) ``` # Interpret Pooled And Differential Outputs Each row in `rowData(se_meta)` is one feature (`ct1`, `ct2`, `radius_um`) with pooled effects and contrasts. For differential terms: - positive `beta_diff`: stronger spatial effect in `group2` vs `group1` - negative `beta_diff`: weaker spatial effect in `group2` - small `fdr_diff`: stronger evidence after multiplicity correction ```{r extract_tables, message=FALSE, warning=FALSE} res <- panoramic_extract_contrast(se_meta) head(res[, c("ct1", "ct2", "radius_um", "beta_diff", "p_diff", "fdr_diff")]) res %>% dplyr::filter(fdr_diff < 0.05) ``` # One-Call Workflow (`panoramic_analyze`) Use `panoramic_analyze()` when you want a single reproducible entry point for end-to-end runs. ```{r one_call, message=FALSE, warning=FALSE} out <- panoramic_analyze( spe_list = spe_list, design = design, cell_type = "cell_type", radii_um = radii_um, nsim = 20, min_cells = 5, window = "concave", BPPARAM = BiocParallel::SerialParam() ) names(out) names(out$tables) ``` # Visualization Tools ## Forest Plot (Feature-Level) ```{r forest_plot, message=FALSE, warning=FALSE} plot_forest( se_meta, ct1 = "A", ct2 = "B", radius_um = 25, group_col = "group" ) ``` ## Volcano Plot (Global Differential Overview) ```{r volcano_plot, message=FALSE, warning=FALSE} plot_volcano(se_meta) ``` ## Network Summary ```{r network_plot, message=FALSE, warning=FALSE} net <- create_spatial_network( se_meta, fdr_threshold = 0.2, directed = FALSE, leiden_resolution = 1.0 ) plot_spatial_network( se_meta, fdr_threshold = 0.2, directed = FALSE, layout = "fr", node_size_by = "degree" ) ``` # Practical Notes - Start with `local_comp_enrichment` unless you have a specific reason to use L/K-function alternatives. - Use at least two biological samples per group for stable contrasts. - Increase `nsim` for final analyses (the vignette uses modest values for speed). - Keep sample metadata (`sample`, `group`, `patient`) explicit and consistent early in your workflow. # Related Resources For complementary spatial analysis workflows, see: - `r BiocStyle::CRANpkg("spatstat")` family (`r BiocStyle::CRANpkg("spatstat.geom")`, `r BiocStyle::CRANpkg("spatstat.explore")`) for point-process spatial statistics. - `r BiocStyle::Biocpkg("SPIAT")` for spatial histology and tumor microenvironment analyses. - `r BiocStyle::Biocpkg("imcRtools")` for imaging mass cytometry spatial workflows. - `r BiocStyle::Biocpkg("spicyR")` for differential spatial localization testing in imaging data. ```{r session_info} sessionInfo() ```