This vignette demonstrates how to use the
SpaceMarkersExperiment class, a
SpatialExperiment subclass that serves as a unified
container for all SpaceMarkers inputs and outputs. Instead of tracking
separate matrices, data frames, and lists across pipeline steps, a
single SpaceMarkersExperiment object holds everything:
expression data, spatial coordinates, pattern features, optimal
parameters, hotspots, interacting genes, and interaction scores.
We will reproduce the same breast cancer analysis from the main vignette using this integrated workflow.
data_dir <- "visiumBrCA"
unlink(file.path(data_dir), recursive = TRUE)
dir.create(data_dir, showWarnings = FALSE)
main_10xlink <- "https://cf.10xgenomics.com/samples/spatial-exp/1.3.0"
counts_folder <- "Visium_Human_Breast_Cancer"
counts_file <- "Visium_Human_Breast_Cancer_filtered_feature_bc_matrix.h5"
counts_url <- paste(c(main_10xlink, counts_folder, counts_file), collapse = "/")
sp_folder <- "Visium_Human_Breast_Cancer"
sp_file <- "Visium_Human_Breast_Cancer_spatial.tar.gz"
sp_url <- paste(c(main_10xlink, sp_folder, sp_file), collapse = "/")Downloads are routed through BiocFileCache so the ~260
MB dataset is fetched once and reused across vignette builds and across
the other SpaceMarkers vignettes.
bfc <- BiocFileCache::BiocFileCache(ask = FALSE)
counts_cache <- BiocFileCache::bfcrpath(bfc, counts_url)
sp_cache <- BiocFileCache::bfcrpath(bfc, sp_url)
file.copy(counts_cache, file.path(data_dir, counts_file), overwrite = TRUE)## [1] TRUE
load10Xload10X() searches the Visium directory for expression
and spatial files using flexible pattern matching — it works regardless
of whether your files are in an outs/, data/,
or flat directory layout. It log-transforms the counts and returns a
SpaceMarkersExperiment directly:
cogaps_result <- readRDS(system.file("extdata", "CoGAPS_result.rds",
package = "SpaceMarkers", mustWork = TRUE))
sme <- load10X(data_dir, features = cogaps_result, method = "CoGAPS")
sme## class: SpaceMarkersExperiment
## dim: 36601 4898
## metadata(0):
## assays(1): logcounts
## rownames(36601): MIR1302-2HG FAM138A ... AC007325.4 AC007325.2
## rowData names(0):
## colnames(4898): AAACAACGAATAGTTC-1 AAACAAGTATCTCCCA-1 ...
## TTGTTTGTATTACACG-1 TTGTTTGTGTAAATTC-1
## colData names(6): Pattern_1 Pattern_2 ... Pattern_5 sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : y x
## imgData names(0):
## spacemarkers analysis: none
## patterns: Pattern_1, Pattern_2, Pattern_3, Pattern_4, Pattern_5
Alternatively, you can load first and add features separately with
add_features():
# Step 1: Load without features
sme <- load10X(data_dir)
# Step 2: Add features later
sme <- add_features(sme, cogaps_result, method = "CoGAPS")If the features cover fewer spots than the expression data (or vice
versa), add_features() will keep only the shared spots and
emit a warning.
The object inherits from SpatialExperiment so all
standard accessors work:
## [1] 36601 4898
## y x
## AAACAACGAATAGTTC-1 4662 14376
## AAACAAGTATCTCCCA-1 16545 25987
## AAACAATCTACTAGCA-1 5392 18039
## AAACACCAATAACTGC-1 18606 14703
## AAACAGAGCGACTCCT-1 8032 24950
## AAACAGCTTTCAGAAG-1 14818 13367
## DataFrame with 6 rows and 5 columns
## Pattern_1 Pattern_2 Pattern_3 Pattern_4 Pattern_5
## <numeric> <numeric> <numeric> <numeric> <numeric>
## AAACAACGAATAGTTC-1 0.467625588 1.04939e-01 2.57606e-01 0.684806228 4.74709e-02
## AAACAAGTATCTCCCA-1 0.269075811 4.39442e-01 2.05647e-01 0.292133719 1.16758e-02
## AAACAATCTACTAGCA-1 0.110593386 1.14852e-02 2.30915e-01 0.411131471 9.50832e-02
## AAACACCAATAACTGC-1 0.000250838 1.68580e-01 1.22360e-01 0.000156279 8.04193e-01
## AAACAGAGCGACTCCT-1 0.284930885 1.10251e-01 9.05316e-08 0.242940620 3.43081e-08
## AAACAGCTTTCAGAAG-1 0.158373639 9.74108e-06 1.72347e-01 0.305995703 7.16760e-01
For this demonstration we will use two patterns and a curated gene list:
data("curated_genes")
good_gene_threshold <- 3
keep_genes <- rownames(sme)[
apply(assay(sme, "logcounts"), 1, function(x) sum(x > 0) >= good_gene_threshold)
]
keep_genes <- intersect(keep_genes, curated_genes)
sme <- sme[keep_genes, ]
# Keep all patterns for now; we will subset per analysis below
sme## class: SpaceMarkersExperiment
## dim: 114 4898
## metadata(0):
## assays(1): logcounts
## rownames(114): C1orf127 CDC42 ... ATRX UBE2A
## rowData names(0):
## colnames(4898): AAACAACGAATAGTTC-1 AAACAAGTATCTCCCA-1 ...
## TTGTTTGTATTACACG-1 TTGTTTGTGTAAATTC-1
## colData names(6): Pattern_1 Pattern_2 ... Pattern_5 sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : y x
## imgData names(0):
## spacemarkers analysis: none
## patterns: Pattern_1, Pattern_2, Pattern_3, Pattern_4, Pattern_5
The SpaceMarkers() function now accepts a
SpaceMarkersExperiment as its first argument. All results
are stored back in the returned object. We keep all five CoGAPS patterns
on the SME (so the overlap heatmap shows every pairwise relationship)
and restrict the gene-level interaction tests to Pattern_1 × Pattern_2 —
our narrative focus — via pattern_pairs, which
SpaceMarkers() forwards to
get_pairwise_interacting_genes().
sme_undirected <- SpaceMarkers(
x = sme,
directed = FALSE,
min.gene.expr = 3,
pattern_pairs = rbind(c("Pattern_1", "Pattern_2"))
)Summary scores are accessed via dedicated accessors:
## [1] "undirected"
## Gene Pattern_1_Pattern_2
## COMP COMP 5.115554
## MAP2K2 MAP2K2 4.808176
## HLA-A HLA-A 4.476319
## CST1 CST1 3.512456
## ARHGEF40 ARHGEF40 3.483833
## SLC12A9 SLC12A9 3.480979
Detailed intermediate results are stored in metadata()
and accessed via convenience functions:
## barcode y x Pattern_1 Pattern_2 Pattern_3
## AAACAACGAATAGTTC-1 AAACAACGAATAGTTC-1 4662 14376 Pattern_1 <NA> Pattern_3
## AAACAAGTATCTCCCA-1 AAACAAGTATCTCCCA-1 16545 25987 Pattern_1 Pattern_2 Pattern_3
## AAACAATCTACTAGCA-1 AAACAATCTACTAGCA-1 5392 18039 <NA> <NA> Pattern_3
## AAACACCAATAACTGC-1 AAACACCAATAACTGC-1 18606 14703 <NA> <NA> <NA>
## AAACAGAGCGACTCCT-1 AAACAGAGCGACTCCT-1 8032 24950 Pattern_1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 AAACAGCTTTCAGAAG-1 14818 13367 <NA> <NA> <NA>
## Pattern_4 Pattern_5
## AAACAACGAATAGTTC-1 Pattern_4 <NA>
## AAACAAGTATCTCCCA-1 Pattern_4 <NA>
## AAACAATCTACTAGCA-1 Pattern_4 <NA>
## AAACACCAATAACTGC-1 <NA> Pattern_5
## AAACAGAGCGACTCCT-1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 Pattern_4 Pattern_5
## [1] "Pattern_1_Pattern_2"
# Detailed stats for the first pair (if any interactions were found)
if (length(pair_names) > 0) {
pair1 <- interactions(sme_undirected, pair = pair_names[1])
if (length(pair1$interacting_genes) > 0) {
head(pair1$interacting_genes[[1]])
} else {
message("No interacting genes found for this pair.")
}
}## Gene Pattern_1 x Pattern_2 KW.obs.tot KW.obs.groups KW.df
## COMP COMP vsBoth 2552 3 2
## MAP2K2 MAP2K2 vsBoth 2552 3 2
## HLA-A HLA-A vsBoth 2552 3 2
## CST1 CST1 vsBoth 2552 3 2
## ARHGEF40 ARHGEF40 vsBoth 2552 3 2
## SLC12A9 SLC12A9 vsBoth 2552 3 2
## KW.statistic KW.pvalue KW.p.adj Dunn.zP1_Int Dunn.zP2_Int
## COMP 46.91225 6.503296e-11 2.471253e-10 -5.034959 -6.686948
## MAP2K2 34.41648 3.361679e-08 1.094947e-07 -5.136201 -5.259909
## HLA-A 181.00503 4.957447e-40 3.532181e-39 -7.453216 -13.450969
## CST1 48.70199 2.657646e-11 1.044730e-10 -3.785701 -6.974570
## ARHGEF40 13.93738 9.408859e-04 1.986315e-03 -2.813471 -3.621020
## SLC12A9 17.76551 1.387617e-04 3.101732e-04 -2.655819 -4.203384
## Dunn.zP2_P1 Dunn.pval_1_Int Dunn.pval_2_Int Dunn.pval_2_1
## COMP -0.8161070 4.779513e-07 2.278734e-11 4.144389e-01
## MAP2K2 0.7093684 2.803481e-07 1.441270e-07 4.780959e-01
## HLA-A -4.7157840 9.109187e-14 3.038637e-41 2.407818e-06
## CST1 -2.5359309 1.532757e-04 3.068081e-12 1.121489e-02
## ARHGEF40 -0.3419204 4.900978e-03 2.934437e-04 7.324108e-01
## SLC12A9 -1.0981652 7.911614e-03 2.629536e-05 2.721324e-01
## Dunn.pval_1_Int.adj Dunn.pval_2_Int.adj Dunn.pval_2_1.adj
## COMP 1.537756e-06 1.053915e-10 5.444109e-01
## MAP2K2 9.429891e-07 5.078761e-07 6.170773e-01
## HLA-A 5.185229e-13 3.212274e-40 7.322404e-06
## CST1 3.739254e-04 1.513586e-11 1.945082e-02
## ARHGEF40 9.628470e-03 6.715928e-04 9.083531e-01
## SLC12A9 1.439654e-02 6.949488e-05 3.683743e-01
## SpaceMarkersMetric
## COMP 5.115554
## MAP2K2 4.808176
## HLA-A 4.476319
## CST1 3.512456
## ARHGEF40 3.483833
## SLC12A9 3.480979
The directed analysis keeps every pattern on the SME so the directed
overlap heatmap covers all source → target pairs, and restricts
gene-level scoring to Pattern_1 → Pattern_5 (immune → cancer) via
pattern_pairs, which SpaceMarkers() forwards
to calculate_gene_scores_directed().
sme_directed <- SpaceMarkers(
x = sme,
directed = TRUE,
min.gene.expr = 3,
pattern_pairs = rbind(c("Pattern_1", "Pattern_5"))
)## number of iterations= 85
## number of iterations= 195
## number of iterations= 77
## number of iterations= 119
## number of iterations= 64
## number of iterations= 85
## number of iterations= 179
## number of iterations= 31
## number of iterations= 113
## number of iterations= 65
## [1] "directed"
## statistic p.value n1 n2 effect_size gene
## C1orf127 1.026674 3.053435e-01 246 922 0.08665055 C1orf127
## CDC42 5.566032 4.968939e-08 246 922 0.40810687 CDC42
## C1QC 6.081864 2.580757e-09 246 922 0.39911510 C1QC
## C1QB 2.676324 7.733878e-03 246 922 0.18059121 C1QB
## ZNF593 7.414034 8.697873e-13 246 922 0.55808811 ZNF593
## YTHDF2 3.475487 5.741073e-04 246 922 0.27064062 YTHDF2
## cell_interaction p.adj
## C1orf127 Pattern_1_near_Pattern_5 3.588573e-01
## CDC42 Pattern_1_near_Pattern_5 1.452459e-07
## C1QC Pattern_1_near_Pattern_5 8.287501e-09
## C1QB Pattern_1_near_Pattern_5 1.366918e-02
## ZNF593 Pattern_1_near_Pattern_5 3.361212e-12
## YTHDF2 Pattern_1_near_Pattern_5 1.223331e-03
Directed analysis also stores hotspot and influence intermediates:
## barcode y x Pattern_1 Pattern_2 Pattern_3
## AAACAACGAATAGTTC-1 AAACAACGAATAGTTC-1 4662 14376 Pattern_1 <NA> Pattern_3
## AAACAAGTATCTCCCA-1 AAACAAGTATCTCCCA-1 16545 25987 Pattern_1 <NA> Pattern_3
## AAACAATCTACTAGCA-1 AAACAATCTACTAGCA-1 5392 18039 <NA> <NA> Pattern_3
## AAACACCAATAACTGC-1 AAACACCAATAACTGC-1 18606 14703 <NA> <NA> <NA>
## AAACAGAGCGACTCCT-1 AAACAGAGCGACTCCT-1 8032 24950 Pattern_1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 AAACAGCTTTCAGAAG-1 14818 13367 <NA> <NA> Pattern_3
## Pattern_4 Pattern_5
## AAACAACGAATAGTTC-1 Pattern_4 <NA>
## AAACAAGTATCTCCCA-1 <NA> <NA>
## AAACAATCTACTAGCA-1 <NA> <NA>
## AAACACCAATAACTGC-1 <NA> Pattern_5
## AAACAGAGCGACTCCT-1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 <NA> Pattern_5
## barcode y x Pattern_1 Pattern_2 Pattern_3
## AAACAACGAATAGTTC-1 AAACAACGAATAGTTC-1 4662 14376 Pattern_1 <NA> Pattern_3
## AAACAAGTATCTCCCA-1 AAACAAGTATCTCCCA-1 16545 25987 <NA> <NA> Pattern_3
## AAACAATCTACTAGCA-1 AAACAATCTACTAGCA-1 5392 18039 <NA> <NA> Pattern_3
## AAACACCAATAACTGC-1 AAACACCAATAACTGC-1 18606 14703 <NA> <NA> <NA>
## AAACAGAGCGACTCCT-1 AAACAGAGCGACTCCT-1 8032 24950 Pattern_1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 AAACAGCTTTCAGAAG-1 14818 13367 <NA> <NA> <NA>
## Pattern_4 Pattern_5
## AAACAACGAATAGTTC-1 Pattern_4 <NA>
## AAACAAGTATCTCCCA-1 <NA> <NA>
## AAACAATCTACTAGCA-1 <NA> <NA>
## AAACACCAATAACTGC-1 <NA> Pattern_5
## AAACAGAGCGACTCCT-1 <NA> <NA>
## AAACAGCTTTCAGAAG-1 <NA> Pattern_5
## barcode x y Pattern_1 Pattern_2
## AAACAACGAATAGTTC-1 AAACAACGAATAGTTC-1 14376 4662 0.34859369 0.10298807
## AAACAAGTATCTCCCA-1 AAACAAGTATCTCCCA-1 25987 16545 0.16456297 0.31619677
## AAACAATCTACTAGCA-1 AAACAATCTACTAGCA-1 18039 5392 0.06432877 0.11121815
## AAACACCAATAACTGC-1 AAACACCAATAACTGC-1 14703 18606 0.05964437 0.15881943
## AAACAGAGCGACTCCT-1 AAACAGAGCGACTCCT-1 24950 8032 0.34172043 0.26491974
## AAACAGCTTTCAGAAG-1 AAACAGCTTTCAGAAG-1 13367 14818 0.08296669 0.06691587
## Pattern_3 Pattern_4 Pattern_5
## AAACAACGAATAGTTC-1 0.17965177 0.62482126 0.022594385
## AAACAAGTATCTCCCA-1 0.23027045 0.35776379 0.039713306
## AAACAATCTACTAGCA-1 0.24245830 0.40456609 0.059994366
## AAACACCAATAACTGC-1 0.06351464 0.07489106 0.566755550
## AAACAGAGCGACTCCT-1 0.02251637 0.17334740 0.009465141
## AAACAGCTTTCAGAAG-1 0.07077331 0.27682185 0.522630852
Since the SME holds all results, every visualization function accepts the object directly and reads from the stored slots.
plot_spatial() automatically extracts coordinates from
the object, finds the tissue image stored in the object’s
visiumDir, and falls back to a blank background if none is
available. The source argument selects where
feature_col is resolved: colData (patterns,
continuous features), assay (gene expression),
hotspots (hotspot labels), influence_map
(per-pattern influence values), or interaction (a three-way
categorical overlay).
Pass source = "hotspots" so plot_spatial()
reads hotspot labels from metadata(sme)$hotspots$undirected
rather than the continuous pattern column in colData (which
shares the same name). colors maps the “hot” level to a
single fill.
plot_spatial(sme_undirected, feature_col = "Pattern_1",
source = "hotspots", hotspot_type = "undirected",
colors = "steelblue", title = "Pattern_1 hotspots")…and Pattern_2, from the same find_all_hotspots()
pass:
plot_spatial(sme_undirected, feature_col = "Pattern_2",
source = "hotspots", hotspot_type = "undirected",
colors = "firebrick", title = "Pattern_2 hotspots")With hotspots in hand we can already answer the broader spatial
question — which pattern pairs actually share territory.
plot_overlap_scores() reads
overlap_scores(sme) directly and auto-detects the
undirected shape (symmetric pattern1 × pattern2 heatmap)
across every pair:
For a chosen pair, a three-way overlay classifies each hotspot spot as Pattern_1 only, interacting (hot in both), or Pattern_2 only. This is the spatial region the per-gene tests below operate on.
plot_spatial(sme_undirected,
source = "interaction", hotspot_type = "undirected",
interaction_patterns = c("Pattern_1", "Pattern_2"),
colors = c(Pattern_1 = "steelblue",
interacting = "purple",
Pattern_2 = "firebrick"))The column names of undirected_scores(sme) are the
pattern-pair keys — use the first pair key robustly instead of
hard-coding:
pair_col <- setdiff(colnames(undirected_scores(sme_undirected)), "Gene")[1]
plot_im_scores(sme_undirected, interaction = pair_col, nGenes = 10)plot_spatial() accepts gene names directly — it falls
through to the expression assay via source = "auto" (the
default):
ims <- undirected_scores(sme_undirected)
top_genes <- head(ims[order(ims[[pair_col]], decreasing = TRUE), "Gene"], 5)
for (g in top_genes) {
print(plot_spatial(sme_undirected, feature_col = g, title = g))
}The directed workflow splits hotspots into two flavors: pattern-value
GMM hotspots (hotspot_type = "pattern") and influence-map
GMM hotspots (hotspot_type = "influence"). Here we focus on
the pattern-level hotspots that define the source/target regions of each
directed pair.
plot_spatial(sme_directed, feature_col = "Pattern_1",
source = "hotspots", hotspot_type = "pattern",
colors = "steelblue", title = "Pattern_1 GMM hotspots")…and Pattern_5, the target region the immune pattern radiates into:
plot_spatial(sme_directed, feature_col = "Pattern_5",
source = "hotspots", hotspot_type = "pattern",
colors = "firebrick", title = "Pattern_5 GMM hotspots")plot_overlap_scores() detects the directed shape
(pattern × influence with a relative-abundance fill)
automatically, again across every source → target pair:
Because a directed analysis runs two independent t-tests per pair — one per source → target direction — the per-spot classification is best rendered one direction at a time. The forward view shows Pattern_1’s hotspot split by whether each spot is also in Pattern_5’s influence zone (the test’s interaction group):
labs_fwd <- overlap_map(sme_directed, c("Pattern_1", "Pattern_5"),
direction = "forward")
table(labs_fwd, useNA = "ifany")## labs_fwd
## Pattern_1 Pattern_1 near Pattern_5 Pattern_5 influence
## 922 246 1451
## <NA>
## 2279
plot_spatial(sme_directed,
source = "interaction",
interaction_patterns = c("Pattern_1", "Pattern_5"),
direction = "forward",
colors = c(Pattern_1 = "steelblue",
"Pattern_1 near Pattern_5" = "purple",
"Pattern_5 influence" = "gray60"))And the reverse (Pattern_5 → Pattern_1):
labs_rev <- overlap_map(sme_directed, c("Pattern_1", "Pattern_5"),
direction = "reverse")
table(labs_rev, useNA = "ifany")## labs_rev
## Pattern_5 Pattern_5 near Pattern_1 Pattern_1 influence
## 1586 246 767
## <NA>
## 2299
plot_spatial(sme_directed,
source = "interaction",
interaction_patterns = c("Pattern_1", "Pattern_5"),
direction = "reverse",
colors = c(Pattern_5 = "firebrick",
"Pattern_5 near Pattern_1" = "goldenrod",
"Pattern_1 influence" = "gray60"))If you already have expression and coordinate data in R, you can
construct a SpaceMarkersExperiment directly:
mat <- matrix(rnorm(100), nrow = 10, ncol = 10)
rownames(mat) <- paste0("gene_", 1:10)
colnames(mat) <- paste0("spot_", 1:10)
coords <- matrix(runif(20) * 100, ncol = 2)
colnames(coords) <- c("y", "x")
rownames(coords) <- colnames(mat)
cd <- S4Vectors::DataFrame(
CellType_A = runif(10),
CellType_B = runif(10),
row.names = colnames(mat)
)
my_sme <- SpaceMarkersExperiment(
assays = list(logcounts = mat),
colData = cd,
spatialCoords = coords,
spaceMarkers = list(
params = list(pattern_names = c("CellType_A", "CellType_B"))
)
)
my_sme## class: SpaceMarkersExperiment
## dim: 10 10
## metadata(0):
## assays(1): logcounts
## rownames(10): gene_1 gene_2 ... gene_9 gene_10
## rowData names(0):
## colnames(10): spot_1 spot_2 ... spot_9 spot_10
## colData names(3): CellType_A CellType_B sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : y x
## imgData names(0):
## spacemarkers analysis: none
## patterns: CellType_A, CellType_B
You can also extract features from any SpatialExperiment
object using get_spatial_features():
spe <- SpatialExperiment::SpatialExperiment(
assays = list(logcounts = mat),
colData = cd,
spatialCoords = coords
)
# Extract numeric colData columns as spatial features
features <- get_spatial_features(spe)
head(features)## CellType_A CellType_B
## spot_1 0.57537205 0.28089681
## spot_2 0.02315201 0.93983437
## spot_3 0.28931453 0.94880442
## spot_4 0.93442737 0.02686118
## spot_5 0.90508234 0.50136822
## spot_6 0.44910650 0.75619456
The legacy file-path interface still works. Set
returnSME = FALSE to get the old-style output:
| Feature | Legacy Workflow | SME Workflow |
|---|---|---|
| Input | File paths | SpaceMarkersExperiment or file paths |
| Output | Scattered data.frames/lists | Single SpaceMarkersExperiment |
| Expression | Separate matrix | assay(sme, "logcounts") |
| Coordinates | Separate data.frame | spatialCoords(sme) |
| Patterns | Merged in spPatterns | spatial_patterns(sme) |
| Kernel Parameters | Separate matrix | spatial_params(sme) |
| All Hyperparameters | Scattered arguments | params(sme) |
| Undirected Scores | Separate data.frame | undirected_scores(sme) |
| Directed Scores | Separate data.frame | directed_scores(sme) |
| LR Scores | Separate matrix | lr_scores(sme) |
| Overlap Scores | Separate data.frame | overlap_scores(sme) |
| Hotspots | Separate data.frame | hotspots(sme, type = ...) |
| Interactions | Nested list | interactions(sme, pair = ...) |
| Influence | Separate data.frame | influence_map(sme) |
| Spatial Plotting | Manual coord extraction +
plot_spatial_data_over_image |
plot_spatial(sme, feature_col = ...) |
## R version 4.6.0 (2026-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] grid stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] ggplot2_4.0.3 viridis_0.6.5 viridisLite_0.4.3 data.table_1.18.4
## [5] dplyr_1.2.1 readbitmap_0.1.5 RColorBrewer_1.1-3 cowplot_1.2.0
## [9] rjson_0.2.23 Matrix_1.7-5 SpaceMarkers_2.3.1 BiocStyle_2.41.0
##
## loaded via a namespace (and not attached):
## [1] sys_3.4.3 jsonlite_2.0.0
## [3] shape_1.4.6.1 magrittr_2.0.5
## [5] spatstat.utils_3.2-2 magick_2.9.1
## [7] farver_2.1.2 rmarkdown_2.31
## [9] GlobalOptions_0.1.4 vctrs_0.7.3
## [11] memoise_2.0.1 spatstat.explore_3.8-0
## [13] matrixTests_0.2.3.1 rstatix_0.7.3
## [15] bmp_0.3.1 htmltools_0.5.9
## [17] S4Arrays_1.13.0 curl_7.1.0
## [19] broom_1.0.12 SparseArray_1.13.2
## [21] Formula_1.2-5 sass_0.4.10
## [23] bslib_0.10.0 htmlwidgets_1.6.4
## [25] plyr_1.8.9 httr2_1.2.2
## [27] plotly_4.12.0 cachem_1.1.0
## [29] buildtools_1.0.0 lifecycle_1.0.5
## [31] pkgconfig_2.0.3 R6_2.6.1
## [33] fastmap_1.2.0 MatrixGenerics_1.25.0
## [35] digest_0.6.39 colorspace_2.1-2
## [37] S4Vectors_0.51.1 tensor_1.5.1
## [39] GenomicRanges_1.65.0 RSQLite_2.4.6
## [41] labeling_0.4.3 filelock_1.0.3
## [43] spatstat.sparse_3.1-0 httr_1.4.8
## [45] polyclip_1.10-7 abind_1.4-8
## [47] compiler_4.6.0 bit64_4.8.0
## [49] withr_3.0.2 S7_0.2.2
## [51] backports_1.5.1 tiff_0.1-12
## [53] BiocParallel_1.47.0 carData_3.0-6
## [55] DBI_1.3.0 MASS_7.3-65
## [57] rappdirs_0.3.4 DelayedArray_0.39.1
## [59] tools_4.6.0 otel_0.2.0
## [61] ape_5.8-1 goftest_1.2-3
## [63] glue_1.8.1 nanoparquet_0.5.1
## [65] nlme_3.1-169 reshape2_1.4.5
## [67] generics_0.1.4 hdf5r_1.3.12
## [69] gtable_0.3.6 spatstat.data_3.1-9
## [71] tidyr_1.3.2 car_3.1-5
## [73] XVector_0.53.0 BiocGenerics_0.59.0
## [75] spatstat.geom_3.7-3 pillar_1.11.1
## [77] stringr_1.6.0 circlize_0.4.18
## [79] splines_4.6.0 BiocFileCache_3.3.0
## [81] lattice_0.22-9 survival_3.8-6
## [83] bit_4.6.0 deldir_2.0-4
## [85] tidyselect_1.2.1 SingleCellExperiment_1.35.0
## [87] maketools_1.3.2 knitr_1.51
## [89] gridExtra_2.3 IRanges_2.47.0
## [91] Seqinfo_1.3.0 SummarizedExperiment_1.43.0
## [93] stats4_4.6.0 xfun_0.57
## [95] Biobase_2.73.1 mixtools_2.0.0.1
## [97] matrixStats_1.5.0 stringi_1.8.7
## [99] lazyeval_0.2.3 yaml_2.3.12
## [101] codetools_0.2-20 evaluate_1.0.5
## [103] kernlab_0.9-33 effsize_0.8.1
## [105] tibble_3.3.1 qvalue_2.45.0
## [107] BiocManager_1.30.27 cli_3.6.6
## [109] segmented_2.2-1 jquerylib_0.1.4
## [111] Rcpp_1.1.1-1.1 spatstat.random_3.4-5
## [113] dbplyr_2.5.2 png_0.1-9
## [115] spatstat.univar_3.1-7 parallel_4.6.0
## [117] blob_1.3.0 jpeg_0.1-11
## [119] SpatialExperiment_1.23.0 scales_1.4.0
## [121] purrr_1.2.2 rlang_1.2.0