StatescopeR is the R version of Statescope, a computational framework designed to discover cell states from cell type-specific gene expression profiles inferred from bulk RNA profiles. StatescopeR works in three steps: 1) Multi-omic single cell RNAseq reference-based deconvolution 2) Refinement of inferred cell type-specific gene expression profiles 3) Cell state discovery In addition to bulk RNA and a single cell RNAseq (scRNAseq) reference, StatescopeR can also integrate prior expectations of cell fractions in the bulk RNA. For example tumor fractions estimated by DNA copy number analysis or immune cell fractions by immunohistochemistry. Below is a manual for running StatescopeR with some pancreas scRNAseq data used for both the reference and pseudobulk which is deconvolved.
StatescopeR is available through Bioconductor, install it using the following commands in R:
if (!requireNamespace("BiocManager", quietly = TRUE)) {
install.packages("BiocManager")
}
BiocManager::install("StatescopeR")
## Check that you have a valid Bioconductor installation
BiocManager::valid()
After the package is installed, we will be loading it and the scRNAseq package, from which we will use the SegerstolpePancreas data set for our example. We will use it for two purposes: 1) Creating a scRNAseq reference, as you would do for your own application. 2) Creating pseudobulk to be deconvolved, which we can conveniently use for evaluation as well. In your own application this would likely be bulk mRNA sequencing.
## Load StatescopeR & scRNAseq (for example data)
suppressMessages(library(StatescopeR))
suppressMessages(library(scRNAseq))
suppressMessages(library(scuttle))
## Load SegerstolpePancreas data set
scRNAseq <- SegerstolpePancreasData()
## remove duplicate genes
scRNAseq <- scRNAseq[!duplicated(rownames(scRNAseq)), ]
## Subset to 1 healthy and type 2 diabetes samples
scRNAseq <- scRNAseq[, scRNAseq$individual %in% c("H2",
"T2D2", "T2D3")]
## remove cells with no cell type label
scRNAseq <- scRNAseq[, !is.na(scRNAseq$`cell type`)]
## remove very rare cell types (<100 cells in total data set)
celltypes_to_remove <-
names(table(scRNAseq$`cell type`)[(table(scRNAseq$`cell type`) < 100)])
scRNAseq <- scRNAseq[, !scRNAseq$`cell type` %in% celltypes_to_remove]
StatescopeR takes the following inputs: 1) scRNAseq reference/signature 2) bulk mRNA to be deconvolved 3) selected genes for deconvolution, by default this uses AutoGeneS (optional) Prior expectation of fractions (optional) Hyperparameters (i.e. cores for parallelization, nrepfinal for maximum number of optimization iteration etc.)
## Normalize (cp10k) and logtransform scRNAseq
cpm(scRNAseq) <- calculateCPM(scRNAseq)
logcounts(scRNAseq) <- log1p(cpm(scRNAseq) / 100)
## Create scRNAseq reference/signature with 5 hvg for quick example
signature <- create_signature(scRNAseq,
hvg_genes = TRUE, n_hvg_genes = 5L,
labels = scRNAseq$`cell type`
)
## select subset of genes for deconvolution (3/5 hvg to make it quick)
selected_genes <- select_genes(scRNAseq, 3L, 5L,
labels = scRNAseq$`cell type`
)
## Create pseudobulk and normalize to cp10k (logging is done within Statescope)
pseudobulk <- aggregateAcrossCells(scRNAseq, ids = scRNAseq$individual)
normcounts(pseudobulk) <- calculateCPM(pseudobulk) / 100
pseudobulk <- as(pseudobulk, "SummarizedExperiment")
rownames(pseudobulk) <- rownames(scRNAseq)
## (optional) Create prior expectation
prior <- gather_true_fractions(scRNAseq, ids = scRNAseq$individual,
label_col = "cell type") # Use True sc fractions for this
prior[rownames(prior) != "ductal cell", ] <- NA # Keep only ductal cells
prior <- t(prior) # Transpose it to nSample x nCelltype
Now that all inputs have been prepared you can deconvolve the (pseudo)bulk mRNA and discover transcriptional states. In short the modules work as follows: BLADE_deconvolution employs Bayesian variational inference, modelling the bulk gene expression level of each gene j in sample i by: \[\begin{equation} y_{ij} = log(\sum_{t} f_{i}^t exp(x_{ij}^t)) + \epsilon_{ij} \end{equation}\] with hidden variables \(f_{i}^t\) and \(x_{ij}^t\) denoting the fraction of cell type t for sample i and the expression of gene j for sample i in cell type j respectively. Refinement further optimizes cell type-specific gene expression profiles, \(x_{ij}^t\), by fixing estimated cell fractions and putting more emphasis on capturing inter-sample heterogeneity over staying close to the single cell prior knowledge. StateDiscovery uses convex-NMF to cluster inferred cell type-specific gene expression profiles in an unsupervised manner. These modules are wrappers around the original Python version Statescope.
## Run Deconvolution module
Statescope <- BLADE_deconvolution(signature, pseudobulk, selected_genes, prior,
cores = 1L, Nrep = 1L)
## Run Refinement module
Statescope <- Refinement(Statescope, signature, pseudobulk, cores = 2L)
## Run State Discovery module (pick k=2 for quick example)
Statescope <- StateDiscovery(Statescope, k = 2L, Ncores = 2L)
As we used pseudobulk from scRNAseq in this manual we are able to check if the output of StatescopeR is correct.
## Count fractions of cells per sample in the pseudobulk
true_fractions <- gather_true_fractions(scRNAseq, ids = scRNAseq$individual,
label_col = "cell type")
## Plot correlation and RMSE with true fractions per celltype
fraction_eval(Statescope, true_fractions)
After you have run Statescope you can do all kinds of analyses, here we show some visualizations to give you an idea, but feel free to try whatever comes to your mind.
## Show predicted fractions
metadata(Statescope)$fractions
#> DataFrame with 2 rows and 3 columns
#> H2 T2D2 T2D3
#> <numeric> <numeric> <numeric>
#> alpha cell 0.819047 0.564368 0.403671
#> ductal cell 0.180953 0.435632 0.596329
## Show predicted ductal cell type specific gene expression profiles
assay(metadata(Statescope)$ct_specific_gep$`ductal cell`, "weighted_gep")
#> DataFrame with 5 rows and 3 columns
#> H2 T2D2 T2D3
#> <numeric> <numeric> <numeric>
#> GCG -0.333965 -0.327920 -0.226353
#> TTR -0.382767 -0.458388 -0.461426
#> SPP1 0.230114 0.190857 0.229792
#> CHGB -0.369818 -0.370663 -0.317205
#> SCG5 -0.187184 -0.212151 -0.213579
## Show ductal cell state scores per sample
metadata(Statescope)$statescores$`ductal cell`
#> DataFrame with 3 rows and 2 columns
#> V1 V2
#> <numeric> <numeric>
#> H2 0.992997 0.00700325
#> T2D2 0.815823 0.18417679
#> T2D3 0.042516 0.95748403
## Show ductal cell state loadings
metadata(Statescope)$stateloadings$`ductal cell`
#> DataFrame with 5 rows and 2 columns
#> V1 V2
#> <numeric> <numeric>
#> GCG -0.300287 -0.207507
#> TTR -0.376014 -0.422140
#> SPP1 0.193745 0.210188
#> CHGB -0.335546 -0.290462
#> SCG5 -0.179287 -0.195407
## Plot heatmap of fractions
fraction_heatmap(Statescope)
## Plot barplot of top loadings
barplot_stateloadings(Statescope)
If you have an expected fraction for a group of cell types instead of for one specific cell (i.e. Lymphoid cells instead of T/B cells), StatescopeR also allows you to give a group prior. Here we show an example where we know the expected fraction of the pancreatic non-duct cells (acinar & alpha cells) together, but not their individual contributions.
## define cell types in one group
grouped_cts <- c("alpha cell", "acinar cell")
## initialize group assignment matrix (nCelltype x nGroup) with 0's
celltype_names <- colnames(signature$mu)
group_names <- c("Group", setdiff(celltype_names, grouped_cts))
nCelltype <- length(celltype_names)
nGroup <- length(group_names)
group <- matrix(0,
nrow = nGroup, ncol = nCelltype,
dimnames = list(group_names, celltype_names)
)
## initialize prior fraction matrix (nSamples x nCelltypes) with NA's
nSamples <- ncol(pseudobulk)
sample_names <- colnames(pseudobulk)
prior <- matrix(
nrow = nSamples, ncol = nGroup,
dimnames = list(sample_names, group_names)
)
## Assign cell types to groups
for (ct in celltype_names) {
if (ct %in% grouped_cts) {
## assign cell type to group
group["Group", ct] <- 1
} else {
## Assign cell type to its own group
group[ct, ct] <- 1
}
}
## Add grouped prior fractions to prior
prior[names(true_fractions), "Group"] <- colSums(as.matrix(
true_fractions[grouped_cts, ]))
## init group prior
group_prior <- list("Group" = group, "Expectation" = prior)
## Now you can just run Statescope as with a normal prior. Note we do not give
## a ductal cell prior this time
Statescope <- BLADE_deconvolution(signature, pseudobulk, selected_genes,
group_prior,
cores = 1L,
Nrep = 1L
)
StatescopeRWe hope that StatescopeR will be useful for your research. Please use the following information to cite the package and the overall approach. Thank you!
## Citation info
citation("StatescopeR")
#> To cite package 'StatescopeR' in publications use:
#>
#> Janssen J, Radonic T, Steketee M, van Asten S, Eijk P, van Maldegem
#> F, Noske D, Bahce I, Vallejo JG, van de Wiel M, Ylstra B, Kim Y
#> (2024). "Hidden RNA profiles of cells in the tumor microenvironment
#> accurately revealed by malignant cell fraction-informed
#> deconvolution." _Research Square_. doi:10.21203/rs.3.rs-4252952/v1
#> <https://doi.org/10.21203/rs.3.rs-4252952/v1>.
#> <https://www.researchgate.net/publication/380850052_Hidden_RNA_profiles_of_cells_in_the_tumor_microenvironment_accurately_revealed_by_malignant_cell_fraction-informed_deconvolution>.
#>
#> A BibTeX entry for LaTeX users is
#>
#> @Article{,
#> title = {Hidden RNA profiles of cells in the tumor microenvironment accurately revealed by malignant cell fraction-informed deconvolution},
#> author = {Jurriaan Janssen and Teodora Radonic and Mischa Steketee and Saskia {van Asten} and Paul Eijk and Febe {van Maldegem} and David Noske and Idris Bahce and Juan Garcia Vallejo and Mark {van de Wiel} and Bauke Ylstra and Yongsoo Kim},
#> year = {2024},
#> journal = {Research Square},
#> doi = {10.21203/rs.3.rs-4252952/v1},
#> url = {https://www.researchgate.net/publication/380850052_Hidden_RNA_profiles_of_cells_in_the_tumor_microenvironment_accurately_revealed_by_malignant_cell_fraction-informed_deconvolution},
#> }
The StatescopeR package (Janssen, Radonic, Steketee, van Asten, Eijk, van Maldegem, Noske, Bahce, Vallejo, van de Wiel, Ylstra, and Kim, 2024) was made possible thanks to:
This package was developed using biocthis.
R session information.
#> ─ Session info ───────────────────────────────────────────────────────────────
#> setting value
#> version R Under development (unstable) (2025-10-20 r88955)
#> os Ubuntu 24.04.3 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate C
#> ctype en_US.UTF-8
#> tz America/New_York
#> date 2025-12-02
#> pandoc 2.7.3 @ /usr/bin/ (via rmarkdown)
#> quarto 1.8.25 @ /usr/local/bin/quarto
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [2] CRAN (R 4.6.0)
#> alabaster.base 1.11.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.matrix 1.11.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.ranges 1.11.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.sce 1.11.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.schemas 1.11.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.se 1.11.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationDbi 1.73.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationFilter 1.35.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationHub 4.1.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> backports 1.5.0 2024-05-23 [2] CRAN (R 4.6.0)
#> basilisk 1.23.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> beachmat 2.27.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> bibtex 0.5.1 2023-01-26 [2] CRAN (R 4.6.0)
#> Biobase * 2.71.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocFileCache 3.1.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocGenerics * 0.57.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocIO 1.21.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocManager 1.30.27 2025-11-14 [2] CRAN (R 4.6.0)
#> BiocNeighbors 2.5.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocParallel 1.45.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocSingular 1.27.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocStyle * 2.39.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocVersion 3.23.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> Biostrings 2.79.2 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> bit 4.6.0 2025-03-06 [2] CRAN (R 4.6.0)
#> bit64 4.6.0-1 2025-01-16 [2] CRAN (R 4.6.0)
#> bitops 1.0-9 2024-10-03 [2] CRAN (R 4.6.0)
#> blob 1.2.4 2023-03-17 [2] CRAN (R 4.6.0)
#> bluster 1.21.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> bookdown 0.45 2025-10-03 [2] CRAN (R 4.6.0)
#> bslib 0.9.0 2025-01-30 [2] CRAN (R 4.6.0)
#> cachem 1.1.0 2024-05-16 [2] CRAN (R 4.6.0)
#> Cairo 1.7-0 2025-10-29 [2] CRAN (R 4.6.0)
#> cigarillo 1.1.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> circlize 0.4.16 2024-02-20 [2] CRAN (R 4.6.0)
#> cli 3.6.5 2025-04-23 [2] CRAN (R 4.6.0)
#> clue 0.3-66 2024-11-13 [2] CRAN (R 4.6.0)
#> cluster 2.1.8.1 2025-03-12 [3] CRAN (R 4.6.0)
#> codetools 0.2-20 2024-03-31 [3] CRAN (R 4.6.0)
#> colorspace 2.1-2 2025-09-22 [2] CRAN (R 4.6.0)
#> ComplexHeatmap 2.27.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> cowplot 1.2.0 2025-07-07 [2] CRAN (R 4.6.0)
#> crayon 1.5.3 2024-06-20 [2] CRAN (R 4.6.0)
#> curl 7.0.0 2025-08-19 [2] CRAN (R 4.6.0)
#> DBI 1.2.3 2024-06-02 [2] CRAN (R 4.6.0)
#> dbplyr 2.5.1 2025-09-10 [2] CRAN (R 4.6.0)
#> DelayedArray 0.37.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> dichromat 2.0-0.1 2022-05-02 [2] CRAN (R 4.6.0)
#> digest 0.6.39 2025-11-19 [2] CRAN (R 4.6.0)
#> dir.expiry 1.19.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> doParallel 1.0.17 2022-02-07 [2] CRAN (R 4.6.0)
#> dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.6.0)
#> dqrng 0.4.1 2024-05-28 [2] CRAN (R 4.6.0)
#> edgeR 4.9.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> ensembldb 2.35.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> evaluate 1.0.5 2025-08-27 [2] CRAN (R 4.6.0)
#> ExperimentHub 3.1.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> farver 2.1.2 2024-05-13 [2] CRAN (R 4.6.0)
#> fastmap 1.2.0 2024-05-15 [2] CRAN (R 4.6.0)
#> filelock 1.0.3 2023-12-11 [2] CRAN (R 4.6.0)
#> foreach 1.5.2 2022-02-02 [2] CRAN (R 4.6.0)
#> generics * 0.1.4 2025-05-09 [2] CRAN (R 4.6.0)
#> GenomeInfoDb 1.47.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicAlignments 1.47.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicFeatures 1.63.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicRanges * 1.63.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> GetoptLong 1.1.0 2025-11-28 [2] CRAN (R 4.6.0)
#> ggplot2 4.0.1 2025-11-14 [2] CRAN (R 4.6.0)
#> GlobalOptions 0.1.3 2025-11-28 [2] CRAN (R 4.6.0)
#> glue 1.8.0 2024-09-30 [2] CRAN (R 4.6.0)
#> gtable 0.3.6 2024-10-25 [2] CRAN (R 4.6.0)
#> gypsum 1.7.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> h5mread 1.3.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> HDF5Array 1.39.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.6.0)
#> httr 1.4.7 2023-08-15 [2] CRAN (R 4.6.0)
#> httr2 1.2.1 2025-07-22 [2] CRAN (R 4.6.0)
#> igraph 2.2.1 2025-10-27 [2] CRAN (R 4.6.0)
#> IRanges * 2.45.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> irlba 2.3.5.1 2022-10-03 [2] CRAN (R 4.6.0)
#> iterators 1.0.14 2022-02-05 [2] CRAN (R 4.6.0)
#> jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.6.0)
#> jsonlite 2.0.0 2025-03-27 [2] CRAN (R 4.6.0)
#> KEGGREST 1.51.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> knitr 1.50 2025-03-16 [2] CRAN (R 4.6.0)
#> labeling 0.4.3 2023-08-29 [2] CRAN (R 4.6.0)
#> lattice 0.22-7 2025-04-02 [3] CRAN (R 4.6.0)
#> lazyeval 0.2.2 2019-03-15 [2] CRAN (R 4.6.0)
#> lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.6.0)
#> limma 3.67.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> locfit 1.5-9.12 2025-03-05 [2] CRAN (R 4.6.0)
#> lubridate 1.9.4 2024-12-08 [2] CRAN (R 4.6.0)
#> magick 2.9.0 2025-09-08 [2] CRAN (R 4.6.0)
#> magrittr 2.0.4 2025-09-12 [2] CRAN (R 4.6.0)
#> Matrix 1.7-4 2025-08-28 [3] CRAN (R 4.6.0)
#> MatrixGenerics * 1.23.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> matrixStats * 1.5.0 2025-01-07 [2] CRAN (R 4.6.0)
#> memoise 2.0.1 2021-11-26 [2] CRAN (R 4.6.0)
#> metapod 1.19.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> pillar 1.11.1 2025-09-17 [2] CRAN (R 4.6.0)
#> pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.6.0)
#> plyr 1.8.9 2023-10-02 [2] CRAN (R 4.6.0)
#> png 0.1-8 2022-11-29 [2] CRAN (R 4.6.0)
#> ProtGenerics 1.43.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> R6 2.6.1 2025-02-15 [2] CRAN (R 4.6.0)
#> rappdirs 0.3.3 2021-01-31 [2] CRAN (R 4.6.0)
#> RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.6.0)
#> Rcpp 1.1.0 2025-07-02 [2] CRAN (R 4.6.0)
#> RCurl 1.98-1.17 2025-03-22 [2] CRAN (R 4.6.0)
#> RefManageR * 1.4.0 2022-09-30 [2] CRAN (R 4.6.0)
#> restfulr 0.0.16 2025-06-27 [2] CRAN (R 4.6.0)
#> reticulate 1.44.1 2025-11-14 [2] CRAN (R 4.6.0)
#> rhdf5 2.55.11 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> rhdf5filters 1.23.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> Rhdf5lib 1.33.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> rjson 0.2.23 2024-09-16 [2] CRAN (R 4.6.0)
#> rlang 1.1.6 2025-04-11 [2] CRAN (R 4.6.0)
#> rmarkdown 2.30 2025-09-28 [2] CRAN (R 4.6.0)
#> Rsamtools 2.27.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> RSQLite 2.4.5 2025-11-30 [2] CRAN (R 4.6.0)
#> rsvd 1.0.5 2021-04-16 [2] CRAN (R 4.6.0)
#> rtracklayer 1.71.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> S4Arrays 1.11.1 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> S4Vectors * 0.49.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> S7 0.2.1 2025-11-14 [2] CRAN (R 4.6.0)
#> sass 0.4.10 2025-04-11 [2] CRAN (R 4.6.0)
#> ScaledMatrix 1.19.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> scales 1.4.0 2025-04-24 [2] CRAN (R 4.6.0)
#> scran 1.39.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> scRNAseq * 2.25.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> scuttle * 1.21.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> Seqinfo * 1.1.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> sessioninfo * 1.2.3 2025-02-05 [2] CRAN (R 4.6.0)
#> shape 1.4.6.1 2024-02-23 [2] CRAN (R 4.6.0)
#> SingleCellExperiment * 1.33.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> SparseArray 1.11.7 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> StatescopeR * 0.99.31 2025-12-02 [1] Bioconductor 3.23 (R 4.6.0)
#> statmod 1.5.1 2025-10-09 [2] CRAN (R 4.6.0)
#> stringi 1.8.7 2025-03-27 [2] CRAN (R 4.6.0)
#> stringr 1.6.0 2025-11-04 [2] CRAN (R 4.6.0)
#> SummarizedExperiment * 1.41.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> tibble 3.3.0 2025-06-08 [2] CRAN (R 4.6.0)
#> tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.6.0)
#> timechange 0.3.0 2024-01-18 [2] CRAN (R 4.6.0)
#> UCSC.utils 1.7.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.6.0)
#> withr 3.0.2 2024-10-28 [2] CRAN (R 4.6.0)
#> xfun 0.54 2025-10-30 [2] CRAN (R 4.6.0)
#> XML 3.99-0.20 2025-11-08 [2] CRAN (R 4.6.0)
#> xml2 1.5.1 2025-12-01 [2] CRAN (R 4.6.0)
#> XVector 0.51.0 2025-12-02 [2] Bioconductor 3.23 (R 4.6.0)
#> yaml 2.3.11 2025-11-28 [2] CRAN (R 4.6.0)
#>
#> [1] /tmp/Rtmp76vXrG/Rinst52a164fd2916
#> [2] /home/biocbuild/bbs-3.23-bioc/R/site-library
#> [3] /home/biocbuild/bbs-3.23-bioc/R/library
#> * ── Packages attached to the search path.
#>
#> ─ Python configuration ───────────────────────────────────────────────────────
#> python: /var/cache/basilisk/1.23.0/StatescopeR/0.99.31/autogenes/bin/python
#> libpython: /home/pkgbuild/.pyenv/versions/3.6.15/lib/libpython3.6m.so
#> pythonhome: /var/cache/basilisk/1.23.0/StatescopeR/0.99.31/autogenes:/var/cache/basilisk/1.23.0/StatescopeR/0.99.31/autogenes
#> version: 3.6.15 (default, Nov 26 2025, 12:41:41) [GCC 13.3.0]
#> numpy: /var/cache/basilisk/1.23.0/StatescopeR/0.99.31/autogenes/lib/python3.6/site-packages/numpy
#> numpy_version: 1.19.5
#>
#> NOTE: Python version was forced by use_python() function
#>
#> ──────────────────────────────────────────────────────────────────────────────