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:
beta_diff,
p_diff, fdr_diff)This vignette focuses on practical usage:
panoramic_analyze()).| Step | Function | Primary output | Why it matters |
|---|---|---|---|
| 1 | panoramic_prepare() |
list of prepared SpatialExperiments 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 |
We simulate two groups of samples with different spatial structure:
group1: mostly random cell mixinggroup2: stronger local colocalization for early cell
typesset.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
#> sample group
#> 1 sample_1 group1
#> 2 sample_2 group1
#> 3 sample_3 group1
#> 4 sample_4 group2
#> 5 sample_5 group2
#> 6 sample_6 group2Quick geometric sanity check:
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")Default PANORAMIC statistic is
stat = "local_comp_enrichment" (percentage-point enrichment
vs local null expectation).
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)
#> [1] 9 6
head(as.data.frame(rowData(se))[, c("ct1", "ct2", "radius_um", "stat")])
#> ct1 ct2 radius_um stat
#> A|A|25 A A 25 local_comp_enrichment
#> B|A|25 B A 25 local_comp_enrichment
#> C|A|25 C A 25 local_comp_enrichment
#> A|B|25 A B 25 local_comp_enrichment
#> B|B|25 B B 25 local_comp_enrichment
#> C|B|25 C B 25 local_comp_enrichmentpanoramic_meta_mv() pools sample-level effects and
computes differential contrasts between groups.
# 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"
)Each row in rowData(se_meta) is one feature
(ct1, ct2, radius_um) with pooled
effects and contrasts.
For differential terms:
beta_diff: stronger spatial effect in
group2 vs group1beta_diff: weaker spatial effect in
group2fdr_diff: stronger evidence after multiplicity
correctionres <- panoramic_extract_contrast(se_meta)
head(res[, c("ct1", "ct2", "radius_um", "beta_diff", "p_diff", "fdr_diff")])
#> ct1 ct2 radius_um beta_diff p_diff fdr_diff
#> 1 A A 25 13.475382 8.399035e-08 9.448914e-08
#> 2 B A 25 15.885177 3.371782e-34 1.011535e-33
#> 3 C A 25 9.730690 1.370860e-16 2.056290e-16
#> 4 A B 25 11.789955 1.187947e-24 2.672881e-24
#> 5 B B 25 12.239678 9.753484e-10 1.254019e-09
#> 6 C B 25 4.292444 5.390817e-04 5.390817e-04
res %>% dplyr::filter(fdr_diff < 0.05)
#> ct1 ct2 radius_um beta_diff se_diff z_diff p_diff fdr_diff
#> 1 A A 25 13.475382 2.5148430 5.358339 8.399035e-08 9.448914e-08
#> 2 B A 25 15.885177 1.3027716 12.193371 3.371782e-34 1.011535e-33
#> 3 C A 25 9.730690 1.1770168 8.267248 1.370860e-16 2.056290e-16
#> 4 A B 25 11.789955 1.1502807 10.249632 1.187947e-24 2.672881e-24
#> 5 B B 25 12.239678 2.0021089 6.113393 9.753484e-10 1.254019e-09
#> 6 C B 25 4.292444 1.2403950 3.460546 5.390817e-04 5.390817e-04
#> 7 A C 25 -10.596653 0.7602971 -13.937516 3.747556e-44 1.686400e-43
#> 8 B C 25 -12.582864 0.7201785 -17.471869 2.346806e-68 2.112126e-67
#> 9 C C 25 -8.802775 0.9407116 -9.357570 8.159165e-21 1.468650e-20
#> coloc_source coloc_target coloc_direction
#> 1 A A A -> A
#> 2 A B A -> B
#> 3 A C A -> C
#> 4 B A B -> A
#> 5 B B B -> B
#> 6 B C B -> C
#> 7 C A C -> A
#> 8 C B C -> B
#> 9 C C C -> Cpanoramic_analyze)Use panoramic_analyze() when you want a single
reproducible entry point for end-to-end runs.
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)
#> [1] "prep" "stats" "pooled" "tables"
names(out$tables)
#> [1] "spatialstats" "meta" "contrast"local_comp_enrichment unless you have a
specific reason to use L/K-function alternatives.nsim for final analyses (the vignette uses
modest values for speed).sample, group,
patient) explicit and consistent early in your
workflow.panoramic_analyze)