| Type: | Package |
| Title: | Detection and Spatial Analysis of Tertiary Lymphoid Structures |
| Version: | 0.3.0 |
| Date: | 2026-04-02 |
| Description: | Fast, reproducible detection and quantitative analysis of tertiary lymphoid structures (TLS) in multiplexed tissue imaging. Implements Independent Component Analysis Trace (ICAT) index, local Ripley's K scanning, automated K Nearest Neighbor (KNN)-based TLS detection, and T-cell clusters identification as described in Amiryousefi et al. (2025) <doi:10.1101/2025.09.21.677465>. |
| License: | MIT + file LICENSE |
| URL: | https://github.com/labsyspharm/tlsR |
| Depends: | R (≥ 4.0.0) |
| Imports: | dbscan (≥ 1.1-10), fastICA (≥ 1.2-3), FNN (≥ 1.1.3), spatstat.explore (≥ 3.0-0), spatstat.geom (≥ 3.0-0), ggplot2 (≥ 3.4.0), rlang (≥ 1.0.0), grDevices, graphics, stats, methods |
| Suggests: | knitr, rmarkdown, testthat (≥ 3.0.0) |
| Encoding: | UTF-8 |
| LazyData: | true |
| RoxygenNote: | 7.3.3 |
| VignetteBuilder: | knitr |
| Config/testthat/edition: | 3 |
| NeedsCompilation: | no |
| Packaged: | 2026-04-23 20:08:52 UTC; ali |
| Author: | Ali Amiryousefi |
| Maintainer: | Ali Amiryousefi <ali_amiryousefi@hms.harvard.edu> |
| Repository: | CRAN |
| Date/Publication: | 2026-04-23 20:40:02 UTC |
tlsR: Detection and Spatial Analysis of Tertiary Lymphoid Structures
Description
Fast, reproducible detection and quantitative analysis of tertiary lymphoid structures (TLS) in multiplexed tissue imaging data.
Typical workflow
Load or prepare a named list of data frames (
ldata), one per tissue sample. Each data frame must contain columnsx,y(spatial coordinates in microns), andphenotype(character:"B cell"/"T cell"/ other).Run
detect_TLSto label B+T co-localised regions.(Optional) Run
scan_clusteringto identify windows of significant immune clustering via local Ripley's L.Run
calc_icatto score the internal linearity/organisation of each detected TLS.Run
detect_ticto identify T-cell clusters outside TLS.Use
summarize_TLSto obtain a tidy summary table.Use
plot_TLSto produce publication-ready spatial plots.
Author(s)
Maintainer: Ali Amiryousefi ali_amiryousefi@hms.harvard.edu (ORCID)
Authors:
Jeremiah Wala jeremiah_wala@dfci.harvard.edu (ORCID)
Other contributors:
Peter Sorger (ORCID) [contributor]
References
Amiryousefi et al. (2025) doi:10.1101/2025.09.21.677465
See Also
Useful links:
Calculate ICAT (Independent Component Analysis Trace) Index for TLS
Description
Quantifies the spatial spread and linear organisation of cells within a detected TLS. FastICA is applied to the (x, y) coordinates of TLS cells to estimate independent components; the mixing matrix is used to reconstruct the data, and the ICAT index is defined as the normalised trace-standard-deviation of the reconstructed coordinates.
The index is always non-negative because it measures the average spatial spread per cell rather than the signed trace of the mixing matrix (which can be negative due to ICA sign ambiguity). Higher values indicate a more spatially extended, structured cluster.
Usage
calc_icat(patientID, tlsID, ldata = NULL)
Arguments
patientID |
Character. Sample name in |
tlsID |
Numeric or integer. TLS identifier (value of |
ldata |
Named list of data frames, or |
Details
The ICAT index is computed as follows:
Centre the (x, y) coordinates of TLS cells.
Run
fastICAwith 2 components.Reconstruct
\hat{X} = S A^T + \mu.Let
v_1, v_2be the marginal variances of\hat{X}.-
\text{ICAT} = 100 \times \frac{\sqrt{v_1 + v_2 + 2\sqrt{v_1 v_2}}}{\text{nrow}(X)}
If the requested TLS contains fewer than 3 cells, or FastICA does not
converge, the function returns NA_real_ with an informative message
rather than throwing an error.
Value
A single non-negative numeric value (the ICAT index), or
NA_real_ if computation is not possible (fewer than 3 cells, or
FastICA did not converge).
Examples
data(toy_ldata)
ldata <- detect_TLS("ToySample", k = 30, ldata = toy_ldata)
if (max(ldata[["ToySample"]]$tls_id_knn, na.rm = TRUE) > 0) {
icat <- calc_icat("ToySample", tlsID = 1, ldata = ldata)
print(icat)
}
Detect Tertiary Lymphoid Structures using a KNN-density approach
Description
Identifies TLS candidates by finding regions of high local B-cell density
that also contain a sufficient number of nearby T cells (B+T
co-localisation). Phenotype labels "B cell" and "B cells"
(and their T-cell equivalents) are both accepted.
Usage
detect_TLS(
LSP,
ldata,
k = 30L,
bcell_density_threshold = 10,
min_B_cells = 50L,
min_T_cells_nearby = 10L,
max_distance_T = 50,
expand_distance = 80
)
Arguments
LSP |
Character. Sample name in |
k |
Integer. Number of nearest neighbours used for density estimation
(default |
bcell_density_threshold |
Numeric. Minimum average 1/k-distance (in
microns) for a B cell to be considered locally dense (default |
min_B_cells |
Integer. Minimum B cells per candidate TLS cluster
(default |
min_T_cells_nearby |
Integer. Minimum T cells within
|
max_distance_T |
Numeric. Search radius (microns) for T-cell proximity
check (default |
expand_distance |
Integer. The extended values from the boundary of the deteced B-cells clusters that the T cells are bieng integrated (default |
ldata |
Named list of data frames, or |
Value
The similarly formatted ldata list, with the data frame for LSP
augmented by three new columns:
tls_id_knnInteger.
0= non-TLS cell; positive integer = TLS cluster ID.tls_center_xNumeric. X coordinate of the TLS centre for TLS cells;
NAotherwise.tls_center_yNumeric. Y coordinate of the TLS centre for TLS cells;
NAotherwise.
Examples
# Use a 70% sample of the data to keep CRAN check time under 10s.
# TLS detection requires sufficient cell density; 70% preserves
# the spatial structure needed for reliable detection.
# For production use, run on the full dataset (see \donttest{} below).
data(toy_ldata)
set.seed(42)
idx <- sample(nrow(toy_ldata[["ToySample"]]),
size = floor(0.7 * nrow(toy_ldata[["ToySample"]])))
sub_ldata <- list(ToySample = toy_ldata[["ToySample"]][idx, ])
ldata <- detect_TLS("ToySample", k = 30, ldata = sub_ldata)
table(ldata[["ToySample"]]$tls_id_knn)
plot(ldata[["ToySample"]]$x, ldata[["ToySample"]]$y,
col = ifelse(ldata[["ToySample"]]$tls_id_knn > 0, "red", "gray"),
pch = 19, cex = 0.5, main = "Detected TLS (70% sample)")
# Full dataset with default settings
data(toy_ldata)
ldata <- detect_TLS("ToySample", k = 30, ldata = toy_ldata)
table(ldata[["ToySample"]]$tls_id_knn)
plot(ldata[["ToySample"]]$x, ldata[["ToySample"]]$y,
col = ifelse(ldata[["ToySample"]]$tls_id_knn > 0, "red", "gray"),
pch = 19, cex = 0.5, main = "Detected TLS in toy data")
Detect T-cell Immune Clusters (TIC) on the Tissue
Description
Applies HDBSCAN to T cells that lie outside of previously detected TLS
regions to identify spatially compact T-cell clusters (TIC). Phenotype
labels "T cell" and "T cells" are both accepted.
Usage
detect_tic(sample, min_pts = 10L, min_cluster_size = 10L, ldata = NULL)
Arguments
sample |
Character. Sample name in |
min_pts |
Integer. HDBSCAN |
min_cluster_size |
Integer. Minimum number of T cells for a HDBSCAN
cluster to be retained; smaller clusters are merged back into noise
(label |
ldata |
Named list of data frames, or |
Value
The input ldata list with the sample data frame augmented by
one new column:
tcell_cluster_hdbscanInteger.
0= noise / not a T-cell cluster; positive integer = TIC cluster ID. Non-T-cell rows receiveNA.
Examples
data(toy_ldata)
ldata <- detect_TLS("ToySample", k = 30, ldata = toy_ldata)
ldata <- detect_tic("ToySample", ldata = ldata)
table(ldata[["ToySample"]]$tcell_cluster_hdbscan, useNA = "ifany")
Plot Spatial Map of TLS and T-cell Clusters
Description
Produces a ggplot2 scatter plot of cell positions, coloured by
TLS membership, T-cell cluster membership, and background phenotype.
Background (non-TLS, non-TIC) cells are rendered with a lower alpha to keep them visually recessive, while TIC cells are drawn slightly larger than TLS cells so they stand out without dominating the plot.
Usage
plot_TLS(
sample,
ldata = NULL,
show_tic = TRUE,
point_size = 0.5,
alpha = 0.7,
bg_alpha = 0.25,
tic_size_mult = 1.8,
tls_palette = c("#0072B2", "#009E73", "#CC79A7", "#D55E00", "#56B4E9", "#F0E442"),
tic_colour = "#E69F00",
bg_colour = "grey80"
)
Arguments
sample |
Character. Sample name in |
ldata |
Named list of data frames, or |
show_tic |
Logical. Colour T-cell clusters (if |
point_size |
Numeric. Base point size for TLS cells and background
cells (default |
alpha |
Numeric. Point transparency for TLS and TIC cells
(default |
bg_alpha |
Numeric. Point transparency for background (non-TLS,
non-TIC) cells (default |
tic_size_mult |
Numeric. Multiplier applied to |
tls_palette |
Character vector of colours for TLS IDs. Recycled if there are more TLS than colours. Default uses a colourblind-friendly palette. |
tic_colour |
Character. Colour for T-cell cluster cells
(default |
bg_colour |
Character. Colour for non-TLS, non-TIC cells
(default |
Value
A ggplot object (invisibly). The plot is also printed unless
the return value is assigned.
Examples
data(toy_ldata)
ldata <- detect_TLS("ToySample", k = 30, ldata = toy_ldata)
p <- plot_TLS("ToySample", ldata = ldata)
Scan Tissue for Local Immune Cell Clustering (Ripley's L Heatmap)
Description
Applies a sliding-window Ripley's L analysis across the tissue to produce
a spatial clustering map. For each window a K-integral index is
computed as the mean positive excess of the observed L function over its
theoretical CSR value. When plot = TRUE a base-graphics spatial
map is drawn with LOESS-smoothed L-excess curves and numeric CI labels
overlaid inside each qualifying window, plus a legend identifying point
and curve colours.
Usage
scan_clustering(
ws = 500,
sample,
phenotype = c("T cells", "B cells", "Both"),
plot = TRUE,
creep = 1L,
min_cells = 10L,
min_phen_cells = 5L,
label_cex = 1.1,
ldata = NULL
)
Arguments
ws |
Numeric. Window side length in microns (default |
sample |
Character. Sample name in |
phenotype |
One of |
plot |
Logical. Draw the spatial clustering map?
(default |
creep |
Integer. Grid density factor; |
min_cells |
Integer. Minimum total cell count required in a window
before it is analysed (default |
min_phen_cells |
Integer. Minimum phenotype-specific cell count per
window (default |
label_cex |
Numeric. Base character expansion for the CI numeric
labels drawn inside each window (default |
ldata |
Named list of data frames, or |
Details
The K-integral clustering index for window w is:
\text{CI}_w = \frac{1}{N_+}\sum_{i:\,L_i > L_{\text{theo},i}}
(L_i - L_{\text{theo},i})
where N_+ is the number of spatial lags where the observed L exceeds
the theoretical CSR value.
When plot = TRUE the map shows:
All cells as small light-grey points.
Phenotype cells (T cells green, B cells red).
Navy dashed grid lines marking window boundaries.
A LOESS-smoothed L-excess curve inside each qualifying window.
A bold numeric CI label centred in the window.
A legend identifying all point and curve colours.
When phenotype = "Both" two side-by-side panels are produced -
one for B cells and one for T cells - so the two clustering maps can be
compared directly on the same spatial layout.
Value
A named list with elements B and/or T (depending on
phenotype), each containing the Lest objects for all
qualifying windows of that phenotype. Returned invisibly when
plot = TRUE.
Examples
data(toy_ldata)
L_models <- scan_clustering(
ws = 200,
sample = "ToySample",
phenotype = "B cells",
plot = TRUE,
ldata = toy_ldata
)
cat("B-cell windows analysed:", length(L_models$B), "\n")
# Side-by-side B and T cell panels
L_both <- scan_clustering(
ws = 200,
sample = "ToySample",
phenotype = "Both",
plot = TRUE,
ldata = toy_ldata
)
Summarize Detected TLS Across Samples
Description
Produces a tidy data.frame with one row per sample summarising the
number of detected TLS, their sizes, and (optionally) ICAT scores.
Usage
summarize_TLS(ldata, calc_icat_scores = FALSE)
Arguments
ldata |
Named list of data frames as returned by |
calc_icat_scores |
Logical. Should ICAT scores be computed for each TLS
and appended as a list-column? Default |
Value
A data.frame with columns:
sampleSample name.
n_TLSNumber of TLS detected.
total_cellsTotal cells in the sample.
TLS_cellsNumber of cells assigned to any TLS.
TLS_fractionFraction of all cells that are TLS cells.
mean_TLS_sizeMean cells per TLS (
NAif n_TLS = 0).n_TICNumber of T-cell clusters detected by
detect_tic(NAif not yet run).icat_scoresList-column of ICAT scores per TLS (only when
calc_icat_scores = TRUE).
Examples
data(toy_ldata)
ldata <- detect_TLS("ToySample", k = 30, ldata = toy_ldata)
summarize_TLS(ldata)
Toy Multiplexed Imaging Data
Description
A small synthetic dataset mimicking multiplexed tissue imaging data, used in
package examples and tests. The list contains one sample named
"ToySample".
Usage
toy_ldata
Format
A named list with one element:
ToySampleA
data.framewith the following columns:xNumeric. X coordinate in microns.
yNumeric. Y coordinate in microns.
phenotypeCharacter. Cell phenotype label. Values are
"B cell","T cell", and"Other".
Source
Synthetically generated for package examples.
References
Amiryousefi et al. (2025) doi:10.1101/2025.09.21.677465
Examples
data(toy_ldata)
str(toy_ldata)
table(toy_ldata[["ToySample"]]$phenotype)