0.1 Introduction

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.

0.2 Installation

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.

0.3 Loading libraries

library(MuData)
library(SingleCellExperiment)
library(MultiAssayExperiment)
library(SingleCellMultiModal)
library(scater)

library(rhdf5)

0.4 Loading data

We will use CITE-seq data accessible with the SingleCellMultiModal Bioconductor package, which was originally described in Stoeckius et al., 2017.

mae <- CITEseq(
    DataType="cord_blood", modes="*", dry.run=FALSE, version="1.0.0"
)
#> Dataset: cord_blood
#> Working on: scADT_Counts
#> Working on: scRNAseq_Counts
#> Working on: coldata_scRNAseq
#> Working on: scADT_clrCounts
#> see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
#> loading from cache
#> see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
#> loading from cache
#> see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
#> loading from cache
#> see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
#> loading from cache
#> Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
#>   potential for errors with mixed data types
#> harmonizing input:
#>   removing 2277 sampleMap rows with 'primary' not in colData

mae
#> A MultiAssayExperiment object of 3 listed
#>  experiments with user-defined names and respective classes.
#>  Containing an ExperimentList class object of length 3:
#>  [1] scADT: matrix with 13 rows and 7858 columns
#>  [2] scADT_clr: matrix with 13 rows and 7858 columns
#>  [3] scRNAseq: matrix with 36280 rows and 7858 columns
#> Functionality:
#>  experiments() - obtain the ExperimentList instance
#>  colData() - the primary/phenotype DataFrame
#>  sampleMap() - the sample coordination DataFrame
#>  `$`, `[`, `[[` - extract colData columns, subset, or experiment
#>  *Format() - convert into a long or wide DataFrame
#>  assays() - convert ExperimentList to a SimpleList of matrices
#>  exportClass() - save data to flat files

We see two modalities in the object — scRNAseq and scADT, the latter providing counts for antibody-derived tags. Notably, each experiment is a matrix.

0.5 Processing ADT data

While CITE-seq analysis workflows such as CiteFuse should be consulted for more details, below we exemplify simple data transformation in order to demonstrate how their output can be saved to an H5MU file later on.

For ADT counts, we will apply CLR transformation following Hao et al., 2020:

# Define CLR transformation as in the Seurat workflow
clr <- function(data) t(
  apply(data, 1, function(x) log1p(
    x / (exp(sum(log1p(x[x > 0]), na.rm = TRUE) / length(x)))
  ))
)

We will make the ADT modality a SingleCellExperiment object and add an assay with CLR-transformed counts:

adt_counts <- mae[["scADT"]]

mae[["scADT"]] <- SingleCellExperiment(adt_counts)
assay(mae[["scADT"]], "clr") <- clr(adt_counts)

We will also generate reduced dimensions taking advantage of the functionality in the scater package:

mae[["scADT"]] <- runPCA(
  mae[["scADT"]], exprs_values = "clr", ncomponents = 20
)
#> Warning in check_numbers(k = k, nu = nu, nv = nv, limit = min(dim(x)) - : more
#> singular values/vectors requested than available
#> Warning in (function (A, nv = 5, nu = nv, maxit = 1000, work = nv + 7, reorth =
#> TRUE, : You're computing too large a percentage of total singular values, use a
#> standard svd instead.
plotReducedDim(mae[["scADT"]], dimred = "PCA",
               by_exprs_values = "clr", colour_by = "CD3")

plotReducedDim(mae[["scADT"]], dimred = "PCA",
               by_exprs_values = "clr", colour_by = "CD14")

0.6 Writing H5MU files

We can write the contents of the MultiAssayExperiment object into an H5MU file:

writeH5MU(mae, "cord_blood_citeseq.h5mu")

We can check that both modalities were written to the file, whether it was a matrix for RNA or SingleCellExperiment for ADT:

h5 <- rhdf5::H5Fopen("cord_blood_citeseq.h5mu")
h5ls(H5Gopen(h5, "mod"), recursive = FALSE)
#>   group      name     otype dclass dim
#> 0     /     scADT H5I_GROUP           
#> 1     / scADT_clr H5I_GROUP           
#> 2     /  scRNAseq H5I_GROUP

