1 Installation

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("SingleCellMultiModal")

1.1 Load packages

library(MultiAssayExperiment)
library(SpatialExperiment)
library(SingleCellMultiModal)

2 seq-FISH dataset

The dataset consists of two data types, seq-FISH data was provided by Zhu et al. (2018), while scRNA-seq data was provided by Tasic et al. (2016).

Data have been retrievedas part of the Hackathon in the Mathematical Frameworks for Integrative Analysis of Emerging Biological DataTypes workshop.

2.1 Downloading datasets

The user can see the available dataset by using the default options

seqFISH(
    DataType="mouse_visual_cortex", modes="*", dry.run=TRUE, version="2.0.0"
)
##    ah_id                mode file_size rdataclass rdatadateadded
## 1 EH3785        scRNA_Counts    0.2 Mb     matrix     2020-09-14
## 2 EH3786        scRNA_Labels      0 Mb data.frame     2020-09-14
## 3 EH3787 seqFISH_Coordinates      0 Mb data.frame     2020-09-14
## 4 EH3788      seqFISH_Counts    0.2 Mb     matrix     2020-09-14
## 5 EH3789      seqFISH_Labels      0 Mb data.frame     2020-09-14
##   rdatadateremoved
## 1             <NA>
## 2             <NA>
## 3             <NA>
## 4             <NA>
## 5             <NA>

Or simply by running:

seqfish <- seqFISH(
    DataType="mouse_visual_cortex", modes="*", dry.run=FALSE, version="2.0.0"
)
## Working on: scRNA_Counts
## Working on: scRNA_Labels
## Working on: seqFISH_Coordinates
## Working on: seqFISH_Counts
## Working on: seqFISH_Labels
## 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
## see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
## loading from cache
## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types

## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types

## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types
## Note: spatialData and spatialDataNames have been deprecated; all columns should be stored in colData and spatialCoords
seqfish
## A MultiAssayExperiment object of 2 listed
##  experiments with user-defined names and respective classes.
##  Containing an ExperimentList class object of length 2:
##  [1] seqFISH: SpatialExperiment with 113 rows and 1597 columns
##  [2] scRNAseq: SingleCellExperiment with 113 rows and 1722 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

Extract the list of experiments without the associated colData.

experiments(seqfish)
## ExperimentList class object of length 2:
##  [1] seqFISH: SpatialExperiment with 113 rows and 1597 columns
##  [2] scRNAseq: SingleCellExperiment with 113 rows and 1722 columns

2.2 Exploring the data structure

Check row annotations for all experiments:

rownames(seqfish)
## CharacterList of length 2
## [["seqFISH"]] abca15 abca9 acta2 adcy4 aldh3b2 ... wrn zfp182 zfp715 zfp90
## [["scRNAseq"]] abca15 abca9 acta2 adcy4 aldh3b2 ... wrn zfp182 zfp715 zfp90

Take a peek at the sampleMap (graph representation of assays, cells, and barcodes):

sampleMap(seqfish)
## DataFrame with 3319 rows and 3 columns
##         assay     primary     colname
##      <factor> <character> <character>
## 1     seqFISH          V2          V2
## 2     seqFISH          V3          V3
## 3     seqFISH          V4          V4
## 4     seqFISH          V5          V5
## 5     seqFISH          V6          V6
## ...       ...         ...         ...
## 3315 scRNAseq       V1719       V1719
## 3316 scRNAseq       V1720       V1720
## 3317 scRNAseq       V1721       V1721
## 3318 scRNAseq       V1722       V1722
## 3319 scRNAseq       V1723       V1723

2.3 Visualize matching cell identifiers across assays

upsetSamples(seqfish)

This shows that about 1597 cells match across both modalities / assays.

2.4 scRNA-seq data

The scRNA-seq data are accessible with $scRNAseq, which returns a SingleCellExperiment class object, with all its associated methods.

seqfish[["scRNAseq"]]
## class: SingleCellExperiment 
## dim: 113 1722 
## metadata(0):
## assays(1): counts
## rownames(113): abca15 abca9 ... zfp715 zfp90
## rowData names(0):
## colnames(1722): V2 V3 ... V1722 V1723
## colData names(3): broad_type sample_name dissection
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):

Otherwhise the assay function can be used to access the scRNAseq assay stored in the seqfish MultiAssayExperiment object.

head(assay(seqfish, "scRNAseq"))[,1:4]
##         V2 V3  V4 V5
## abca15  11 42  17 42
## abca9   22 46  22 46
## acta2   15 47  15 42
## adcy4   12 45  12 45
## aldh3b2 27 49  27 49
## amigo2  23 43 101 43

2.5 seq-FISH data

The seq-FISH data are accessible with $seqFISH, which returns a SpatialExperiment class object.

