stPipe
packagestPipe::Run_ST
stPipe::Run_Loc_Match
stPipe::Run_QC
stPipe::Run_Clustering
stPipe::Run_Interactive
stPipe::Run_Visualization
stPipe
packagestPipe
is a package designed to process sequencing-based spatial transcriptomics (sST) data generated by different platforms, including 10X Visium, BGI Stereo-seq, and Slide-seq/Curio-seeker. stPipe
consists of two major components. The first is data pre-processing with paired-end FASTQ files as input and a gene count matrix with matched spatial information as output. The second component takes the pre-processing results as input and performs quality control, downstream object construction and generates an HTML report for overall summary.
10X Genomics’ Visium platform is one of the most widely adopted methods for sST, offering high-throughput capabilities and compatibility with diverse tissue types. Visium is developed based on the original spatial transcriptomics technology via positioning histological sections on arrayed reverse transcription primers with unique positional barcodes. Visium supports both fresh frozen (FF) and formalin-fixed, paraffin-embedded (FFPE) tissues. Tissue processing can be performed manually or via the automated CytAssist system. Two distinct mRNA capture methodologies are available: probe-based and poly-A-based methods. Visium offers two capture area configurations: 6.5 x 6.5 mm area comprising 4,992 spots and an 11 x 11 mm area containing 14,336 spots. Each spot is 55 µm in diameter, with a center-to-center distance of 100 µm. Recent available Visium HD provides an enhanced spatial resolution of 2 µm without any gaps between.
Slide-seq uses a method whereby RNA from tissue sections is transferred onto a surface embedded with DNA-barcoded 10 µm beads. The random distribution of beads across the surface results in variable distances between barcodes. Notably, the spatial barcode in Slide-seq is divided into two segments by a linker sequence. Slide-seq V2 modifies the library generation, bead synthesis, and array indexing steps to improve the mRNA capture sensitivity by approximately 10-fold over the original Slide-seq V1 method. Slide-seq has been commercialized as Curio Seeker by Curio Bioscience.
Based on DNA nanoball (DNB) technology, BGI Stereo-seq (SpaTial Enhanced REsolution Omics-sequencing) resolves the sequencing data via 1 x 1 cm RNA capture unit on the STOmics-GeneExpression Chip. The original DNA nanoball sequencing amplifies small fragments of genomic DNA into DNA nanoballs via rolling circle replication with fluorescent nucleotides binding to complementary nucleotides that are polymerized to anchor sequences bound to known sequences on the DNA template. In the refined Stereo-seq protocol, DNBs are generated by rolling circle amplification of circularized DNA oligonucleotides, forming spherical structures. Each DNB carries a 25-nt coordinate identity (CID) used to identify its spatial position during RNA capture and sequencing. Then the spatial location of each DNB is recorded with diameter of 0.2um and a center-to-center distance of 0.5um which allows for sub-cellular single-cell resolution data.
Protocol | Resolution (micron) | Spot distance (micron) | Typical array size |
---|---|---|---|
Visium | 50 | 100 | 0.65cm x 0.65cm / 1.1cm x 1.1cm |
Visium HD | 2 | 0 (no gaps) | 0.65cm x 0.65cm / 1.1cm x 1.1cm |
Slide-seq/Curio-seeker | 10 | 10 | / |
Stereo-seq | 10 | 0.5 | 1cm x 1cm; (max: 13.2cm x 13.2cm) |
Curio-seeker | 10 | 10 | 0.3cm x 0.3cm / 1cm x 1cm |
The sample data used for this study is from a FFPE mouse spleen performed using 10x Visium’s probe based assay. It has been downsampled for computing efficiency.
All files can be found in the extdata
folder of the stPipe
package. The original data comes from the SpatialBenchVisium
project and is available under GEO accession number GSE254652, detailed in Spotlight on 10x Visium: a multi-sample protocol comparison of spatial technologies.
We begin by filling in specific details of this particular sample in the config file. Note it is important to read the explanation inside the sample config file as we will be using it for the majority of stPipe
implemented functions.
stPipe
implements different strategies to cope with corresponding platform features, hence the config file needs to be filled in correctly before running stPipe
. This will be further explained in the following sections.
library(stPipe)
config <- system.file("config", "config_stPipe.yml", package = "stPipe")
yaml_content <- readLines(config)
cat(yaml_content, sep = "\n")
## ########## Config File ##########
## ### 'Run_ST' ###
## # Directory stores paired FASTQ files, for multiple samples should be separated by ','
## data_directory: /path/to/FASTQ
## # Directory stores outputs, for multiple samples should be separated by ','
## output_directory: /path/to/output
## species: mouse # "human"; "mouse"
## index_fa: /path/to/genome.fa # for polyA-based protocol only
## index_gff3: /path/to/anno.gff3 # for polyA-based protocol only
## # read_structure design, positions are 0-indexed so the first base is considered base 0
## bs1: -1 # barcode start position in fq_R1, -1 indicates no barcode
## bl1: 0 # barcode length in fq_R1, 0 since no barcode present
## bs2: 0 # barcode start position in fq_R2
## bl2: 16 # barcode length in fq_R2
## us: 16 # UMI start position in fq_R2
## ul: 12 # UMI length
## ll: 0 # for Curio-seeker only - linker sequence length
## h5_mapping_path: /path/to/barcodeToPos.h5 # for Stereoseq only - path to mapping h5 file
## bin_size: 1 # for Stereoseq only - binning size of n * n, 1 means no binning performed
## # Rsubread::align
## scpipe_nthreads: 4 # used thread
## # scPipe::sc_detect_bc
## max_reads: 1000000000 # process first N reads
## min_count: 100 # discard barcodes with few than 100 hits
## number_of_locations: 4992 # number of mapped spatial location
## technology_version: "Visium_probe_v1" # “Visium_probe_vn”; “Visium_polyA”; “Stereoseq”; “Slideseq”; "Curio-seeker"
##
## ### 'Run_loc_match' ###
## visium_coordination: "V4" # for Visium only - chip chemistry version
## image_path: /path/to/Image.tif # for Visium only - path to image file
## bead_location: /path/to/bead_locations.csv # for Slideseq/Curio-seeker only - path to sample coordination file
##
## ### 'Run_QC' ###
## qc_filter: EmptyDropletUtils # 'slope_max' [exact raw UMI count filter] or 'EmptyDropletUtils' [DropletUtils::emptyDrops]
## qc_per: 0.1_0.8 # threshold for calculating max slope [i.e., 0.4 - first 40%]
stPipe::Run_ST
The implemented Run_ST
function streamlines multiple steps into one cohesive process, including FASTQ or BAM file reformatting where barcodes and UMI information are extracted and incorporated into the read headers, read alignment, barcode demultiplexing, and gene counting. stPipe
performs different strategies for probe-based and polyA-based protocols, the former one utilizes probe set version to construct annotation files such as fa and GFF3 files while the latter one relies the annotation files specified by the user.
The read structure for 10X Visium has the 8 bp long spatial barcode starting at position 6 and the 6 bp long UMI sequence starting at the first position. So the read structure will be : list(bs1=-1, bl1=0, bs2=6, bl2=8, us=0, ul=6)
. bs1=-1, bl1=0
means we don’t have an index in read 1 so we set a negative value as its start position and give it zero length. bs2=6, bl2=8
means we have an index in read two which starts at position 6 in the read and is 8 bases in length. us=0, ul=6
means we have a UMI at the start of read 2 which is 6 bases long. NOTE: we use a zero based index system, so the indexing of the sequence starts at zero.
We begin by creating a temporary folder to download the demo data.
# FASTQ file path:
data_directory <- file.path(tempdir(), "stPipe_demo_data")
if (!dir.exists(data_directory)) {
dir.create(data_directory, recursive = TRUE)
}
options(timeout = 600)
R1_url <- "https://zenodo.org/records/14920583/files/stPipe_demo_R1.fq.gz?download=1"
R2_url <- "https://zenodo.org/records/14920583/files/stPipe_demo_R2.fq.gz?download=1"
R1_path <- file.path(data_directory, "stPipe_demo_R1.fq.gz")
R2_path <- file.path(data_directory, "stPipe_demo_R2.fq.gz")
download.file(R1_url, R1_path)
download.file(R2_url, R2_path)
config <- system.file("config", "config_stPipe.yml", package = "stPipe")
config <- yaml::read_yaml(config)
config$data_directory <- data_directory
We then create a temporary directory to store results before we’re ready to run Run_ST and adjust some config information for this sample. NOTE: In most functions there will be a flag called show.config
which simply prints all config information and helps you make sure that everything is configured correctly and reduces the risk of encountering issues later.
output_directory <- file.path(tempdir(), "stPipe_output")
if (!dir.exists(output_directory)) {
dir.create(output_directory, recursive = TRUE)
}
config$output_directory <- output_directory
demo_config_path <- tempfile(fileext = ".yml")
yaml::write_yaml(config, demo_config_path)
Now stPipe
is ready to run.
Run_ST(config = demo_config_path, show.config = FALSE)
## Using the following FASTQ files:
## /tmp/RtmpPl7sg7/stPipe_demo_data/stPipe_demo_R1.fq.gz
## /tmp/RtmpPl7sg7/stPipe_demo_data/stPipe_demo_R2.fq.gz
##
##
## trimming fastq file...
## pass QC: 557592
## removed_have_N: 12408
## removed_low_qual: 0
## time elapsed: 635 milliseconds
## FASTA file generated based on the provided probe set:/tmp/RtmpPl7sg7/stPipe_output/probe_set.fa
## GFF3 file generated based on the provided probe set:/tmp/RtmpPl7sg7/stPipe_output/probe_set.gff3
##
## ========== _____ _ _ ____ _____ ______ _____
## ===== / ____| | | | _ \| __ \| ____| /\ | __ \
## ===== | (___ | | | | |_) | |__) | |__ / \ | | | |
## ==== \___ \| | | | _ <| _ /| __| / /\ \ | | | |
## ==== ____) | |__| | |_) | | \ \| |____ / ____ \| |__| |
## ========== |_____/ \____/|____/|_| \_\______/_/ \_\_____/
## Rsubread 2.23.4
##
## //================================= setting ==================================\\
## || ||
## || Index name : stPipe_output ||
## || Index space : base space ||
## || Index split : no-split ||
## || Repeat threshold : 100 repeats ||
## || Gapped index : no ||
## || ||
## || Free / total memory : 83.8GB / 125.4GB ||
## || ||
## || Input files : 1 file in total ||
## || o probe_set.fa ||
## || ||
## \\============================================================================//
##
## //================================= Running ==================================\\
## || ||
## || Check the integrity of provided reference sequences ... ||
## || No format issues were found ||
## || Scan uninformative subreads in reference sequences ... ||
## || Estimate the index size... ||
## || 3.0 GB of memory is needed for index building. ||
## || Build the index... ||
## || Save current index block... ||
## || [ 0.0% finished ] ||
## || [ 10.0% finished ] ||
## || [ 20.0% finished ] ||
## || [ 30.0% finished ] ||
## || [ 40.0% finished ] ||
## || [ 50.0% finished ] ||
## || [ 60.0% finished ] ||
## || [ 70.0% finished ] ||
## || [ 80.0% finished ] ||
## || [ 90.0% finished ] ||
## || [ 100.0% finished ] ||
## || ||
## || Total running time: 0.2 minutes. ||
## || Index /tmp/RtmpPl7sg7/stPipe_output was successfully built. ||
## || ||
## \\============================================================================//
##
##
## ========== _____ _ _ ____ _____ ______ _____
## ===== / ____| | | | _ \| __ \| ____| /\ | __ \
## ===== | (___ | | | | |_) | |__) | |__ / \ | | | |
## ==== \___ \| | | | _ <| _ /| __| / /\ \ | | | |
## ==== ____) | |__| | |_) | | \ \| |____ / ____ \| |__| |
## ========== |_____/ \____/|____/|_| \_\______/_/ \_\_____/
## Rsubread 2.23.4
##
## //================================= setting ==================================\\
## || ||
## || Function : Read alignment (RNA-Seq) ||
## || Input file : trimmed.fastq ||
## || Output file : out.aln.bam (BAM) ||
## || Index name : stPipe_output ||
## || ||
## || ------------------------------------ ||
## || ||
## || Threads : 4 ||
## || Phred offset : 33 ||
## || Min votes : 3 / 10 ||
## || Max mismatches : 3 ||
## || Max indel length : 5 ||
## || Report multi-mapping reads : yes ||
## || Max alignments per multi-mapping read : 1 ||
## || ||
## \\============================================================================//
##
## //=============== Running (12-Sep-2025 20:20:21, pid=2136283) ================\\
## || ||
## || Check the input reads. ||
## || The input file contains base space reads. ||
## || Initialise the memory objects. ||
## || Estimate the mean read length. ||
## || The range of Phred scores observed in the data is [2,34] ||
## || Create the output BAM file. ||
## || Check the index. ||
## || Init the voting space. ||
## || Global environment is initialised. ||
## || Load the 1-th index block... ||
## || The index block has been loaded. ||
## || Start read mapping in chunk. ||
## || 4% completed, 0.6 mins elapsed, rate=94.7k reads per second ||
## || 11% completed, 0.6 mins elapsed, rate=84.1k reads per second ||
## || 18% completed, 0.6 mins elapsed, rate=80.1k reads per second ||
## || 24% completed, 0.6 mins elapsed, rate=76.6k reads per second ||
## || 31% completed, 0.6 mins elapsed, rate=79.5k reads per second ||
## || 37% completed, 0.6 mins elapsed, rate=78.4k reads per second ||
## || 43% completed, 0.7 mins elapsed, rate=78.4k reads per second ||
## || 50% completed, 0.7 mins elapsed, rate=78.4k reads per second ||
## || 57% completed, 0.7 mins elapsed, rate=79.6k reads per second ||
## || 63% completed, 0.7 mins elapsed, rate=80.2k reads per second ||
## || 70% completed, 0.7 mins elapsed, rate=9.4k reads per second ||
## || 73% completed, 0.7 mins elapsed, rate=9.7k reads per second ||
## || 76% completed, 0.7 mins elapsed, rate=10.1k reads per second ||
## || 79% completed, 0.7 mins elapsed, rate=10.4k reads per second ||
## || 82% completed, 0.7 mins elapsed, rate=10.7k reads per second ||
## || 86% completed, 0.7 mins elapsed, rate=11.0k reads per second ||
## || 89% completed, 0.7 mins elapsed, rate=11.3k reads per second ||
## || 92% completed, 0.7 mins elapsed, rate=11.6k reads per second ||
## || 95% completed, 0.7 mins elapsed, rate=11.9k reads per second ||
## || 98% completed, 0.8 mins elapsed, rate=12.2k reads per second ||
## || ||
## || Completed successfully. ||
## || ||
## \\==================================== ====================================//
##
## //================================ Summary =================================\\
## || ||
## || Total reads : 557,592 ||
## || Mapped : 471,616 (84.6%) ||
## || Uniquely mapped : 471,598 ||
## || Multi-mapping : 18 ||
## || ||
## || Unmapped : 85,976 ||
## || ||
## || Indels : 824 ||
## || ||
## || Running time : 0.8 minutes ||
## || ||
## \\============================================================================//
##
## adding annotation files...
## time elapsed: 79 milliseconds
##
## annotating exon features...
## updating progress every 3 minutes...
## 557592 reads processed, 185k reads/sec
## number of read processed: 557592
## unique map to exon: 471616 (84.58%)
## ambiguous map to multiple exon: 0 (0.00%)
## map to intron: 0 (0.00%)
## not mapped: 0 (0.00%)
## unaligned: 85976 (15.42%)
## adding annotation files...
## time elapsed: 77 milliseconds
##
## annotating exon features...
## updating progress every 3 minutes...
## 557592 reads processed, 69k reads/sec
## number of read processed: 557592
## unique map to exon: 471616 (84.58%)
## ambiguous map to multiple exon: 0 (0.00%)
## map to intron: 0 (0.00%)
## not mapped: 0 (0.00%)
## unaligned: 85976 (15.42%)
## demultiplexing reads by barcode...
## Warning: mitochondrial chromosome not found using chromosome name `MT`.
## time elapsed: 2905 milliseconds
##
## summarising gene counts...
## time elapsed: 2627 milliseconds
## [[1]]
## NULL
stPipe::Run_Loc_Match
The Run_Loc_Match
function is responsible for pairing the gene expression matrix with match spatial locations. The key input in the config file is technology_version which specifies the used platform. In further segmentation, each platform has its own argument: Visium needs visium_coordinates to match internally stored spatial coordinates system and image_path is optional for pixel computation; Slide-seq/Curio-seeker needs bead_location which should be a .csv
file with three columns, where the first column has the spatial barcode sequence, the second column contains the X coordinate, and the third column contains the Y coordinate to match bead location across the spatial dimension. Stereo-seq only requires technology_version.
st.result <- Run_Loc_Match(config = demo_config_path, pixel = FALSE, show.config = FALSE)
## Heading Coordinates for Visium:
## c("AACACCTACTATCGAA", "AACACGTGCATCGCAC", "AACACTTGGCAAGGAA", "AACAGGAAGAGCATAG", "AACAGGATTCATAGTT", "AACAGGCCAACGATTA")c(123, 23, 72, 8, 44, 128)c(1, 77, 48, 70, 50, 72)
## Start mapping coordination...
head(st.result$matched_data)
## barcode_sequence spatial_name count X_coordinate Y_coordinate UMI_count
## 1 AACACTTGGCAAGGAA SPATIAL_00559 316 72 48 285
## 2 AACAGGATTCATAGTT SPATIAL_01115 227 44 50 200
## 3 AACAGGTTCACCGAAG SPATIAL_01408 196 42 52 170
## 4 AACAGTCAGGCTCCGC SPATIAL_00282 391 7 25 336
## 5 AACAGTCCACGCGGTG SPATIAL_01114 227 11 13 186
## 6 AACATCTTAAGGCTCA SPATIAL_00720 282 9 17 229
stPipe::Run_QC
At this stage, the pre-processing work is almost done. But under specific case, you might want to do some QC filtering, for example, to clean up background noise in 10X Visium scenario. The Run_QC
function helps with this purpose.
# adjust some parameter and update the config file
config$qc_per <- c("0.4_0.6")
yaml::write_yaml(config, demo_config_path)
Under this filtering criteria, more than 65% of the original spots are lost which might be concerning.
qc.results <- Run_QC(config = demo_config_path, matched.data = st.result$matched_data, gene.matrix = st.result$gene_count, show.config = FALSE)
dim(qc.results$filtered_spatial_coords)[1]/dim(st.result$matched_data)[1]
## [1] 0.3329923
stPipe::Run_Clustering
Once you are satisfied with the data, the Run_Clustering
function can be utilized to perform some preliminary clustering analysis, including UMAP and T-SNE. The clustering function itself prints message every 50 iteration, to avoid long printout only 15 are shown.
output <- capture.output(
Cluster.results <- Run_Clustering(
gene.count = st.result$gene_count,
matched.data = st.result$matched_data,
num_clusters = 5
)
)
cat(output[seq_along(1:15)], sep = "\n")
## Performing PCA
## Read the 1955 x 50 data matrix successfully!
## OpenMP is working. 1 threads.
## Using no_dims = 2, perplexity = 30.000000, and theta = 0.500000
## Computing input similarities...
## Building tree...
## Done in 0.34 seconds (sparsity = 0.069427)!
## Learning embedding...
## Iteration 50: error is 73.857527 (50 iterations in 0.37 seconds)
## Iteration 100: error is 66.439156 (50 iterations in 0.37 seconds)
## Iteration 150: error is 66.211237 (50 iterations in 0.39 seconds)
## Iteration 200: error is 66.190734 (50 iterations in 0.39 seconds)
## Iteration 250: error is 66.190540 (50 iterations in 0.39 seconds)
## Iteration 300: error is 1.718403 (50 iterations in 0.31 seconds)
## Iteration 350: error is 1.524888 (50 iterations in 0.30 seconds)
stPipe::Run_Interactive
Previously the Run_QC
function serves as a step to get rid of spatial locations with low quality reads, but sometimes you might want to select a specific spatial region and perform further analysis. In this case, consider to use the Run_Interactive
function which provides you with an interactive R shiny interface for personalized selection.
Run_Interactive(
matched_data = st.result$matched_data,
clustering_result = Cluster.results$tsne_df,
background_img = NULL,
reduction_method = "tsne",
point_size = 1
)
stPipe::Run_Visualization
Run_Visualization
is the step where you want to check how the data looks like after a series of pre-processing steps. This function contains two flags for visualization, Vis.spatial
and Vis.read
corresponds to spatial visualization and sequencing visualization respectively. The former one contains UMI heatmap in raw value and log format. The latter includes barcode demultiplexing statistics, distribution of UMI duplication numbers, and alignment statistics.
Vis.results <- Run_Visualization(
matched.data = st.result$matched_data,
config = demo_config_path,
Vis.spatial = TRUE,
Vis.read = TRUE,
show.config = FALSE
)
## organism/gene_id_type not provided. Make a guess:mmusculus_gene_ensembl/ensembl_gene_id
## no ERCC spike-in. Skip `non_ERCC_percent`
## Warning: `aes_string()` was deprecated in ggplot2 3.0.0.
## ℹ Please use tidy evaluation idioms with `aes()`.
## ℹ See also `vignette("ggplot2-in-packages")` for more information.
## ℹ The deprecated feature was likely used in the scPipe package.
## Please report the issue at <https://github.com/LuyiTian/scPipe>.
## This warning is displayed once every 8 hours.
## Call `lifecycle::last_lifecycle_warnings()` to see where this warning was
## generated.
Barcode demultiplexing statistics this can suggest the percentage of reads that uniquely match to the spatial barcodes, as well as the unmatched proportion of reads and their alignment rates to introns and exons.
Vis.results$demultiplex_plot
Mapping statistics this plot shows the mapping quality of sample. In this demo case, the 10X Visium data is generated as probe-based, hence either mapped to exon or unaligned will be observed. For polyA-based protocols, mixed variable should be expected.
Vis.results$mapping_plot
Spatial UMI Plot this is plotted as both raw and log UMI count format, suggesting the UMI spatial distribution across the captured tissue.
After all the previous work, a HTML report in RMarkDown format can be generated to summarize your data and save it for further documentation archive purposes. The Run_HTML
function can help with this. NOTE: You will need to include ALL previous outputs to ensure the successful operation of this function.
Upstream pre-processing is never the end of exploring data in research, the Run_Create_Obj
function is designed for further downstream analysis, catering for different pipelines such as Seurat, squidpy, and SpatialExperiment supported bioconductor packages.
library(Seurat)
## Loading required package: SeuratObject
## Loading required package: sp
##
## Attaching package: 'SeuratObject'
## The following objects are masked from 'package:base':
##
## intersect, t
library(SeuratObject)
seu <- Run_Create_Obj(
gene.matrix = st.result$gene_count,
matched.data = st.result$matched_data,
obj.type = 'Seurat',
tech = 'Visium'
)
## Warning: Feature names cannot have underscores ('_'), replacing with dashes
## ('-')
## Warning: Data is of class matrix. Coercing to dgCMatrix.
seu
## An object of class Seurat
## 14669 features across 1955 samples within 1 assay
## Active assay: Spatial (14669 features, 0 variable features)
## 1 layer present: counts
## 1 image present: image
Session Info
sessionInfo()
## R version 4.5.1 Patched (2025-08-23 r88802)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.22-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=en_US.UTF-8
## [3] LC_TIME=en_GB LC_COLLATE=C
## [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: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] Seurat_5.3.0 SeuratObject_5.2.0 sp_2.2-0 stPipe_0.99.994
## [5] BiocStyle_2.37.1
##
## loaded via a namespace (and not attached):
## [1] spatstat.sparse_3.1-0 matrixStats_1.5.0
## [3] bitops_1.0-9 fontawesome_0.5.3
## [5] httr_1.4.7 RColorBrewer_1.1-3
## [7] sctransform_0.4.2 tools_4.5.1
## [9] R6_2.6.1 HDF5Array_1.37.0
## [11] uwot_0.2.3 lazyeval_0.2.2
## [13] rhdf5filters_1.21.0 withr_3.0.2
## [15] gridExtra_2.3 prettyunits_1.2.0
## [17] progressr_0.15.1 cli_3.6.5
## [19] Biobase_2.69.0 spatstat.explore_3.5-2
## [21] fastDummies_1.7.5 labeling_0.4.3
## [23] sass_0.4.10 spatstat.data_3.1-8
## [25] S7_0.2.0 robustbase_0.99-6
## [27] pbapply_1.7-4 ggridges_0.5.7
## [29] askpass_1.2.1 Rsamtools_2.25.3
## [31] R.utils_2.13.0 Rsubread_2.23.4
## [33] dichromat_2.0-0.1 parallelly_1.45.1
## [35] limma_3.65.4 RSQLite_2.4.3
## [37] generics_0.1.4 BiocIO_1.19.0
## [39] spatstat.random_3.4-1 ica_1.0-3
## [41] crosstalk_1.2.2 dplyr_1.1.4
## [43] Matrix_1.7-4 S4Vectors_0.47.0
## [45] abind_1.4-8 R.methodsS3_1.8.2
## [47] lifecycle_1.0.4 yaml_2.3.10
## [49] edgeR_4.7.4 SummarizedExperiment_1.39.1
## [51] rhdf5_2.53.4 SparseArray_1.9.1
## [53] BiocFileCache_2.99.6 Rtsne_0.17
## [55] grid_4.5.1 blob_1.2.4
## [57] promises_1.3.3 dqrng_0.4.1
## [59] crayon_1.5.3 dir.expiry_1.17.0
## [61] miniUI_0.1.2 lattice_0.22-7
## [63] beachmat_2.25.5 cowplot_1.2.0
## [65] KEGGREST_1.49.1 magick_2.9.0
## [67] pillar_1.11.0 knitr_1.50
## [69] GenomicRanges_1.61.2 rjson_0.2.23
## [71] future.apply_1.20.0 codetools_0.2-20
## [73] glue_1.8.0 spatstat.univar_3.1-4
## [75] data.table_1.17.8 vctrs_0.6.5
## [77] png_0.1-8 spam_2.11-1
## [79] org.Mm.eg.db_3.21.0 gtable_0.3.6
## [81] cachem_1.1.0 xfun_0.53
## [83] S4Arrays_1.9.1 mime_0.13
## [85] DropletUtils_1.29.7 Seqinfo_0.99.2
## [87] survival_3.8-3 SingleCellExperiment_1.31.1
## [89] pbmcapply_1.5.1 tinytex_0.57
## [91] statmod_1.5.0 fitdistrplus_1.2-4
## [93] ROCR_1.0-11 nlme_3.1-168
## [95] bit64_4.6.0-1 progress_1.2.3
## [97] filelock_1.0.3 RcppAnnoy_0.0.22
## [99] bslib_0.9.0 irlba_2.3.5.1
## [101] KernSmooth_2.23-26 BiocGenerics_0.55.1
## [103] DBI_1.2.3 tidyselect_1.2.1
## [105] scPipe_2.9.0 bit_4.6.0
## [107] compiler_4.5.1 curl_7.0.0
## [109] httr2_1.2.1 h5mread_1.1.1
## [111] xml2_1.4.0 DelayedArray_0.35.2
## [113] plotly_4.11.0 bookdown_0.44
## [115] rtracklayer_1.69.1 scales_1.4.0
## [117] DEoptimR_1.1-4 lmtest_0.9-40
## [119] rappdirs_0.3.3 goftest_1.2-3
## [121] stringr_1.5.2 SpatialExperiment_1.19.1
## [123] digest_0.6.37 spatstat.utils_3.1-5
## [125] rmarkdown_2.29 basilisk_1.21.5
## [127] XVector_0.49.0 htmltools_0.5.8.1
## [129] pkgconfig_2.0.3 umap_0.2.10.0
## [131] sparseMatrixStats_1.21.0 MatrixGenerics_1.21.0
## [133] dbplyr_2.5.1 fastmap_1.2.0
## [135] rlang_1.1.6 htmlwidgets_1.6.4
## [137] shiny_1.11.1 DelayedMatrixStats_1.31.0
## [139] farver_2.1.2 jquerylib_0.1.4
## [141] zoo_1.8-14 jsonlite_2.0.0
## [143] BiocParallel_1.43.4 mclust_6.1.1
## [145] R.oo_1.27.1 RCurl_1.98-1.17
## [147] magrittr_2.0.4 scuttle_1.19.0
## [149] dotCall64_1.2 patchwork_1.3.2
## [151] Rhdf5lib_1.31.0 Rcpp_1.1.0
## [153] reticulate_1.43.0 stringi_1.8.7
## [155] MASS_7.3-65 plyr_1.8.9
## [157] org.Hs.eg.db_3.21.0 parallel_4.5.1
## [159] listenv_0.9.1 ggrepel_0.9.6
## [161] deldir_2.0-4 Biostrings_2.77.2
## [163] splines_4.5.1 tensor_1.5.1
## [165] hms_1.1.3 locfit_1.5-9.12
## [167] igraph_2.1.4 spatstat.geom_3.5-0
## [169] RcppHNSW_0.6.0 reshape2_1.4.4
## [171] biomaRt_2.65.7 stats4_4.5.1
## [173] XML_3.99-0.19 evaluate_1.0.5
## [175] BiocManager_1.30.26 httpuv_1.6.16
## [177] polyclip_1.10-7 RANN_2.6.2
## [179] tidyr_1.3.1 openssl_2.3.3
## [181] purrr_1.1.0 scattermore_1.2
## [183] future_1.67.0 reshape_0.8.10
## [185] ggplot2_4.0.0 xtable_1.8-4
## [187] restfulr_0.0.16 RSpectra_0.16-2
## [189] later_1.4.4 viridisLite_0.4.2
## [191] tibble_3.3.0 memoise_2.0.1
## [193] AnnotationDbi_1.71.1 GenomicAlignments_1.45.3
## [195] IRanges_2.43.0 cluster_2.1.8.1
## [197] globals_0.18.0 Rhtslib_3.5.0