StatescopeR is an R wrapper around Statescope, a computational framework designed to discover cell states from cell type-specific gene expression profiles inferred from bulk RNAseq profiles. StatescopeR executes the three following steps consecutively: 1) Single cell RNAseq reference-based deconvolution 2) Refinement of inferred cell type-specific gene expression profiles 3) Cell state discovery. In addition to bulk RNAseq and a single cell RNAseq (scRNAseq) reference, StatescopeR can integrate prior expectations of (groups of) cell fractions in the bulk RNAseq. For example tumor fractions estimated by DNA copy number analysis or immune cell fractions by immunohistochemistry. Below is a manual for running StatescopeR with scRNAseq data from human pancreas cells 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 installing StatescopeR, the next step is to load it along with the scuttle and the scRNAseq packages. From the latter, the SegerstolpePancreas dataset will be used as an example. This dataset will serve two purposes: 1) Creating a scRNAseq reference. Alternatively, premade signatures can be obtained from https://github.com/tgac-vumc/StatescopeData using the fetch_signature function. 2) Creating pseudobulk for deconvolution, which can also be used for evaluation, as the true cell type fractions and expression are known. Usually, bulk RNA sequencing data, with unknown cell type fractions and expression, would be used for deconvolution.
## Load StatescopeR & scRNAseq (for example data)
suppressMessages(library(StatescopeR))
suppressMessages(library(scRNAseq))
suppressMessages(library(scuttle))
## Load the SegerstolpePancreas data set
scRNAseq <- SegerstolpePancreasData()
## remove duplicate genes
scRNAseq <- scRNAseq[!duplicated(rownames(scRNAseq)), ]
## 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 RNAseq to be deconvolved 3) selected genes for deconvolution, by default this uses AutoGeneS (optional) Prior expectation of fractions. See Group expectations for expectations of groups of cell types. (optional) Hyperparameters (i.e. cores for parallelization , nrepfinal for maximum number of optimization iterations etc.)
## Normalize (cp10k) and logtransform scRNAseq
cpm(scRNAseq) <- calculateCPM(scRNAseq)
logcounts(scRNAseq) <- log1p(cpm(scRNAseq) / 100)
## Create scRNAseq reference/signature with 50 hvg for quick example
signature <- create_signature(scRNAseq,
hvg_genes = TRUE, n_hvg_genes = 50L,
labels = scRNAseq$`cell type`
)
## select subset of genes for deconvolution (30/50 hvg to make it quick)
selected_genes <- select_genes(scRNAseq, 30L, 50L,
labels = scRNAseq$`cell type`
)
#> Using Python: /home/biocbuild/.pyenv/versions/3.7.17/bin/python3.7
#> Creating virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes' ...
#> + /home/biocbuild/.pyenv/versions/3.7.17/bin/python3.7 -m venv /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes
#> Done!
#> Installing packages: pip, wheel, setuptools
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes/bin/python -m pip install --upgrade pip wheel setuptools
#> Installing packages: 'autogenes==1.0.4', 'pandas==1.3.5', 'anndata==0.8.0', 'numpy==1.21.6', 'dill==0.3.7', 'deap==1.4.3', 'scipy==1.7.3', 'cachetools==5.5.2', 'scikit-learn==1.0.2', 'matplotlib==3.5.3'
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes/bin/python -m pip install --upgrade --no-user 'autogenes==1.0.4' 'pandas==1.3.5' 'anndata==0.8.0' 'numpy==1.21.6' 'dill==0.3.7' 'deap==1.4.3' 'scipy==1.7.3' 'cachetools==5.5.2' 'scikit-learn==1.0.2' 'matplotlib==3.5.3'
#> Virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes' successfully created.
## 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 for ductal cells
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
Once all input data have been prepared, the (pseudo)bulk RNAseq can be deconvolved to identify 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 functions Statescope.
## Run Deconvolution module
Statescope <- BLADE_deconvolution(signature, pseudobulk, selected_genes, prior)
#> Using Python: /home/biocbuild/.pyenv/versions/3.11.14/bin/python3.11
#> Creating virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/deconvolution' ...
#> + /home/biocbuild/.pyenv/versions/3.11.14/bin/python3.11 -m venv /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/deconvolution
#> Done!
#> Installing packages: pip, wheel, setuptools
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/deconvolution/bin/python -m pip install --upgrade pip wheel setuptools
#> Installing packages: 'numba==0.62.1', 'pandas==1.5.3', 'joblib==1.5.2', 'numpy==1.26.4', 'scipy==1.16.3', 'scikit-learn==1.5.2', 'matplotlib==3.6.3', 'mkl==2025.3.0', 'torch==2.9.1', 'dill==0.3.4'
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/deconvolution/bin/python -m pip install --upgrade --no-user 'numba==0.62.1' 'pandas==1.5.3' 'joblib==1.5.2' 'numpy==1.26.4' 'scipy==1.16.3' 'scikit-learn==1.5.2' 'matplotlib==3.6.3' 'mkl==2025.3.0' 'torch==2.9.1' 'dill==0.3.4'
#> Virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/deconvolution' successfully created.
## 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)
#> Using Python: /home/biocbuild/.pyenv/versions/3.11.14/bin/python3.11
#> Creating virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/statescope' ...
#> + /home/biocbuild/.pyenv/versions/3.11.14/bin/python3.11 -m venv /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/statescope
#> Done!
#> Installing packages: pip, wheel, setuptools
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/statescope/bin/python -m pip install --upgrade pip wheel setuptools
#> Installing packages: 'numba==0.62.1', 'pandas==1.5.3', 'joblib==1.5.2', 'numpy==1.23.5', 'scipy==1.15.3', 'matplotlib==3.6.3', 'seaborn==0.13.2', 'scikit-learn==1.5.2'
#> + /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/statescope/bin/python -m pip install --upgrade --no-user 'numba==0.62.1' 'pandas==1.5.3' 'joblib==1.5.2' 'numpy==1.23.5' 'scipy==1.15.3' 'matplotlib==3.6.3' 'seaborn==0.13.2' 'scikit-learn==1.5.2'
#> Virtual environment '/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/statescope' successfully created.
Since pseudobulk data derived from scRNAseq are used in this manual, it is possible to verify the accuracy of the StatescopeR output.
## 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 running StatescopeR, a wide range of analyses can be performed. In this section, we describe how to access the results and show some examples of visualizations with them.
## Show predicted fractions
metadata(Statescope)$fractions
#> DataFrame with 6 rows and 10 columns
#> H1 H2 H3 H4 H5 H6
#> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
#> gamma cell 0.1557854 0.147386 0.116859 0.1475245 0.1256517 0.1641532
#> alpha cell 0.3005740 0.309393 0.185046 0.3402188 0.2774545 0.3100585
#> beta cell 0.1627064 0.170174 0.120897 0.1696483 0.1252289 0.1767926
#> acinar cell 0.1747952 0.143390 0.195112 0.1374128 0.0967196 0.1265545
#> ductal cell 0.0467147 0.082192 0.285222 0.0449931 0.2353011 0.0744429
#> delta cell 0.1594244 0.147465 0.096864 0.1602024 0.1396442 0.1479983
#> T2D1 T2D2 T2D3 T2D4
#> <numeric> <numeric> <numeric> <numeric>
#> gamma cell 0.1880705 0.1131061 0.0852324 0.132946
#> alpha cell 0.3493517 0.2261601 0.1868113 0.271618
#> beta cell 0.1659803 0.0985992 0.0832572 0.154354
#> acinar cell 0.1232199 0.1453575 0.0971355 0.127659
#> ductal cell 0.0109012 0.3052454 0.4657599 0.170218
#> delta cell 0.1624762 0.1115317 0.0818036 0.143204
## Show predicted ductal cell type specific gene expression profiles
assay(metadata(Statescope)$ct_specific_gep$`ductal cell`, "weighted_gep")
#> DataFrame with 50 rows and 10 columns
#> H1 H2 H3 H4 H5 H6
#> <numeric> <numeric> <numeric> <numeric> <numeric> <numeric>
#> GCG -0.209496 -0.209491 -0.197886 -0.209661 -0.208636 -0.209491
#> TTR -0.444675 -0.444587 -0.440269 -0.444894 -0.443817 -0.444571
#> INS -0.320821 -0.320540 -0.318760 -0.320757 -0.319193 -0.320696
#> PPY -0.280438 -0.280413 -0.279423 -0.280036 -0.279870 -0.280474
#> REG1A 0.677720 0.706514 0.683367 0.740765 0.496637 0.770259
#> ... ... ... ... ... ... ...
#> PNLIP -0.1859452 -0.1849785 -0.1858142 -0.1856118 -0.1922385 -0.1849097
#> SERPING1 0.0605501 0.1899599 0.1452673 0.1487681 0.1968191 0.1427350
#> CFTR -0.0792104 0.1428833 0.0438710 0.0851317 0.0509284 0.1764416
#> ABCC8 -0.0971111 -0.0971003 -0.0969924 -0.0971142 -0.0970550 -0.0971016
#> CTSD 0.2807088 0.3094438 0.5529918 0.3981679 0.4716538 -0.0178438
#> T2D1 T2D2 T2D3 T2D4
#> <numeric> <numeric> <numeric> <numeric>
#> GCG -0.209770 -0.207699 -0.205704 -0.208810
#> TTR -0.445024 -0.441274 -0.437029 -0.443693
#> INS -0.321006 -0.305197 -0.281531 -0.320110
#> PPY -0.280584 -0.276617 -0.271524 -0.280023
#> REG1A 0.688349 0.781995 0.999680 0.812959
#> ... ... ... ... ...
#> PNLIP -0.1861397 -0.1813423 -0.1745991 -0.183668
#> SERPING1 0.0658480 0.1771080 0.2021813 0.133801
#> CFTR 0.0234802 0.1296970 0.1960245 0.169654
#> ABCC8 -0.0971272 -0.0969133 -0.0966557 -0.097064
#> CTSD 0.3082869 0.4522108 0.3825589 0.300022
## Show ductal cell state scores per sample
metadata(Statescope)$statescores$`ductal cell`
#> DataFrame with 10 rows and 2 columns
#> V1 V2
#> <numeric> <numeric>
#> H1 3.26191e-03 9.96738e-01
#> H2 8.09262e-01 1.90738e-01
#> H3 2.19928e-07 1.00000e+00
#> H4 3.59346e-01 6.40654e-01
#> H5 8.06776e-02 9.19322e-01
#> H6 1.00000e+00 3.25940e-09
#> T2D1 2.65436e-01 7.34564e-01
#> T2D2 6.78607e-01 3.21393e-01
#> T2D3 9.98755e-01 1.24452e-03
#> T2D4 8.76205e-01 1.23795e-01
## Show ductal cell state loadings
metadata(Statescope)$stateloadings$`ductal cell`
#> DataFrame with 50 rows and 2 columns
#> V1 V2
#> <numeric> <numeric>
#> GCG -0.190864 -0.183576
#> TTR -0.405330 -0.395745
#> INS -0.281834 -0.285590
#> PPY -0.254393 -0.250040
#> REG1A 0.765300 0.571998
#> ... ... ...
#> PNLIP -0.1663051 -0.1674517
#> SERPING1 0.1561813 0.1166340
#> CFTR 0.1546821 0.0163687
#> ABCC8 -0.0889382 -0.0866923
#> CTSD 0.2465490 0.3845163
## Plot heatmap of fractions
fraction_heatmap(Statescope)
## Plot barplot of top state loadings
barplot_stateloadings(Statescope)
When an expected fraction is available for a group of cell types rather than for individual cell types (i.e. Lymphoid cells instead of T/B cells), StatescopeR allows the use of a group prior. Here we show an example where the expected fraction of the pancreatic islet cells (alpha, beta, gamma & delta cells) together is known, but not their individual contributions.
## define cell types in one group
grouped_cts <- c("alpha cell", "beta cell", "gamma cell", "delta 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
)
StatescopeRIf you made use of StatescopeR for your research, please use the following information to cite the package and the underlying methodology.
## 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) (2026-01-15 r89304)
#> 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 2026-01-28
#> 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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.matrix 1.11.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.ranges 1.11.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.sce 1.11.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.schemas 1.11.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> alabaster.se 1.11.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationDbi 1.73.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationFilter 1.35.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> AnnotationHub 4.1.0 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> beachmat 2.27.2 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocFileCache 3.1.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocGenerics * 0.57.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocIO 1.21.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocManager 1.30.27 2025-11-14 [2] CRAN (R 4.6.0)
#> BiocNeighbors 2.5.2 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocParallel 1.45.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocSingular 1.27.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocStyle * 2.39.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> BiocVersion 3.23.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> Biostrings 2.79.4 2026-01-28 [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.3.0 2026-01-14 [2] CRAN (R 4.6.0)
#> bluster 1.21.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> bookdown 0.46 2025-12-05 [2] CRAN (R 4.6.0)
#> bslib 0.10.0 2026-01-26 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> circlize 0.4.17 2025-12-08 [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 2026-01-28 [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 2026-01-28 [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 2026-01-28 [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.2 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> ensembldb 2.35.0 2026-01-28 [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 2026-01-28 [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.2 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicAlignments 1.47.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicFeatures 1.63.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> GenomicRanges * 1.63.1 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> h5mread 1.3.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> HDF5Array 1.39.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> htmltools 0.5.9 2025-12-04 [2] CRAN (R 4.6.0)
#> httr 1.4.7 2023-08-15 [2] CRAN (R 4.6.0)
#> httr2 1.2.2 2025-12-08 [2] CRAN (R 4.6.0)
#> igraph 2.2.1 2025-10-27 [2] CRAN (R 4.6.0)
#> IRanges * 2.45.0 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> knitr 1.51 2025-12-20 [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.5 2026-01-08 [2] CRAN (R 4.6.0)
#> limma 3.67.0 2026-01-28 [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 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> otel 0.2.0 2025-08-29 [2] CRAN (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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> R6 2.6.1 2025-02-15 [2] CRAN (R 4.6.0)
#> rappdirs 0.3.4 2026-01-17 [2] CRAN (R 4.6.0)
#> RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.6.0)
#> Rcpp 1.1.1 2026-01-10 [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.12 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> rhdf5filters 1.23.3 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> Rhdf5lib 1.33.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> rjson 0.2.23 2024-09-16 [2] CRAN (R 4.6.0)
#> rlang 1.1.7 2026-01-09 [2] CRAN (R 4.6.0)
#> rmarkdown 2.30 2025-09-28 [2] CRAN (R 4.6.0)
#> Rsamtools 2.27.0 2026-01-28 [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.3 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> S4Arrays 1.11.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> S4Vectors * 0.49.0 2026-01-28 [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 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> scRNAseq * 2.25.0 2026-01-27 [2] Bioconductor 3.23 (R 4.6.0)
#> scuttle * 1.21.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> Seqinfo * 1.1.0 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> SparseArray 1.11.10 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> StatescopeR * 0.99.34 2026-01-28 [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 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> tibble 3.3.1 2026-01-11 [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.1 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> vctrs 0.7.1 2026-01-23 [2] CRAN (R 4.6.0)
#> withr 3.0.2 2024-10-28 [2] CRAN (R 4.6.0)
#> xfun 0.56 2026-01-18 [2] CRAN (R 4.6.0)
#> XML 3.99-0.20 2025-11-08 [2] CRAN (R 4.6.0)
#> xml2 1.5.2 2026-01-17 [2] CRAN (R 4.6.0)
#> XVector 0.51.0 2026-01-28 [2] Bioconductor 3.23 (R 4.6.0)
#> yaml 2.3.12 2025-12-10 [2] CRAN (R 4.6.0)
#>
#> [1] /tmp/RtmpgGDme8/Rinst35ecb653d1bd0
#> [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.34/autogenes/bin/python
#> libpython: /home/biocbuild/.pyenv/versions/3.7.17/lib/libpython3.7m.so
#> pythonhome: /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes:/var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes
#> version: 3.7.17 (default, Dec 24 2025, 20:54:15) [GCC 13.3.0]
#> numpy: /var/cache/basilisk/1.23.0/StatescopeR/0.99.34/autogenes/lib/python3.7/site-packages/numpy
#> numpy_version: 1.21.6
#>
#> NOTE: Python version was forced by use_python() function
#>
#> ──────────────────────────────────────────────────────────────────────────────