… both assays for ADT — raw counts are stored in X and CLR-transformed counts are in the corresponding layer:

h5ls(H5Gopen(h5, "mod/scADT"), recursive = FALSE)
#>   group   name       otype  dclass       dim
#> 0     /      X H5I_DATASET INTEGER 13 x 7858
#> 1     / layers   H5I_GROUP                  
#> 2     /    obs   H5I_GROUP                  
#> 3     /   obsm   H5I_GROUP                  
#> 4     /    var   H5I_GROUP
h5ls(H5Gopen(h5, "mod/scADT/layers"), recursive = FALSE)
#>   group name       otype dclass       dim
#> 0     /  clr H5I_DATASET  FLOAT 13 x 7858

… as well as reduced dimensions (PCA):

h5ls(H5Gopen(h5, "mod/scADT/obsm"), recursive = FALSE)
#>   group name       otype dclass       dim
#> 0     /  PCA H5I_DATASET  FLOAT 12 x 7858
# There is an alternative way to access groups:
# h5&'mod'&'scADT'&'obsm'
rhdf5::H5close()

0.7 References

  • Muon: multimodal omics analysis framework preprint

  • mudata (Python) documentation

  • muon documentation and web page

  • Stoeckius, M., Hafemeister, C., Stephenson, W., Houck-Loomis, B., Chattopadhyay, P.K., Swerdlow, H., Satija, R. and Smibert, P., 2017. Simultaneous epitope and transcriptome measurement in single cells. Nature methods, 14(9), pp.865-868.

  • Hao, Y., Hao, S., Andersen-Nissen, E., Mauck III, W.M., Zheng, S., Butler, A., Lee, M.J., Wilk, A.J., Darby, C., Zager, M. and Hoffman, P., 2021. Integrated analysis of multimodal single-cell data. Cell.

0.8 Session Info

