library(cellmig)
library(ggplot2)
library(ggforce)
ggplot2::theme_set(new = theme_bw(base_size = 10))High-throughput tracking of cells with time-lapse microscopy enables the analysis of cell migration across many wells treated under different conditions. Such experiments generate substantial technical and biological variability, making technical and biological replicates necessary. This leads to hierarchically structured datasets: cells are nested within technical replicates, which in turn are nested within biological replicates.
Most current statistical analyses ignore the hierarchical structure and do not explicitly quantify uncertainty from technical or biological variability. The Bioconductor package cellmig addresses this gap by implementing Bayesian hierarchical models for cell migration analysis. It quantifies condition-specific changes in migration speed (e.g., drug effects) while modeling nested data structures, producing uncertainty-aware estimates through credible intervals.
Currently, there are no other Bioconductor packages specialized for hierarchical high-throughput cell migration data analysis. cellmig addresses this gap and integrates naturally into the Bioconductor ecosystem.
To install this package, start R (version “4.5”) and enter:
This is how a typical cell migration data looks like \(\rightarrow\) a table.
Each row is a cell with the following features:
FALSE 'data.frame': 7560 obs. of 6 variables:
FALSE $ well : chr "1" "1" "1" "1" ...
FALSE $ plate : chr "1" "1" "1" "1" ...
FALSE $ compound: chr "C1" "C1" "C1" "C1" ...
FALSE $ dose : chr "D1" "D1" "D1" "D1" ...
FALSE $ v : num 21.905 0.535 3.348 5.351 1.194 ...
FALSE $ offset : num 1 1 1 1 1 1 1 1 1 1 ...
In this vignette we will use simulated data from:
vEach dot represents a cell; y-axis is velocity. Facets represent compounds, x-axis represents dose. Plate is indicated by color. Technical replicates (wells) are stacked next to each other and have the same color.
ggplot(data = d)+
facet_wrap(facets = ~paste0("compound=", compound),
scales = "free_y", ncol = 2)+
geom_sina(aes(x = as.factor(dose), col = plate, y = v, group = well),
size = 0.5)+
theme_bw()+
theme(legend.position = "top",
strip.text.x = element_text(margin = margin(0.03,0,0.03,0, "cm")))+
ylab(label = "migration speed")+
xlab(label = '')+
scale_color_grey()+
guides(color = guide_legend(override.aes = list(size = 3)))+
guides(shape = guide_legend(override.aes = list(size = 3)))+
scale_y_log10()+
annotation_logticks(base = 10, sides = "l")VIsualizing mean speed within wells highlights plate-specific batch effects.
dm <- aggregate(v~well+plate+compound+dose, data = d, FUN = mean)
ggplot(data = dm)+
facet_wrap(facets = ~paste0("compound=", compound),
scales = "free_y", ncol = 2)+
geom_sina(aes(x = as.factor(dose), col = plate, y = v, group = well),
size = 1.5, alpha = 0.7)+
theme_bw()+
theme(legend.position = "top",
strip.text.x = element_text(margin = margin(0.03,0,0.03,0, "cm")))+
ylab(label = "migration speed")+
xlab(label = '')+
scale_color_grey()+
guides(color = guide_legend(override.aes = list(size = 3)))+
guides(shape = guide_legend(override.aes = list(size = 3)))+
scale_y_log10()+
annotation_logticks(base = 10, sides = "l")cellmig analysisWe will use this data to infer the overall treatment effects (parameter \(\delta_t\)), relative to a control treatment (the offset) to correct for plate-specific batch effects. At the same time, cellmig will quantify many different features of the data using its model parameters (e.g., variability between technical or biological replicates; or plate-specific treatment effects (\(\gamma_{pt}\))).
We fit the Stan model employed by cellmig
with the control parameters defined in the list control.
There are many other input parameters in control, check the
cellmig function documentation.
To extract the means, medians, and 95% Highest Density Intervals
(HDIs, quantifying parameter value uncertainty) of \(\delta_t\), we have to access the
data.frame delta_t in the output object
posteriors:
FALSE 'data.frame': 35 obs. of 16 variables:
FALSE $ group_id: int 1 2 3 4 5 6 7 8 9 10 ...
FALSE $ mean : num -0.8629 -0.454 -0.2045 0.0638 0.1166 ...
FALSE $ se_mean : num 0.00194 0.00165 0.002 0.00176 0.00171 ...
FALSE $ sd : num 0.071 0.0715 0.0729 0.0741 0.0686 ...
FALSE $ X2.5. : num -1.006 -0.5902 -0.3471 -0.0724 -0.0199 ...
FALSE $ X25. : num -0.90792 -0.50277 -0.25422 0.00997 0.0702 ...
FALSE $ X50. : num -0.8638 -0.4562 -0.2023 0.0612 0.1168 ...
FALSE $ X75. : num -0.813 -0.406 -0.158 0.116 0.16 ...
FALSE $ X97.5. : num -0.7258 -0.3162 -0.0627 0.2076 0.2511 ...
FALSE $ n_eff : num 1336 1870 1334 1777 1613 ...
FALSE $ Rhat : num 1.002 1.002 1.003 0.999 1.001 ...
FALSE $ group : chr "C2|D1" "C2|D2" "C2|D3" "C2|D4" ...
FALSE $ compound: chr "C2" "C2" "C2" "C2" ...
FALSE $ dose : chr "D1" "D2" "D3" "D4" ...
FALSE $ plate_id: num 1 1 1 1 1 1 1 1 1 1 ...
FALSE $ plate : chr "1" "1" "1" "1" ...
It is better to visualize the mean \(\delta_t\)s and their 95% HDIs
As compound t=1 was selected as control (by setting offset=1), the treatment effects of this compounds are not shown.
ggplot(data = o$posteriors$delta_t)+
geom_line(aes(x = dose, y = mean, col = compound, group = compound))+
geom_point(aes(x = dose, y = mean, col = compound))+
geom_errorbar(aes(x = dose, y = mean, ymin = X2.5., ymax = X97.5.,
col = compound), width = 0.1)+
ylab(label = expression("Overall treatment effect ("*delta*")"))+
theme(legend.position = "top")profilesFor “rectangular datasets”, i.e. datasets with multiple compounds and overlapping doses, we can study the treatment dose-response profiles by hierarchical clustering based on the complete posteriors of \(\delta_t\), account for uncertainty in this parameter.
Panel A: dendrogram constructed by hierarchical clustering with average linkage, based on euclidean distances between vectors of \(\delta_t\) (shown in panel B) of each compound (leaf) across doses. Branch support values show branch robustness (label = 1000 implies this branch was encountered in each of the 1000 dendrograms constructed from the posterior of \(\delta_t\)). Plate-specific treatment dose-responses based on parameters \(\gamma_pt\).
Pairwise dot-plot comparison \(\rightarrow\) x minus y axis
(Left panel) Differences in overall treatment effects. Log fold change (LFC; described by parameter \(\rho_{ij}\)) between overall treatments effects (\(\delta_t\)) of row (\(i\)) vs. column (\(j\)) treatment groups. Tile colors and labels represent \(\rho_{ij}\). (Right panel) Probability of differential treatment effect described by parameter \(\pi_{ij}\). Tile colors and labels represent \(\pi_{ij}\).
FALSE NULL
from_groups vs. to_group).FALSE 'data.frame': 35 obs. of 4 variables:
FALSE $ group_id: int 1 2 3 4 5 6 7 8 9 10 ...
FALSE $ group : chr "C2|D1" "C2|D2" "C2|D3" "C2|D4" ...
FALSE $ compound: chr "C2" "C2" "C2" "C2" ...
FALSE $ dose : chr "D1" "D2" "D3" "D4" ...
u <- get_violins(x = o,
from_groups = get_groups(x = o)$group,
to_group = "C2|D1",
exponentiate = FALSE)
u$plotTo assess model validity, we performed posterior predictive checks, which showed that the simulated data (pink violin) were consistent with the observed data (black violins). Each dot is a cell.
g_alpha_p <- ggplot(data = o$posteriors$alpha_p)+
geom_errorbarh(aes(y = plate_id, x = mean, xmin = X2.5., xmax = X97.5.),
height = 0.1)+
geom_point(aes(y = plate_id, x = mean))
g_sigma <- ggplot()+
geom_errorbarh(data = o$posteriors$sigma_bio,
aes(y = "sigma_bio",
x = mean, xmin = X2.5., xmax = X97.5.),
height = 0.1)+
geom_errorbarh(data = o$posteriors$sigma_tech,
aes(y = "sigma_tech",
x = mean, xmin = X2.5., xmax = X97.5.),
height = 0.1)+
geom_point(data = o$posteriors$sigma_bio,
aes(y = "sigma_bio", x = mean))+
geom_point(data = o$posteriors$sigma_tech,
aes(y = "sigma_tech", x = mean))+
ylab(label = '')
g_alpha_p|g_sigmaFALSE R version 4.5.2 (2025-10-31)
FALSE Platform: x86_64-pc-linux-gnu
FALSE Running under: Ubuntu 24.04.3 LTS
FALSE
FALSE Matrix products: default
FALSE BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
FALSE LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
FALSE
FALSE locale:
FALSE [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
FALSE [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
FALSE [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
FALSE [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
FALSE [9] LC_ADDRESS=C LC_TELEPHONE=C
FALSE [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
FALSE
FALSE time zone: Etc/UTC
FALSE tzcode source: system (glibc)
FALSE
FALSE attached base packages:
FALSE [1] stats graphics grDevices utils datasets methods base
FALSE
FALSE other attached packages:
FALSE [1] ggforce_0.5.0 ggplot2_4.0.1 cellmig_1.1.6 BiocStyle_2.39.0
FALSE
FALSE loaded via a namespace (and not attached):
FALSE [1] ggiraph_0.9.3 tidyselect_1.2.1 dplyr_1.1.4
FALSE [4] farver_2.1.2 loo_2.9.0 S7_0.2.1
FALSE [7] fastmap_1.2.0 lazyeval_0.2.2 tweenr_2.0.3
FALSE [10] fontquiver_0.2.1 digest_0.6.39 lifecycle_1.0.5
FALSE [13] StanHeaders_2.32.10 tidytree_0.4.7 magrittr_2.0.4
FALSE [16] compiler_4.5.2 rlang_1.1.7 sass_0.4.10
FALSE [19] tools_4.5.2 yaml_2.3.12 knitr_1.51
FALSE [22] labeling_0.4.3 htmlwidgets_1.6.4 pkgbuild_1.4.8
FALSE [25] curl_7.0.0 plyr_1.8.9 RColorBrewer_1.1-3
FALSE [28] aplot_0.2.9 withr_3.0.2 purrr_1.2.1
FALSE [31] sys_3.4.3 grid_4.5.2 polyclip_1.10-7
FALSE [34] stats4_4.5.2 gdtools_0.4.4 inline_0.3.21
FALSE [37] scales_1.4.0 MASS_7.3-65 cli_3.6.5
FALSE [40] rmarkdown_2.30 treeio_1.35.0 generics_0.1.4
FALSE [43] RcppParallel_5.1.11-1 ggtree_4.1.1 reshape2_1.4.5
FALSE [46] ape_5.8-1 cachem_1.1.0 rstan_2.32.7
FALSE [49] stringr_1.6.0 parallel_4.5.2 ggplotify_0.1.3
FALSE [52] BiocManager_1.30.27 matrixStats_1.5.0 vctrs_0.7.1
FALSE [55] yulab.utils_0.2.3 V8_8.0.1 jsonlite_2.0.0
FALSE [58] fontBitstreamVera_0.1.1 gridGraphics_0.5-1 patchwork_1.3.2
FALSE [61] systemfonts_1.3.1 maketools_1.3.2 tidyr_1.3.2
FALSE [64] jquerylib_0.1.4 glue_1.8.0 codetools_0.2-20
FALSE [67] stringi_1.8.7 gtable_0.3.6 QuickJSR_1.9.0
FALSE [70] tibble_3.3.1 pillar_1.11.1 rappdirs_0.3.4
FALSE [73] htmltools_0.5.9 R6_2.6.1 evaluate_1.0.5
FALSE [76] lattice_0.22-7 ggfun_0.2.0 fontLiberation_0.1.0
FALSE [79] bslib_0.10.0 rstantools_2.6.0 Rcpp_1.1.1
FALSE [82] gridExtra_2.3 nlme_3.1-168 xfun_0.56
FALSE [85] fs_1.6.6 buildtools_1.0.0 pkgconfig_2.0.3