seqfish[["seqFISH"]]
## class: SpatialExperiment 
## dim: 113 1597 
## metadata(0):
## assays(1): counts
## rownames(113): abca15 abca9 ... zfp715 zfp90
## rowData names(1): X
## colnames(1597): V2 V3 ... V1597 V1598
## colData names(7): Cell_ID cluster ... Prob sample_id
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
## spatialCoords names(2) : x y
## imgData names(0):

Otherwhise the assay function can be used to access the seqFISH assay stored in the seqfish MultiAssayExperiment object.

head(assay(seqfish, "seqFISH"))[,1:4]
##          V2 V3 V4 V5
## abca15   68 49 50 39
## abca9    41 42 38 36
## acta2    25 23 16 21
## adcy4    39 54 37 18
## aldh3b2 101 47 41 52
## amigo2   93 64 93 93

Spatial data can be retrieved with spatialData function on the SpatialExperiment object.

(sd <- spatialData(seqfish[["seqFISH"]]))
## Note: spatialData and spatialDataNames have been deprecated; all columns should be stored in colData and spatialCoords
## DataFrame with 1597 rows and 2 columns
##         Cell_ID Irrelevant
##       <integer>  <integer>
## V2            1        100
## V3            2        100
## V4            3        100
## V5            4        100
## V6            5        100
## ...         ...        ...
## V1594      1593        100
## V1595      1594        100
## V1596      1595        100
## V1597      1596        100
## V1598      1597        100

Spatial coordinates within the spatial data can be retrieved in matrix form with spatialCoords function on the SpatialExperiment object.

head(sc <- spatialCoords(seqfish[["seqFISH"]]))
##           x       y
## [1,] 265.76 -231.14
## [2,] 290.48 -261.52
## [3,] 257.12 -133.35
## [4,] 753.46 -261.14
## [5,] 700.01 -169.05
## [6,] 415.63 -252.45

Direct access to the colnames of the spacial coordinates with spatialCoordsNames function.

spatialCoordsNames(seqfish[["seqFISH"]])
## [1] "x" "y"

2.6 Other data version

The provided seqFISH dataset comes out in two different versions:

  • 1.0.0 - provides the same seqFISH data as shown in the rest of this vignette, but it returns the full normalized scRNA-seq data matrix (with labels), as released from the original authors on the GEO database.
  • 2.0.0 - provides the same seqFISH data as shown in the rest of this vignette, but it returns a processed subset of the original scRNA-seq data, providing only the same genes present in the seqFISH data matrix.

2.6.1 Data version 1.0.0

The full scRNA-seq data matrix is 24057 rows x 1809 columns.

To access the v1.0.0 simply run

seqFISH(
    DataType="mouse_visual_cortex", modes="*", dry.run=FALSE, version="1.0.0"
)
## Working on: scRNA_Full_Counts
## Working on: scRNA_Full_Labels
## Working on: seqFISH_Coordinates
## Working on: seqFISH_Counts
## Working on: seqFISH_Labels
## 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
## see ?SingleCellMultiModal and browseVignettes('SingleCellMultiModal') for documentation
## loading from cache
## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types

## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types

## Warning: 'ExperimentList' contains 'data.frame' or 'DataFrame',
##   potential for errors with mixed data types
## Note: spatialData and spatialDataNames have been deprecated; all columns should be stored in colData and spatialCoords
## A MultiAssayExperiment object of 2 listed
##  experiments with user-defined names and respective classes.
##  Containing an ExperimentList class object of length 2:
##  [1] seqFISH: SpatialExperiment with 113 rows and 1597 columns
##  [2] scRNAseq: SingleCellExperiment with 24057 rows and 1809 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