sessionInfo()
#> R version 4.3.1 (2023-06-16)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.18-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] SingleCellMultiModal_1.13.15 scater_1.30.0               
#>  [3] ggplot2_3.4.4                scuttle_1.12.0              
#>  [5] CiteFuse_1.14.0              MultiAssayExperiment_1.28.0 
#>  [7] SingleCellExperiment_1.24.0  SummarizedExperiment_1.32.0 
#>  [9] Biobase_2.62.0               GenomicRanges_1.54.0        
#> [11] GenomeInfoDb_1.38.0          IRanges_2.36.0              
#> [13] MatrixGenerics_1.14.0        matrixStats_1.0.0           
#> [15] MuData_1.6.0                 rhdf5_2.46.0                
#> [17] S4Vectors_0.40.0             BiocGenerics_0.48.0         
#> [19] Matrix_1.6-1.1               BiocStyle_2.30.0            
#> 
#> loaded via a namespace (and not attached):
#>   [1] later_1.3.1                   splines_4.3.1                
#>   [3] bitops_1.0-7                  filelock_1.0.2               
#>   [5] tibble_3.2.1                  polyclip_1.10-6              
#>   [7] lifecycle_1.0.3               edgeR_4.0.0                  
#>   [9] lattice_0.22-5                MASS_7.3-60                  
#>  [11] magrittr_2.0.3                limma_3.58.0                 
#>  [13] plotly_4.10.3                 sass_0.4.7                   
#>  [15] rmarkdown_2.25                jquerylib_0.1.4              
#>  [17] yaml_2.3.7                    metapod_1.10.0               
#>  [19] httpuv_1.6.12                 cowplot_1.1.1                
#>  [21] bayesm_3.1-6                  DBI_1.1.3                    
#>  [23] RColorBrewer_1.1-3            abind_1.4-5                  
#>  [25] zlibbioc_1.48.0               Rtsne_0.16                   
#>  [27] purrr_1.0.2                   mixtools_2.0.0               
#>  [29] ggraph_2.1.0                  RCurl_1.98-1.12              
#>  [31] tensorA_0.36.2                tweenr_2.0.2                 
#>  [33] rappdirs_0.3.3                GenomeInfoDbData_1.2.11      
#>  [35] ggrepel_0.9.4                 irlba_2.3.5.1                
#>  [37] pheatmap_1.0.12               dqrng_0.3.1                  
#>  [39] DelayedMatrixStats_1.24.0     codetools_0.2-19             
#>  [41] DelayedArray_0.28.0           ggforce_0.4.1                
#>  [43] tidyselect_1.2.0              farver_2.1.1                 
#>  [45] ScaledMatrix_1.10.0           viridis_0.6.4                
#>  [47] BiocFileCache_2.10.0          jsonlite_1.8.7               
#>  [49] BiocNeighbors_1.20.0          ellipsis_0.3.2               
#>  [51] tidygraph_1.2.3               ggridges_0.5.4               
#>  [53] survival_3.5-7                dbscan_1.1-11                
#>  [55] segmented_1.6-4               tools_4.3.1                  
#>  [57] Rcpp_1.0.11                   glue_1.6.2                   
#>  [59] BiocBaseUtils_1.4.0           gridExtra_2.3                
#>  [61] SparseArray_1.2.0             xfun_0.40                    
#>  [63] dplyr_1.1.3                   withr_2.5.1                  
#>  [65] BiocManager_1.30.22           fastmap_1.1.1                
#>  [67] rhdf5filters_1.14.0           bluster_1.12.0               
#>  [69] fansi_1.0.5                   digest_0.6.33                
#>  [71] rsvd_1.0.5                    mime_0.12                    
#>  [73] R6_2.5.1                      colorspace_2.1-0             
#>  [75] RSQLite_2.3.1                 utf8_1.2.4                   
#>  [77] tidyr_1.3.0                   generics_0.1.3               
#>  [79] data.table_1.14.8             robustbase_0.99-0            
#>  [81] graphlayouts_1.0.1            httr_1.4.7                   
#>  [83] htmlwidgets_1.6.2             S4Arrays_1.2.0               
#>  [85] uwot_0.1.16                   pkgconfig_2.0.3              
#>  [87] gtable_0.3.4                  blob_1.2.4                   
#>  [89] XVector_0.42.0                htmltools_0.5.6.1            
#>  [91] bookdown_0.36                 scales_1.2.1                 
#>  [93] png_0.1-8                     SpatialExperiment_1.12.0     
#>  [95] scran_1.30.0                  knitr_1.44                   
#>  [97] rjson_0.2.21                  reshape2_1.4.4               
#>  [99] nlme_3.1-163                  curl_5.1.0                   
#> [101] cachem_1.0.8                  stringr_1.5.0                
#> [103] BiocVersion_3.18.0            parallel_4.3.1               
#> [105] vipor_0.4.5                   AnnotationDbi_1.64.0         
#> [107] pillar_1.9.0                  grid_4.3.1                   
#> [109] vctrs_0.6.4                   promises_1.2.1               
#> [111] randomForest_4.7-1.1          BiocSingular_1.18.0          
#> [113] dbplyr_2.3.4                  beachmat_2.18.0              
#> [115] xtable_1.8-4                  cluster_2.1.4                
#> [117] beeswarm_0.4.0                evaluate_0.22                
#> [119] magick_2.8.1                  cli_3.6.1                    
#> [121] locfit_1.5-9.8                compiler_4.3.1               
#> [123] rlang_1.1.1                   crayon_1.5.2                 
#> [125] labeling_0.4.3                plyr_1.8.9                   
#> [127] ggbeeswarm_0.7.2              stringi_1.7.12               
#> [129] viridisLite_0.4.2             BiocParallel_1.36.0          
#> [131] Biostrings_2.70.0             munsell_0.5.0                
#> [133] lazyeval_0.2.2                compositions_2.0-6           
#> [135] ExperimentHub_2.10.0          sparseMatrixStats_1.14.0     
#> [137] bit64_4.0.5                   Rhdf5lib_1.24.0              
#> [139] KEGGREST_1.42.0               statmod_1.5.0                
#> [141] shiny_1.7.5.1                 interactiveDisplayBase_1.40.0
#> [143] AnnotationHub_3.10.0          kernlab_0.9-32               
#> [145] igraph_1.5.1                  memoise_2.0.1                
#> [147] bslib_0.5.1                   DEoptimR_1.1-3               
#> [149] bit_4.0.5