CITE-seq data provide RNA and surface protein counts for the same cells. This tutorial shows how MuData can be integrated into with Bioconductor workflows to analyse CITE-seq data.
The most recent dev build can be installed from GitHub:
library(remotes)
remotes::install_github("ilia-kats/MuData")
Stable version of MuData
will be available in future bioconductor versions.
library(MuData)
library(SingleCellExperiment)
library(MultiAssayExperiment)
library(CiteFuse)
library(scater)
library(rhdf5)
We will use CITE-seq data available within
CiteFuse
Bioconductor package.
data("CITEseq_example", package = "CiteFuse")
lapply(CITEseq_example, dim)
#> $RNA
#> [1] 19521 500
#>
#> $ADT
#> [1] 49 500
#>
#> $HTO
#> [1] 4 500
This dataset contains three matrices — one with RNA
counts, one with
antibody-derived tags (ADT
) counts and one with hashtag oligonucleotide
(HTO
) counts.
While CITE-seq analysis workflows such as CiteFuse should be consulted for more details, below we exemplify simple data transformations in order to demonstrate how their output can be saved to an H5MU file later on.
Following the CiteFuse tutorial, we start with creating a SingleCellExperiment
object with the three matrices:
sce_citeseq <- preprocessing(CITEseq_example)
sce_citeseq
#> class: SingleCellExperiment
#> dim: 19521 500
#> metadata(0):
#> assays(1): counts
#> rownames(19521): hg19_AL627309.1 hg19_AL669831.5 ... hg19_MT-ND6
#> hg19_MT-CYB
#> rowData names(0):
#> colnames(500): AAGCCGCGTTGTCTTT GATCGCGGTTATCGGT ... TTGGCAACACTAGTAC
#> GCTGCGAGTTGTGGCC
#> colData names(0):
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(2): ADT HTO
We will add a new assay with normalised RNA counts:
sce_citeseq <- scater::logNormCounts(sce_citeseq)
sce_citeseq # new assay: logcounts
#> class: SingleCellExperiment
#> dim: 19521 500
#> metadata(0):
#> assays(2): counts logcounts
#> rownames(19521): hg19_AL627309.1 hg19_AL669831.5 ... hg19_MT-ND6
#> hg19_MT-CYB
#> rowData names(0):
#> colnames(500): AAGCCGCGTTGTCTTT GATCGCGGTTATCGGT ... TTGGCAACACTAGTAC
#> GCTGCGAGTTGTGGCC
#> colData names(1): sizeFactor
#> reducedDimNames(0):
#> mainExpName: NULL
#> altExpNames(2): ADT HTO
To the ADT modality, we will add an assay with normalised counts:
sce_citeseq <- CiteFuse::normaliseExprs(
sce_citeseq, altExp_name = "ADT", transform = "log"
)
altExp(sce_citeseq, "ADT") # new assay: logcounts
#> class: SummarizedExperiment
#> dim: 49 500
#> metadata(0):
#> assays(2): counts logcounts
#> rownames(49): B220 (CD45R) B7-H1 (PD-L1) ... TCRb TCRg
#> rowData names(0):
#> colnames(500): AAGCCGCGTTGTCTTT GATCGCGGTTATCGGT ... TTGGCAACACTAGTAC
#> GCTGCGAGTTGTGGCC
#> colData names(0):
We will also generate reduced dimensions:
sce_citeseq <- scater::runPCA(
sce_citeseq, exprs_values = "logcounts", ncomponents = 20
)
scater::plotReducedDim(sce_citeseq, dimred = "PCA",
by_exprs_values = "logcounts", colour_by = "CD27")
An appropriate structure for multimodal datasets is
MultiAssayExperiment
.
We will make a respective MultiAssayExperiment object from sce_citeseq
:
experiments <- list(
ADT = altExp(sce_citeseq, "ADT"),
HTO = altExp(sce_citeseq, "HTO")
)
# Drop other modalities from sce_citeseq
altExp(sce_citeseq) <- NULL
experiments[["RNA"]] <- sce_citeseq
mae <- MultiAssayExperiment(experiments)
We can write the contents of the MultiAssayExperiment object into an H5MU file:
writeH5MU(mae, "citefuse_example.h5mu")
We can check that all the modalities were written to the file:
h5 <- rhdf5::H5Fopen("citefuse_example.h5mu")
h5ls(H5Gopen(h5, "mod"), recursive = FALSE)
#> group name otype dclass dim
#> 0 / ADT H5I_GROUP
#> 1 / HTO H5I_GROUP
#> 2 / RNA H5I_GROUP
… both assays for ADT — raw counts are stored in X
and normalised counts
are in the corresponding layer:
h5ls(H5Gopen(h5, "mod/ADT"), FALSE)
#> group name otype dclass dim
#> 0 / X H5I_GROUP
#> 1 / layers H5I_GROUP
#> 2 / obs H5I_GROUP
#> 3 / var H5I_GROUP
h5ls(H5Gopen(h5, "mod/ADT/layers"), FALSE)
#> group name otype dclass dim
#> 0 / logcounts H5I_DATASET FLOAT 49 x 500
… as well as reduced dimensions (PCA):
h5ls(H5Gopen(h5, "mod/RNA/obsm"), FALSE)
#> group name otype dclass dim
#> 0 / PCA H5I_DATASET FLOAT 20 x 500
# There is an alternative way to access groups:
# h5&'mod'&'RNA'&'obsm'
rhdf5::H5close()
mudata (Python) documentation
muon documentation and web page
Kim HJ, Lin Y, Geddes TA, Yang P, Yang JYH (2020). “CiteFuse enables multi-modal analysis of CITE-seq data.” Bioinformatics, 36(14), 4137–4143. https://doi.org/10.1093/bioinformatics/btaa282.
Ramos M, Schiffer L, Re A, Azhar R, Basunia A, Cabrera CR, Chan T, Chapman P, Davis S, Gomez-Cabrero D, Culhane AC, Haibe-Kains B, Hansen K, Kodali H, Louis MS, Mer AS, Reister M, Morgan M, Carey V, Waldron L (2017). “Software For The Integration Of Multi-Omics Experiments In Bioconductor.” Cancer Research, 77(21); e39-42.
sessionInfo()
#> R Under development (unstable) (2024-10-21 r87258)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [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] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] scater_1.35.0 ggplot2_3.5.1
#> [3] scuttle_1.17.0 CiteFuse_1.19.0
#> [5] MultiAssayExperiment_1.33.0 SingleCellExperiment_1.29.0
#> [7] SummarizedExperiment_1.37.0 Biobase_2.67.0
#> [9] GenomicRanges_1.59.0 GenomeInfoDb_1.43.0
#> [11] IRanges_2.41.0 MatrixGenerics_1.19.0
#> [13] matrixStats_1.4.1 MuData_1.11.0
#> [15] rhdf5_2.51.0 S4Vectors_0.45.0
#> [17] BiocGenerics_0.53.0 Matrix_1.7-1
#> [19] BiocStyle_2.35.0
#>
#> loaded via a namespace (and not attached):
#> [1] RColorBrewer_1.1-3 tensorA_0.36.2.1 jsonlite_1.8.9
#> [4] magrittr_2.0.3 magick_2.8.5 ggbeeswarm_0.7.2
#> [7] farver_2.1.2 rmarkdown_2.28 zlibbioc_1.53.0
#> [10] vctrs_0.6.5 memoise_2.0.1 tinytex_0.53
#> [13] htmltools_0.5.8.1 S4Arrays_1.7.0 BiocNeighbors_2.1.0
#> [16] Rhdf5lib_1.29.0 SparseArray_1.7.0 sass_0.4.9
#> [19] bslib_0.8.0 htmlwidgets_1.6.4 plyr_1.8.9
#> [22] plotly_4.10.4 cachem_1.1.0 igraph_2.1.1
#> [25] lifecycle_1.0.4 pkgconfig_2.0.3 rsvd_1.0.5
#> [28] R6_2.5.1 fastmap_1.2.0 GenomeInfoDbData_1.2.13
#> [31] digest_0.6.37 colorspace_2.1-1 dqrng_0.4.1
#> [34] irlba_2.3.5.1 beachmat_2.23.0 labeling_0.4.3
#> [37] randomForest_4.7-1.2 fansi_1.0.6 httr_1.4.7
#> [40] polyclip_1.10-7 abind_1.4-8 compiler_4.5.0
#> [43] withr_3.0.2 BiocParallel_1.41.0 viridis_0.6.5
#> [46] highr_0.11 ggforce_0.4.2 MASS_7.3-61
#> [49] bayesm_3.1-6 DelayedArray_0.33.0 bluster_1.17.0
#> [52] tools_4.5.0 vipor_0.4.7 beeswarm_0.4.0
#> [55] glue_1.8.0 dbscan_1.2-0 nlme_3.1-166
#> [58] rhdf5filters_1.19.0 grid_4.5.0 Rtsne_0.17
#> [61] cluster_2.1.6 reshape2_1.4.4 generics_0.1.3
#> [64] gtable_0.3.6 tidyr_1.3.1 data.table_1.16.2
#> [67] BiocSingular_1.23.0 tidygraph_1.3.1 ScaledMatrix_1.15.0
#> [70] metapod_1.15.0 utf8_1.2.4 XVector_0.47.0
#> [73] ggrepel_0.9.6 pillar_1.9.0 stringr_1.5.1
#> [76] limma_3.63.0 robustbase_0.99-4-1 splines_4.5.0
#> [79] dplyr_1.1.4 tweenr_2.0.3 lattice_0.22-6
#> [82] survival_3.7-0 compositions_2.0-8 tidyselect_1.2.1
#> [85] locfit_1.5-9.10 knitr_1.48 gridExtra_2.3
#> [88] bookdown_0.41 edgeR_4.5.0 xfun_0.48
#> [91] graphlayouts_1.2.0 mixtools_2.0.0 statmod_1.5.0
#> [94] DEoptimR_1.1-3 pheatmap_1.0.12 stringi_1.8.4
#> [97] UCSC.utils_1.3.0 lazyeval_0.2.2 yaml_2.3.10
#> [100] evaluate_1.0.1 codetools_0.2-20 kernlab_0.9-33
#> [103] ggraph_2.2.1 tibble_3.2.1 BiocManager_1.30.25
#> [106] cli_3.6.3 uwot_0.2.2 segmented_2.1-3
#> [109] munsell_0.5.1 jquerylib_0.1.4 Rcpp_1.0.13
#> [112] parallel_4.5.0 scran_1.35.0 viridisLite_0.4.2
#> [115] scales_1.3.0 ggridges_0.5.6 purrr_1.0.2
#> [118] crayon_1.5.3 rlang_1.1.4 cowplot_1.1.3