3 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] SpatialExperiment_1.12.0    scater_1.30.0              
##  [3] ggplot2_3.4.4               scran_1.30.0               
##  [5] scuttle_1.12.0              rhdf5_2.46.0               
##  [7] HDF5Array_1.30.0            RaggedExperiment_1.26.0    
##  [9] Matrix_1.6-1.1              SingleCellExperiment_1.24.0
## [11] SingleCellMultiModal_1.14.0 MultiAssayExperiment_1.28.0
## [13] SummarizedExperiment_1.32.0 Biobase_2.62.0             
## [15] GenomicRanges_1.54.0        GenomeInfoDb_1.38.0        
## [17] IRanges_2.36.0              S4Vectors_0.40.0           
## [19] BiocGenerics_0.48.0         MatrixGenerics_1.14.0      
## [21] matrixStats_1.0.0           BiocStyle_2.30.0           
## 
## loaded via a namespace (and not attached):
##   [1] jsonlite_1.8.7                magrittr_2.0.3               
##   [3] ggbeeswarm_0.7.2              magick_2.8.1                 
##   [5] farver_2.1.1                  rmarkdown_2.25               
##   [7] zlibbioc_1.48.0               vctrs_0.6.4                  
##   [9] memoise_2.0.1                 DelayedMatrixStats_1.24.0    
##  [11] RCurl_1.98-1.12               htmltools_0.5.6.1            
##  [13] S4Arrays_1.2.0                BiocBaseUtils_1.4.0          
##  [15] AnnotationHub_3.10.0          curl_5.1.0                   
##  [17] BiocNeighbors_1.20.0          Rhdf5lib_1.24.0              
##  [19] SparseArray_1.2.0             sass_0.4.7                   
##  [21] bslib_0.5.1                   plyr_1.8.9                   
##  [23] cachem_1.0.8                  igraph_1.5.1                 
##  [25] mime_0.12                     lifecycle_1.0.3              
##  [27] pkgconfig_2.0.3               rsvd_1.0.5                   
##  [29] R6_2.5.1                      fastmap_1.1.1                
##  [31] GenomeInfoDbData_1.2.11       shiny_1.7.5.1                
##  [33] digest_0.6.33                 colorspace_2.1-0             
##  [35] AnnotationDbi_1.64.0          dqrng_0.3.1                  
##  [37] irlba_2.3.5.1                 ExperimentHub_2.10.0         
##  [39] RSQLite_2.3.1                 beachmat_2.18.0              
##  [41] labeling_0.4.3                filelock_1.0.2               
##  [43] fansi_1.0.5                   httr_1.4.7                   
##  [45] abind_1.4-5                   compiler_4.3.1               
##  [47] bit64_4.0.5                   withr_2.5.1                  
##  [49] BiocParallel_1.36.0           viridis_0.6.4                
##  [51] DBI_1.1.3                     UpSetR_1.4.0                 
##  [53] rappdirs_0.3.3                DelayedArray_0.28.0          
##  [55] rjson_0.2.21                  bluster_1.12.0               
##  [57] tools_4.3.1                   vipor_0.4.5                  
##  [59] beeswarm_0.4.0                interactiveDisplayBase_1.40.0
##  [61] httpuv_1.6.12                 glue_1.6.2                   
##  [63] rhdf5filters_1.14.0           promises_1.2.1               
##  [65] grid_4.3.1                    cluster_2.1.4                
##  [67] generics_0.1.3                gtable_0.3.4                 
##  [69] BiocSingular_1.18.0           ScaledMatrix_1.10.0          
##  [71] metapod_1.10.0                utf8_1.2.4                   
##  [73] XVector_0.42.0                RcppAnnoy_0.0.21             
##  [75] ggrepel_0.9.4                 BiocVersion_3.18.0           
##  [77] pillar_1.9.0                  limma_3.58.0                 
##  [79] later_1.3.1                   dplyr_1.1.3                  
##  [81] BiocFileCache_2.10.0          lattice_0.22-5               
##  [83] bit_4.0.5                     tidyselect_1.2.0             
##  [85] locfit_1.5-9.8                Biostrings_2.70.1            
##  [87] knitr_1.44                    gridExtra_2.3                
##  [89] bookdown_0.36                 edgeR_4.0.0                  
##  [91] xfun_0.40                     statmod_1.5.0                
##  [93] yaml_2.3.7                    evaluate_0.22                
##  [95] codetools_0.2-19              tibble_3.2.1                 
##  [97] BiocManager_1.30.22           cli_3.6.1                    
##  [99] uwot_0.1.16                   xtable_1.8-4                 
## [101] munsell_0.5.0                 jquerylib_0.1.4              
## [103] Rcpp_1.0.11                   dbplyr_2.3.4                 
## [105] png_0.1-8                     parallel_4.3.1               
## [107] ellipsis_0.3.2                blob_1.2.4                   
## [109] sparseMatrixStats_1.14.0      bitops_1.0-7                 
## [111] viridisLite_0.4.2             scales_1.2.1                 
## [113] purrr_1.0.2                   crayon_1.5.2                 
## [115] rlang_1.1.1                   cowplot_1.1.1                
## [117] KEGGREST_1.42.0               formatR_1.14

Tasic, Bosiljka, Vilas Menon, Thuc Nghi Nguyen, Tae Kyung Kim, Tim Jarsky, Zizhen Yao, Boaz Levi, et al. 2016. “Adult Mouse Cortical Cell Taxonomy Revealed by Single Cell Transcriptomics.” Nature Neuroscience 19 (2): 335.

Zhu, Qian, Sheel Shah, Ruben Dries, Long Cai, and Guo-Cheng Yuan. 2018. “Identification of Spatially Associated Subpopulations by Combining scRNAseq and Sequential Fluorescence in Situ Hybridization Data.” Nature Biotechnology 36 (12): 1183.