COTAN 2.2.4
library(COTAN)
library(zeallot)
library(data.table)
library(Rtsne)
library(qpdf)
library(GEOquery)
options(parallelly.fork.enable = TRUE)
This tutorial contains the same functionalities as the first release of the COTAN tutorial but done using the new and updated functions.
Download the data-set for "mouse cortex E17.5"
.
dataDir <- tempdir()
GEO <- "GSM2861514"
fName <- "GSM2861514_E175_Only_Cortical_Cells_DGE.txt.gz"
dataSetFile <- file.path(dataDir, GEO, fName)
if (!file.exists(dataSetFile)) {
getGEOSuppFiles(GEO, makeDirectory = TRUE,
baseDir = dataDir, fetch_files = TRUE,
filter_regex = fName)
sample.dataset <- read.csv(dataSetFile, sep = "\t", row.names = 1L)
}
Define a directory where the output will be stored.
outDir <- tempdir()
# Log-level 2 was chosen to showcase better how the package works
# In normal usage a level of 0 or 1 is more appropriate
setLoggingLevel(2L)
#> Setting new log level to 2
# This file will contain all the logs produced by the package
# as if at the highest logging level
setLoggingFile(file.path(outDir, "vignette_v2.log"))
#> Setting log file to be: /tmp/RtmpAwv7mp/vignette_v2.log
Initialize the COTAN
object with the row count table and
the metadata for the experiment.
cond <- "mouse_cortex_E17.5"
#cond <- "test"
#obj = COTAN(raw = sampled.dataset)
obj <- COTAN(raw = sample.dataset)
obj <- initializeMetaDataset(obj,
GEO = GEO,
sequencingMethod = "Drop_seq",
sampleCondition = cond)
#> Initializing `COTAN` meta-data
logThis(paste0("Condition ", getMetadataElement(obj, datasetTags()[["cond"]])),
logLevel = 1L)
#> Condition mouse_cortex_E17.5
Before we proceed to the analysis, we need to clean the data. The analysis will use a matrix of raw UMI counts as the input. To obtain this matrix, we have to remove any potential cell doublets or multiplets, as well as any low quality or dying cells.
We can check the library size (UMI number) with an empirical cumulative distribution function
ECDPlot(obj, yCut = 700L)
cellSizePlot(obj)
genesSizePlot(obj)
#> Warning: Removed 1 row containing missing values or values outside the scale range
#> (`geom_point()`).
mit <- mitochondrialPercentagePlot(obj, genePrefix = "^Mt")
mit[["plot"]]
During the cleaning, every time we want to remove cells or genes
we can use the dropGenesCells()
function.
To drop cells by cell library size:
cells_to_rem <- getCells(obj)[getCellsSize(obj) > 6000L]
obj <- dropGenesCells(obj, cells = cells_to_rem)
To drop cells by gene number: high genes count might indicate doublets…
cells_to_rem <- getCells(obj)[getNumExpressedGenes(obj) > 3000L]
obj <- dropGenesCells(obj, cells = cells_to_rem)
genesSizePlot(obj)
#> Warning: Removed 1 row containing missing values or values outside the scale range
#> (`geom_point()`).
To drop cells by mitochondrial percentage:
to_rem <- mit[["sizes"]][["mit.percentage"]] > 1.5
cells_to_rem <- rownames(mit[["sizes"]])[to_rem]
obj <- dropGenesCells(obj, cells = cells_to_rem)
mit <- mitochondrialPercentagePlot(obj, genePrefix = "^Mt")
mit[["plot"]]
If we do not want to consider the mitochondrial genes we can remove them before starting the analysis.
genes_to_rem <- getGenes(obj)[grep("^Mt", getGenes(obj))]
cells_to_rem <- getCells(obj)[which(getCellsSize(obj) == 0L)]
obj <- dropGenesCells(obj, genes_to_rem, cells_to_rem)
We want also to log the current status.
logThis(paste("n cells", getNumCells(obj)), logLevel = 1L)
#> n cells 873
The clean()
function estimates all the parameters for the data. Therefore, we have to run it again every time we remove any genes or cells from the data.
n_it <- 1L
obj <- clean(obj)
#> Genes/cells selection done: dropped [4744] genes and [0] cells
#> Working on [12278] genes and [873] cells
c(pcaCellsPlot, pcaCellsData, genesPlot, UDEPlot, nuPlot, zoomedNuPlot) %<-%
cleanPlots(obj)
#> PCA: START
#> PCA: DONE
#> Hierarchical clustering: START
#> Hierarchical clustering: DONE
pcaCellsPlot
genesPlot
We can observe here that the red cells are really enriched in hemoglobin genes so we prefer to remove them. They can be extracted from the pcaCellsData
object and removed.
cells_to_rem <- rownames(pcaCellsData)[pcaCellsData[["groups"]] == "B"]
obj <- dropGenesCells(obj, cells = cells_to_rem)
n_it <- 2L
obj <- clean(obj)
#> Genes/cells selection done: dropped [8] genes and [0] cells
#> Working on [12270] genes and [868] cells
c(pcaCellsPlot, pcaCellsData, genesPlot, UDEPlot, nuPlot, zoomedNuPlot) %<-%
cleanPlots(obj)
#> PCA: START
#> PCA: DONE
#> Hierarchical clustering: START
#> Hierarchical clustering: DONE
pcaCellsPlot
To color the PCA based on nu
(so the cells’ efficiency)
UDEPlot
UDE (color) should not correlate with principal components! This is very important.
The next part is used to remove the cells with efficiency too low.
plot(nuPlot)
We can zoom on the smallest values and, if COTAN detects a clear elbow, we can decide to remove the cells.
plot(zoomedNuPlot)
We also save the defined threshold in the metadata and re-run the estimation
yset <- 0.28
obj <- addElementToMetaDataset(obj, "Threshold low UDE cells:", yset)
cells_to_rem <- getCells(obj)[getNu(obj) < yset]
obj <- dropGenesCells(obj, cells = cells_to_rem)
Repeat the estimation after the cells are removed
n_it <- 3L
obj <- clean(obj)
#> Genes/cells selection done: dropped [11] genes and [0] cells
#> Working on [12259] genes and [859] cells
c(pcaCellsPlot, pcaCellsData, genesPlot, UDEPlot, nuPlot, zoomedNuPlot) %<-%
cleanPlots(obj)
#> PCA: START
#> PCA: DONE
#> Hierarchical clustering: START
#> Hierarchical clustering: DONE
pcaCellsPlot
logThis(paste("n cells", getNumCells(obj)), logLevel = 1L)
#> n cells 859
In this part, all the contingency tables are computed
and used to get the statistics necessary to COEX
evaluation and storing
obj <- proceedToCoex(obj, calcCoex = TRUE, cores = 10L,
saveObj = TRUE, outDir = outDir)
#> Cotan analysis functions started
#> Genes/cells selection done: dropped [0] genes and [0] cells
#> Working on [12259] genes and [859] cells
#> PCA: START
#> PCA: DONE
#> Hierarchical clustering: START
#> Hierarchical clustering: DONE
#> Estimate dispersion: START
#> Estimate dispersion: DONE
#> dispersion | min: -0.05609130859375 | max: 375.5 | % negative: 19.406150583245
#> Cotan genes' coex estimation started
#> Calculate genes' coex: START
#> Warning in asMethod(object): sparse->dense coercion: allocating vector of size
#> 1.1 GiB
#> Calculate genes' coex: DONE
#> Saving elaborated data locally at: /tmp/RtmpAwv7mp/mouse_cortex_E17.5.cotan.RDS
# saving the structure
saveRDS(obj, file = file.path(outDir, paste0(cond, ".cotan.RDS")))
It is also possible to run directly a single function if we don’t want to clean anything.
obj2 <- automaticCOTANObjectCreation(
raw = sample.dataset,
GEO = GEO,
sequencingMethod = "Drop_seq",
sampleCondition = cond,
calcCoex = TRUE, saveObj = TRUE, outDir = outDir, cores = 10L)
To calculate the GDI
we can run:
quantP <- calculateGDI(obj)
#> Calculating S: START
#> Calculating S: DONE
#> Calculate GDI dataframe: START
#> Calculate GDI dataframe: DONE
head(quantP)
#> sum.raw.norm GDI exp.cells
#> 0610007N19Rik 3.330038 1.620300 2.211874
#> 0610007P14Rik 5.307461 1.374303 18.975553
#> 0610009B22Rik 4.730827 1.281160 11.525029
#> 0610009D07Rik 6.311724 1.371413 43.189756
#> 0610009E02Rik 2.587426 1.015026 1.164144
#> 0610009L18Rik 3.506976 1.255881 3.026775
The next function can either plot the GDI
for the dataset directly or
use the pre-computed dataframe.
It marks a 1.5
threshold (in red) and
the two highest quantiles (in blue) by default.
We can also specify some gene sets (three in this case) that
we want to label explicitly in the GDI
plot.
genesList <- list(
"NPGs" = c("Nes", "Vim", "Sox2", "Sox1", "Notch1", "Hes1", "Hes5", "Pax6"),
"PNGs" = c("Map2", "Tubb3", "Neurod1", "Nefm", "Nefl", "Dcx", "Tbr1"),
"hk" = c("Calm1", "Cox6b1", "Ppia", "Rpl18", "Cox7c", "Erh", "H3f3a",
"Taf1", "Taf2", "Gapdh", "Actb", "Golph3", "Zfr", "Sub1",
"Tars", "Amacr")
)
# needs to be recalculated after the changes in nu/dispersion above
#obj <- calculateCoex(obj, actOnCells = FALSE)
GDIPlot(obj, GDIIn = quantP, cond = cond, genes = genesList)
#> GDI plot
#> Removed 0 low GDI genes (such as the fully-expressed) in GDI plot
The percentage of cells expressing the gene in the third column of this data-frame is reported.
To perform the Gene Pair Analysis, we can plot a heatmap of the COEX
values
between two gene sets.
We have to define the different gene sets (list.genes
) in a list.
Then we can choose which sets to use in the function parameter sets
(for example, from 1 to 3).
We also have to provide an array of the file name prefixes for each condition
(for example, “mouse_cortex_E17.5”).
In fact, this function can plot genes relationships across many different
conditions to get a complete overview.
print(cond)
#> [1] "mouse_cortex_E17.5"
heatmapPlot(genesLists = genesList, sets = (1L:3L),
conditions = c(cond), dir = outDir)
#> heatmap plot: START
#> Calculating S: START
#> Calculating S: DONE
#> calculating PValues: START
#> Get p-values on a set of genes on columns and on a set of genes on rows
#> calculating PValues: DONE
#> min coex: -0.456342887197441 max coex: 0.46754881344819
We can also plot a general heatmap of COEX
values based on some markers like
the following one.
genesHeatmapPlot(obj, primaryMarkers = c("Satb2", "Bcl11b", "Vim", "Hes1"),
pValueThreshold = 0.001, symmetric = TRUE)
#> Calculating S: START
#> Calculating S: DONE
#> calculating PValues: START
#> Get p-values genome wide on columns and genome wide on rows
#> calculating PValues: DONE
genesHeatmapPlot(obj, primaryMarkers = c("Satb2", "Bcl11b", "Fezf2"),
secondaryMarkers = c("Gabra3", "Meg3", "Cux1", "Neurod6"),
pValueThreshold = 0.001, symmetric = FALSE)
#> Calculating S: START
#> Calculating S: DONE
#> calculating PValues: START
#> Get p-values genome wide on columns and genome wide on rows
#> calculating PValues: DONE
Sometimes we can also be interested in the numbers present directly in the contingency tables for two specific genes. To get them we can use two functions:
contingencyTables()
to produce the observed and expected data
c(observedCT, expectedCT) %<-% contingencyTables(obj, g1 = "Satb2",
g2 = "Bcl11b")
print("Observed CT")
#> [1] "Observed CT"
observedCT
#> Satb2.yes Satb2.no
#> Bcl11b.yes 47 150
#> Bcl11b.no 290 372
print("Expected CT")
#> [1] "Expected CT"
expectedCT
#> Satb2.yes Satb2.no
#> Bcl11b.yes 83.36639 113.6328
#> Bcl11b.no 253.63326 408.3675
Another useful function is getGenesCoex()
. This can be used to extract
the whole or a partial COEX
matrix from a COTAN
object.
# For the whole matrix
coex <- getGenesCoex(obj, zeroDiagonal = FALSE)
coex[1L:5L, 1L:5L]
#> 5 x 5 Matrix of class "dspMatrix"
#> 0610007N19Rik 0610007P14Rik 0610009B22Rik 0610009D07Rik
#> 0610007N19Rik 0.665772858 -0.021430294 0.034923825 0.01382780
#> 0610007P14Rik -0.021430294 0.923212383 0.009594258 -0.01433096
#> 0610009B22Rik 0.034923825 0.009594258 0.922069336 0.00438200
#> 0610009D07Rik 0.013827801 -0.014330962 0.004382000 0.92255549
#> 0610009E02Rik -0.009300606 -0.057914365 0.052396761 0.02412341
#> 0610009E02Rik
#> 0610007N19Rik -0.009300606
#> 0610007P14Rik -0.057914365
#> 0610009B22Rik 0.052396761
#> 0610009D07Rik 0.024123410
#> 0610009E02Rik 0.369268813
# For a partial matrix
coex <- getGenesCoex(obj, genes = c("Satb2", "Bcl11b", "Fezf2"))
head(coex)
#> 6 x 3 Matrix of class "dgeMatrix"
#> Bcl11b Fezf2 Satb2
#> 0610007N19Rik -0.033460958 -0.001391290 -0.11373396
#> 0610007P14Rik 0.026674964 0.011215580 -0.01377812
#> 0610009B22Rik -0.013929787 0.083759974 -0.03546084
#> 0610009D07Rik 0.007327526 0.001390828 -0.03446913
#> 0610009E02Rik 0.059653952 0.031358876 -0.02901426
#> 0610009L18Rik 0.019133761 0.046478092 -0.01750999
COTAN
provides a way to establish genes’ clusters given some lists of markers
layersGenes <- list(
"L1" = c("Reln", "Lhx5"),
"L2/3" = c("Satb2", "Cux1"),
"L4" = c("Rorb", "Sox5"),
"L5/6" = c("Bcl11b", "Fezf2"),
"Prog" = c("Vim", "Hes1")
)
c(gSpace, eigPlot, pcaClustersDF, treePlot) %<-%
establishGenesClusters(obj, groupMarkers = layersGenes,
numGenesPerMarker = 25L, kCuts = 6L)
#> Establishing gene clusters - START
#> Calculating gene co-expression space - START
#> Calculating S: START
#> Calculating S: DONE
#> calculating PValues: START
#> Get p-values on a set of genes on columns and genome wide on rows
#> calculating PValues: DONE
#> Number of selected secondary markers: 185
#> Calculating S: START
#> Calculating S: DONE
#> Calculating gene co-expression space - DONE
#> Establishing gene clusters - DONE
plot(eigPlot)
plot(treePlot)
UMAPPlot(pcaClustersDF[, 1L:10L],
clusters = pcaClustersDF[["hclust"]],
elements = layersGenes,
title = "Genes' clusters UMAP Plot")
#> UMAP plot
#> Found more than one class "dist" in cache; using the first, from namespace 'BiocGenerics'
#> Also defined by 'spam'
It is possible to obtain a cell clusterization based on the concept of
uniformity of expression of the genes across the cells. That is the cluster
satisfies the null hypothesis of the COTAN
model:
the genes expression is not dependent on the cell in consideration.
fineClusters <- cellsUniformClustering(obj, GDIThreshold = 1.4,
initialResolution = 0.8,
cores = 10L, saveObj = TRUE,
outDir = outDir)[["clusters"]]
obj <- addClusterization(obj, clName = "FineClusters", clusters = fineClusters)
coexDF <- DEAOnClusters(obj, clusters = fineClusters)
obj <- addClusterizationCoex(obj, clName = "FineClusters",
coexDF = coexDF)
c(mergedClusters, coexDF) %<-%
mergeUniformCellsClusters(obj, GDIThreshold = 1.4, cores = 10L,
saveObj = TRUE, outDir = outDir)
obj <- addClusterization(obj, clName = "MergedClusters",
clusters = mergedClusters, coexDF = coexDF)
mergedUMAPPlot <- UMAPPlot(coexDF, elements = layersGenes,
title = "Merged Clusters UMAP Plot")
plot(mergedUMAPPlot)
The next few lines are just to clean.
if (file.exists(file.path(outDir, paste0(cond, ".cotan.RDS")))) {
#Delete file if it exists
file.remove(file.path(outDir, paste0(cond, ".cotan.RDS")))
}
#> [1] TRUE
unlink(file.path(outDir, cond), recursive = TRUE)
file.remove(dataSetFile)
#> [1] TRUE
# stop logging to file
setLoggingFile("")
#> Closing previous log file - Setting log file to be:
file.remove(file.path(outDir, "vignette_v1.log"))
#> Warning in file.remove(file.path(outDir, "vignette_v1.log")): cannot remove
#> file '/tmp/RtmpAwv7mp/vignette_v1.log', reason 'No such file or directory'
#> [1] FALSE
options(parallelly.fork.enable = FALSE)
sessionInfo()
#> R version 4.3.3 (2024-02-29)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.4 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] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] GEOquery_2.70.0 Biobase_2.62.0 BiocGenerics_0.48.1
#> [4] qpdf_1.3.3 Rtsne_0.17 data.table_1.15.4
#> [7] zeallot_0.1.0 COTAN_2.2.4 BiocStyle_2.30.0
#>
#> loaded via a namespace (and not attached):
#> [1] RcppAnnoy_0.0.22 splines_4.3.3
#> [3] later_1.3.2 tibble_3.2.1
#> [5] polyclip_1.10-6 fastDummies_1.7.3
#> [7] lifecycle_1.0.4 doParallel_1.0.17
#> [9] globals_0.16.3 lattice_0.22-6
#> [11] MASS_7.3-60.0.1 dendextend_1.17.1
#> [13] magrittr_2.0.3 limma_3.58.1
#> [15] plotly_4.10.4 sass_0.4.9
#> [17] rmarkdown_2.26 jquerylib_0.1.4
#> [19] yaml_2.3.8 httpuv_1.6.15
#> [21] Seurat_5.0.3 sctransform_0.4.1
#> [23] askpass_1.2.0 spam_2.10-0
#> [25] sp_2.1-3 spatstat.sparse_3.0-3
#> [27] reticulate_1.35.0 cowplot_1.1.3
#> [29] pbapply_1.7-2 RColorBrewer_1.1-3
#> [31] abind_1.4-5 zlibbioc_1.48.2
#> [33] purrr_1.0.2 circlize_0.4.16
#> [35] IRanges_2.36.0 S4Vectors_0.40.2
#> [37] ggrepel_0.9.5 irlba_2.3.5.1
#> [39] listenv_0.9.1 spatstat.utils_3.0-4
#> [41] umap_0.2.10.0 goftest_1.2-3
#> [43] RSpectra_0.16-1 spatstat.random_3.2-3
#> [45] dqrng_0.3.2 fitdistrplus_1.1-11
#> [47] parallelly_1.37.1 DelayedMatrixStats_1.24.0
#> [49] leiden_0.4.3.1 codetools_0.2-20
#> [51] DelayedArray_0.28.0 xml2_1.3.6
#> [53] tidyselect_1.2.1 shape_1.4.6.1
#> [55] farver_2.1.1 ScaledMatrix_1.10.0
#> [57] viridis_0.6.5 matrixStats_1.2.0
#> [59] stats4_4.3.3 spatstat.explore_3.2-7
#> [61] jsonlite_1.8.8 GetoptLong_1.0.5
#> [63] progressr_0.14.0 ggridges_0.5.6
#> [65] survival_3.5-8 iterators_1.0.14
#> [67] foreach_1.5.2 tools_4.3.3
#> [69] ica_1.0-3 Rcpp_1.0.12
#> [71] glue_1.7.0 gridExtra_2.3
#> [73] SparseArray_1.2.4 xfun_0.43
#> [75] MatrixGenerics_1.14.0 ggthemes_5.1.0
#> [77] dplyr_1.1.4 withr_3.0.0
#> [79] BiocManager_1.30.22 fastmap_1.1.1
#> [81] fansi_1.0.6 openssl_2.1.1
#> [83] digest_0.6.35 rsvd_1.0.5
#> [85] parallelDist_0.2.6 R6_2.5.1
#> [87] mime_0.12 colorspace_2.1-0
#> [89] Cairo_1.6-2 scattermore_1.2
#> [91] tensor_1.5 spatstat.data_3.0-4
#> [93] utf8_1.2.4 tidyr_1.3.1
#> [95] generics_0.1.3 httr_1.4.7
#> [97] htmlwidgets_1.6.4 S4Arrays_1.2.1
#> [99] uwot_0.1.16 pkgconfig_2.0.3
#> [101] gtable_0.3.4 ComplexHeatmap_2.18.0
#> [103] lmtest_0.9-40 XVector_0.42.0
#> [105] htmltools_0.5.8.1 dotCall64_1.1-1
#> [107] bookdown_0.38 clue_0.3-65
#> [109] SeuratObject_5.0.1 scales_1.3.0
#> [111] png_0.1-8 knitr_1.45
#> [113] tzdb_0.4.0 reshape2_1.4.4
#> [115] rjson_0.2.21 curl_5.2.1
#> [117] nlme_3.1-164 cachem_1.0.8
#> [119] zoo_1.8-12 GlobalOptions_0.1.2
#> [121] stringr_1.5.1 KernSmooth_2.23-22
#> [123] parallel_4.3.3 miniUI_0.1.1.1
#> [125] RcppZiggurat_0.1.6 pillar_1.9.0
#> [127] grid_4.3.3 vctrs_0.6.5
#> [129] RANN_2.6.1 promises_1.3.0
#> [131] BiocSingular_1.18.0 beachmat_2.18.1
#> [133] xtable_1.8-4 cluster_2.1.6
#> [135] evaluate_0.23 magick_2.8.3
#> [137] readr_2.1.5 cli_3.6.2
#> [139] compiler_4.3.3 rlang_1.1.3
#> [141] crayon_1.5.2 future.apply_1.11.2
#> [143] labeling_0.4.3 plyr_1.8.9
#> [145] stringi_1.8.3 viridisLite_0.4.2
#> [147] deldir_2.0-4 BiocParallel_1.36.0
#> [149] assertthat_0.2.1 munsell_0.5.1
#> [151] lazyeval_0.2.2 spatstat.geom_3.2-9
#> [153] PCAtools_2.14.0 Matrix_1.6-5
#> [155] RcppHNSW_0.6.0 hms_1.1.3
#> [157] patchwork_1.2.0 sparseMatrixStats_1.14.0
#> [159] future_1.33.2 ggplot2_3.5.0
#> [161] statmod_1.5.0 shiny_1.8.1.1
#> [163] highr_0.10 ROCR_1.0-11
#> [165] Rfast_2.1.0 igraph_2.0.3
#> [167] RcppParallel_5.1.7 bslib_0.7.0