SuperCellCyto 1.0.0
This vignette describes the steps to generate supercells for cytometry data using SuperCellCyto R package.
Briefly, supercells are “mini” clusters of cells that are similar in their marker expressions. The motivation behind supercells is that instead of analysing millions of individual cells, you can analyse thousands of supercells, making downstream analysis much faster while maintaining biological interpretability.
See other vignettes for how to:
You can install stable version of SuperCellCyto from Bioconductor using:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("SuperCellCyto")
For the latest development version, you can install it from GitHub using pak:
if (!requireNamespace("pak", quietly = TRUE))
install.packages("pak")
pak::install_github("phipsonlab/SuperCellCyto")
The function which creates supercells is called runSuperCellCyto, and it
operates on a data.table object, an enhanced version of R native
data.frame.
In addition to needing the data stored in a data.table object it also
requires:
runSuperCellCyto does not perform any data transformation or scaling.If you are not sure how to import CSV or FCS files into data.table
object, and/or how to subsequently prepare the object ready for
SuperCellCyto, please consult this vignette.
In that vignette, we also provide an explanation behind why we need to have the
cell ID and sample column.
For this vignette, we will simulate some toy data using the simCytoData function.
Specifically, we will simulate 15 markers and 3 samples,
with each sample containing 10,000 cells.
Hence in total, we will have a toy dataset containing 15 markers and
30,000 cells.
n_markers <- 15
n_samples <- 3
dat <- simCytoData(nmarkers = n_markers, ncells = rep(10000, n_samples))
head(dat)
#> Marker_1 Marker_2 Marker_3 Marker_4 Marker_5 Marker_6 Marker_7 Marker_8
#> <num> <num> <num> <num> <num> <num> <num> <num>
#> 1: 12.19003 4.559720 20.55929 19.51386 18.90056 15.92570 6.534928 5.519177
#> 2: 12.19352 6.697852 19.97353 20.03694 18.59536 16.99655 7.681645 5.437464
#> 3: 13.59126 5.724511 18.17800 20.05466 17.52945 17.94238 7.651896 5.984866
#> 4: 12.92135 5.527831 20.21993 19.96520 18.90484 15.78041 7.726421 6.368659
#> 5: 12.46926 4.352631 19.79060 18.66338 18.27634 17.59175 6.501828 5.132806
#> 6: 12.61249 5.853810 20.09856 18.06383 17.58630 18.84149 7.180943 5.923827
#> Marker_9 Marker_10 Marker_11 Marker_12 Marker_13 Marker_14 Marker_15
#> <num> <num> <num> <num> <num> <num> <num>
#> 1: 9.99194 9.555644 17.64321 13.06032 15.63808 15.31595 10.79193
#> 2: 11.35668 10.581524 18.89986 13.99543 14.73649 15.37497 10.30666
#> 3: 10.09858 9.941482 17.72623 13.55976 14.46394 14.14847 12.42011
#> 4: 10.13852 11.955057 16.76198 14.67172 14.32644 17.04551 10.75823
#> 5: 11.42957 7.583365 19.20930 11.48967 15.26172 14.17766 11.15529
#> 6: 12.58768 12.580274 15.14297 12.52631 14.03529 15.34384 11.54611
#> Sample Cell_Id
#> <char> <char>
#> 1: Sample_1 Cell_1
#> 2: Sample_1 Cell_2
#> 3: Sample_1 Cell_3
#> 4: Sample_1 Cell_4
#> 5: Sample_1 Cell_5
#> 6: Sample_1 Cell_6
For our toy dataset, we will transform our data using arcsinh transformation.
We will use the base R asinh function to do this:
# Specify which columns are the markers to transform
marker_cols <- paste0("Marker_", seq_len(n_markers))
# The co-factor for arc-sinh
cofactor <- 5
# Do the transformation
dat_asinh <- asinh(dat[, marker_cols, with = FALSE] / cofactor)
# Rename the new columns
marker_cols_asinh <- paste0(marker_cols, "_asinh")
names(dat_asinh) <- marker_cols_asinh
# Add them our previously loaded data
dat <- cbind(dat, dat_asinh)
head(dat[, marker_cols_asinh, with = FALSE])
#> Marker_1_asinh Marker_2_asinh Marker_3_asinh Marker_4_asinh Marker_5_asinh
#> <num> <num> <num> <num> <num>
#> 1: 1.623958 0.8177185 2.121491 2.070858 2.039954
#> 2: 1.624223 1.1023489 2.093428 2.096503 2.024225
#> 3: 1.725373 0.9802169 2.002320 2.097361 1.967337
#> 4: 1.678081 0.9540860 2.105326 2.093023 2.040173
#> 5: 1.644945 0.7867981 2.084505 2.027751 2.007523
#> 6: 1.655555 0.9971201 2.099482 1.996247 1.970451
#> Marker_6_asinh Marker_7_asinh Marker_8_asinh Marker_9_asinh Marker_10_asinh
#> <num> <num> <num> <num> <num>
#> 1: 1.875422 1.082703 0.9529245 1.442914 1.403170
#> 2: 1.937685 1.214747 0.9419075 1.558789 1.494471
#> 3: 1.989747 1.211497 1.0140340 1.452418 1.438389
#> 4: 1.866682 1.219623 1.0623328 1.455957 1.605975
#> 5: 1.970749 1.078674 0.9000310 1.564648 1.203976
#> 6: 2.036929 1.158838 1.0061835 1.653724 1.653177
#> Marker_11_asinh Marker_12_asinh Marker_13_asinh Marker_14_asinh
#> <num> <num> <num> <num>
#> 1: 1.973560 1.688065 1.858048 1.838242
#> 2: 2.039918 1.752922 1.801647 1.841899
#> 3: 1.978077 1.723195 1.783986 1.763169
#> 4: 1.924360 1.797477 1.774963 1.940444
#> 5: 2.055627 1.569455 1.834871 1.765113
#> 6: 1.827450 1.649184 1.755600 1.839972
#> Marker_15_asinh
#> <num>
#> 1: 1.512305
#> 2: 1.470732
#> 3: 1.641281
#> 4: 1.509467
#> 5: 1.542437
#> 6: 1.573949
We will also create a column Cell_id_dummy which uniquely identify each cell.
It will have values such as Cell_1, Cell_2, all the way until Cell_x
where x is the number of cells in the dataset.
dat$Cell_id_dummy <- paste0("Cell_", seq_len(nrow(dat)))
head(dat$Cell_id_dummy, n = 10)
#> [1] "Cell_1" "Cell_2" "Cell_3" "Cell_4" "Cell_5" "Cell_6" "Cell_7"
#> [8] "Cell_8" "Cell_9" "Cell_10"
By default, the simCytoData function will generate cells for multiple samples,
and that the resulting data.table object will already have a column
called Sample that denotes the sample the cells come from.
unique(dat$Sample)
#> [1] "Sample_1" "Sample_2" "Sample_3"
Let’s take note of the sample and cell id column for later.
sample_col <- "Sample"
cell_id_col <- "Cell_id_dummy"
Now that we have our data, let’s create some supercells.
To do this, we will use runSuperCellCyto function and pass the markers,
sample and cell ID columns as parameters.
The reason why we need to specify the markers is because the function will
create supercells based on only the expression of those markers.
We highly recommend creating supercells using all markers in your data, let
that be cell type or cell state markers.
However, if for any reason you only want to only use a subset of the markers in
your data, then make sure you specify them in a vector that you later pass to
runSuperCellCyto function.
For this tutorial, we will use all the arcsinh transformed markers in the toy data.
supercells <- runSuperCellCyto(
dt = dat,
markers = marker_cols_asinh,
sample_colname = sample_col,
cell_id_colname = cell_id_col
)
Let’s dig deeper into the object it created:
class(supercells)
#> [1] "list"
It is a list containing 3 elements:
names(supercells)
#> [1] "supercell_expression_matrix" "supercell_cell_map"
#> [3] "supercell_object"
The supercell_object contains the metadata used to create the supercells.
It is a list, and each element contains the metadata used to create the
supercells for a sample.
This will come in handy if we need to either regenerate the supercells using
different gamma values (so we get more or less supercells) or do some
debugging later down the line.
More on regenerating supercells on
Controlling supercells granularity
section below.
The supercell_expression_matrix contains the marker expression of each
supercell.
These are calculated by taking the average of the marker expression of
all the cells contained within a supercell.
head(supercells$supercell_expression_matrix)
#> Marker_1_asinh Marker_2_asinh Marker_3_asinh Marker_4_asinh Marker_5_asinh
#> <num> <num> <num> <num> <num>
#> 1: 1.595206 0.7453904 2.094246 2.067504 1.987042
#> 2: 1.676257 0.7352226 2.085970 2.084532 2.028468
#> 3: 1.648110 0.9267569 2.079317 2.068491 1.993542
#> 4: 1.623450 0.7370731 2.096514 2.055442 1.977871
#> 5: 1.628832 1.1021997 2.103388 2.052827 2.001398
#> 6: 1.585950 0.7882386 2.090255 2.058435 1.996948
#> Marker_6_asinh Marker_7_asinh Marker_8_asinh Marker_9_asinh Marker_10_asinh
#> <num> <num> <num> <num> <num>
#> 1: 1.956772 1.100796 1.0472439 1.525352 1.413901
#> 2: 1.971433 1.195208 0.9904803 1.501270 1.567923
#> 3: 1.952632 1.049286 0.8963540 1.481892 1.384651
#> 4: 1.957064 1.128284 0.9425430 1.594824 1.555330
#> 5: 1.940896 1.196791 0.9779183 1.456972 1.437334
#> 6: 1.941216 1.120207 0.7479158 1.476046 1.482651
#> Marker_11_asinh Marker_12_asinh Marker_13_asinh Marker_14_asinh
#> <num> <num> <num> <num>
#> 1: 1.956247 1.664558 1.735851 1.847143
#> 2: 1.982839 1.659290 1.792077 1.837896
#> 3: 1.955532 1.590391 1.754894 1.832268
#> 4: 1.957976 1.651476 1.714221 1.851620
#> 5: 1.943848 1.639942 1.817014 1.847960
#> 6: 1.942685 1.747602 1.688746 1.873493
#> Marker_15_asinh Sample SuperCellId
#> <num> <char> <char>
#> 1: 1.417910 Sample_1 SuperCell_1_Sample_Sample_1
#> 2: 1.550315 Sample_1 SuperCell_2_Sample_Sample_1
#> 3: 1.475592 Sample_1 SuperCell_3_Sample_Sample_1
#> 4: 1.482331 Sample_1 SuperCell_4_Sample_Sample_1
#> 5: 1.524245 Sample_1 SuperCell_5_Sample_Sample_1
#> 6: 1.543158 Sample_1 SuperCell_6_Sample_Sample_1
Therein, we will have the following columns:
markers_col variable.
In this example, they are the arcsinh transformed markers in our toy data.Sample in this case) denoting which sample a supercell
belongs to, (note the column name is the same as what is stored in sample_col
variable).SuperCellId column denoting the unique ID of the supercell.Let’s have a look at SuperCellId:
head(unique(supercells$supercell_expression_matrix$SuperCellId))
#> [1] "SuperCell_1_Sample_Sample_1" "SuperCell_2_Sample_Sample_1"
#> [3] "SuperCell_3_Sample_Sample_1" "SuperCell_4_Sample_Sample_1"
#> [5] "SuperCell_5_Sample_Sample_1" "SuperCell_6_Sample_Sample_1"
Let’s break down one of them, SuperCell_1_Sample_Sample_1.
SuperCell_1 is a numbering (1 to however many supercells there are in
a sample) used to uniquely identify each supercell in a sample.
Notably, you may encounter this (SuperCell_1, SuperCell_2) being repeated
across different samples, e.g.,
supercell_ids <- unique(supercells$supercell_expression_matrix$SuperCellId)
supercell_ids[grep("SuperCell_1_", supercell_ids)]
#> [1] "SuperCell_1_Sample_Sample_1" "SuperCell_1_Sample_Sample_2"
#> [3] "SuperCell_1_Sample_Sample_3"
While these 3 supercells’ id are pre-fixed with SuperCell_1, it does
not make them equal to one another!
SuperCell_1_Sample_Sample_1 will only contain cells from Sample_1 while
SuperCell_1_Sample_Sample_2 will only contain cells from Sample_2.
By now, you may have noticed that we appended the sample name into each supercell id. This aids in differentiating the supercells in different samples.
supercell_cell_map maps each cell in our dataset to the supercell it
belongs to.
head(supercells$supercell_cell_map)
#> SuperCellID CellId Sample
#> <char> <char> <char>
#> 1: SuperCell_289_Sample_Sample_1 Cell_1 Sample_1
#> 2: SuperCell_31_Sample_Sample_1 Cell_2 Sample_1
#> 3: SuperCell_90_Sample_Sample_1 Cell_3 Sample_1
#> 4: SuperCell_190_Sample_Sample_1 Cell_4 Sample_1
#> 5: SuperCell_196_Sample_Sample_1 Cell_5 Sample_1
#> 6: SuperCell_223_Sample_Sample_1 Cell_6 Sample_1
This map is very useful if we later need to expand the supercells out. Additionally, this is also the reason why we need to have a column in the dataset which uniquely identify each cell.
runSuperCellCyto in parallelBy default, runSuperCellCyto will process each sample one after the other.
As each sample is processed independent of one another, strictly speaking, we
can process all of them in parallel.
To do this, we need to:
BiocParallelParam object from the BiocParallel package.
This object can either be of type MulticoreParamor SnowParam.
We highly recommend consulting their vignette for more information.BiocParallelParam object to the number of
samples we have in the dataset.load_balancing parameter for runSuperCellCyto function to TRUE.
This is to ensure even distribution of the supercell creation jobs.
As each sample will be processed by a parallel job, we don’t want a job that
processs large sample to also be assigned other smaller samples if possible.
If you want to know more how this feature works, please refer to our manuscript.supercell_par <- runSuperCellCyto(
dt = dat,
markers = marker_cols_asinh,
sample_colname = sample_col,
cell_id_colname = cell_id_col,
BPPARAM = MulticoreParam(tasks = n_samples),
load_balancing = TRUE
)
This is described in the runSuperCellCyto function’s documentation, but let’s
briefly go through it here.
The runSuperCellCyto function is equipped with various parameters which
can be customised to alter the composition of the supercells.
The one that is very likely to be used the most is the gamma parameter,
denoted as gam in the function.
By default, the value for gam is set to 20, which we found work well for
most cases.
The gamma parameter controls how many supercells to generate, and
indirectly, how many cells are captured within each supercell.
This parameter is resolved into the following formula
gamma=n_cells/n_supercells where n_cell denotes the number of cells and
n_supercells denotes the number of supercells.
In general, the larger gamma parameter is set to, the less supercells we will get. Say for instance we have 10,000 cells. If gamma is set to 10, we will end up with about 1,000 supercells, whereas if gamma is set to 50, we will end up with about 200 supercells.
You may have noticed, after reading the sections above, runSuperCellCyto
is ran on each sample independent of each other, and that we can only set
1 value as the gamma parameter.
Indeed, for now, the same gamma value will be used across all samples,
and that depending on how many cells we have in each sample, we will end up
with different number of supercells for each sample.
For instance, say we have 10,000 cells for sample 1, and 100,000 cells for
sample 2.
If gamma is set to 10, for sample 1, we will get 1,000 supercells (10,000/10)
while for sample 2, we will get 10,000 supercells (100,000/10).
Do note: whatever gamma value you chose, you should not expect each supercell to contain exactly the same number of cells. This behaviour is intentional to ensure rare cell types are not intermixed with non-rare cell types in a supercell.
If you have run runSuperCellCyto once and have not discarded the
SuperCell object it generated (no serious, please don’t!),
you can use the object to quickly
regenerate supercells using different gamma values.
As an example, using the SuperCell object we have generated for our
toy dataset, we will regenerate the supercells using gamma of 10 and 50.
The function to do this is recomputeSupercells.
We will store the output in a list, one element per gamma value.
addt_gamma_vals <- c(10, 50)
supercells_addt_gamma <- lapply(addt_gamma_vals, function(gam) {
recomputeSupercells(
dt = dat,
sc_objects = supercells$supercell_object,
markers = marker_cols_asinh,
sample_colname = sample_col,
cell_id_colname = cell_id_col,
gam = gam
)
})
We should end up with a list containing 2 elements. The 1st element contains supercells generated using gamma = 10, and the 2nd contains supercells generated using gamma = 50.
supercells_addt_gamma[[1]]
#> $supercell_expression_matrix
#> Marker_1_asinh Marker_2_asinh Marker_3_asinh Marker_4_asinh
#> <num> <num> <num> <num>
#> 1: 1.648545 0.8163302 2.091286 2.0528084
#> 2: 1.578355 0.9342016 2.101614 2.0707506
#> 3: 1.660627 0.6878140 2.099027 2.0797250
#> 4: 1.655574 0.9926162 2.088318 2.0754894
#> 5: 1.629174 1.0521980 2.080175 2.0559476
#> ---
#> 2996: 1.104590 1.8356044 1.706120 0.9802269
#> 2997: 1.151995 1.7035979 1.785026 0.8164285
#> 2998: 1.176526 1.7073983 1.788784 1.2843753
#> 2999: 1.179636 1.6677399 1.688395 0.7298266
#> 3000: 1.256157 1.7925856 1.701202 1.0127502
#> Marker_5_asinh Marker_6_asinh Marker_7_asinh Marker_8_asinh
#> <num> <num> <num> <num>
#> 1: 1.979651 1.943302 1.154513 1.1094316
#> 2: 2.007129 1.963738 1.065537 0.9125006
#> 3: 1.988715 1.950827 1.262522 1.1331555
#> 4: 1.993947 1.954081 1.150662 1.0077978
#> 5: 2.000412 1.961344 1.210783 0.9630321
#> ---
#> 2996: 2.050686 1.715387 1.555122 1.9664148
#> 2997: 1.993853 1.793681 1.358905 1.8631798
#> 2998: 2.033584 1.763948 1.612057 1.9195981
#> 2999: 2.021140 1.740146 1.595320 1.9590282
#> 3000: 2.007490 1.618998 1.307202 1.8714763
#> Marker_9_asinh Marker_10_asinh Marker_11_asinh Marker_12_asinh
#> <num> <num> <num> <num>
#> 1: 1.4689459 1.477175 1.942654 1.603274
#> 2: 1.4228915 1.326562 1.949072 1.626912
#> 3: 1.4841234 1.439633 1.958429 1.639451
#> 4: 1.5140737 1.454889 1.953806 1.669801
#> 5: 1.4366355 1.464401 1.946253 1.578554
#> ---
#> 2996: 1.2434630 2.088405 2.130623 2.061341
#> 2997: 1.1472844 2.050717 2.121512 1.971556
#> 2998: 1.2026242 1.966518 2.063505 2.042360
#> 2999: 0.9872839 2.034254 2.122262 2.066022
#> 3000: 0.8787599 1.974959 2.094728 1.973078
#> Marker_13_asinh Marker_14_asinh Marker_15_asinh Sample
#> <num> <num> <num> <char>
#> 1: 1.7258550 1.8211772 1.5094502 Sample_1
#> 2: 1.7838380 1.8515884 1.4872883 Sample_1
#> 3: 1.7467142 1.8478311 1.4517317 Sample_1
#> 4: 1.7655717 1.8455646 1.5035926 Sample_1
#> 5: 1.7113733 1.8587909 1.3655797 Sample_1
#> ---
#> 2996: 1.2592436 0.9414288 1.0298262 Sample_3
#> 2997: 0.9093374 1.1022433 0.8564922 Sample_3
#> 2998: 1.2870053 0.8414334 0.8732482 Sample_3
#> 2999: 1.1089350 0.9437904 1.1485098 Sample_3
#> 3000: 1.0930684 1.0466298 1.0368329 Sample_3
#> SuperCellId
#> <char>
#> 1: SuperCell_1_Sample_Sample_1
#> 2: SuperCell_2_Sample_Sample_1
#> 3: SuperCell_3_Sample_Sample_1
#> 4: SuperCell_4_Sample_Sample_1
#> 5: SuperCell_5_Sample_Sample_1
#> ---
#> 2996: SuperCell_996_Sample_Sample_3
#> 2997: SuperCell_997_Sample_Sample_3
#> 2998: SuperCell_998_Sample_Sample_3
#> 2999: SuperCell_999_Sample_Sample_3
#> 3000: SuperCell_1000_Sample_Sample_3
#>
#> $supercell_cell_map
#> SuperCellID CellId Sample
#> <char> <char> <char>
#> 1: SuperCell_618_Sample_Sample_1 Cell_1 Sample_1
#> 2: SuperCell_661_Sample_Sample_1 Cell_2 Sample_1
#> 3: SuperCell_395_Sample_Sample_1 Cell_3 Sample_1
#> 4: SuperCell_621_Sample_Sample_1 Cell_4 Sample_1
#> 5: SuperCell_894_Sample_Sample_1 Cell_5 Sample_1
#> ---
#> 29996: SuperCell_84_Sample_Sample_3 Cell_29996 Sample_3
#> 29997: SuperCell_475_Sample_Sample_3 Cell_29997 Sample_3
#> 29998: SuperCell_204_Sample_Sample_3 Cell_29998 Sample_3
#> 29999: SuperCell_279_Sample_Sample_3 Cell_29999 Sample_3
#> 30000: SuperCell_494_Sample_Sample_3 Cell_30000 Sample_3
The output generated by recomputeSupercells is essentially a list:
supercell_expression_matrix: A data.table object that contains the marker
expression for each supercell.supercell_cell_map: A data.table that maps each cell to its
corresponding supercell.As mentioned before, gamma dictates the granularity of supercells. Compared to the previous run where gamma was set to 20, we should get more supercells for gamma = 10, and less for gamma = 50. Let’s see if that’s the case.
n_supercells_gamma20 <- nrow(supercells$supercell_expression_matrix)
n_supercells_gamma10 <- nrow(
supercells_addt_gamma[[1]]$supercell_expression_matrix
)
n_supercells_gamma50 <- nrow(
supercells_addt_gamma[[2]]$supercell_expression_matrix
)
n_supercells_gamma10 > n_supercells_gamma20
#> [1] TRUE
n_supercells_gamma50 < n_supercells_gamma20
#> [1] TRUE
In the future, we may add the ability to specify different gam
value for different samples.
For now, if we want to do this, we will need to break down our data
into multiple data.table objects, each containing data from 1 sample,
and run runSuperCellCyto
function on each of them with different gam parameter value.
Something like the following:
n_markers <- 10
dat <- simCytoData(nmarkers = n_markers)
markers_col <- paste0("Marker_", seq_len(n_markers))
sample_col <- "Sample"
cell_id_col <- "Cell_Id"
samples <- unique(dat[[sample_col]])
gam_values <- c(10, 20, 10)
supercells_diff_gam <- lapply(seq_len(length(samples)), function(i) {
sample <- samples[i]
gam <- gam_values[i]
dat_samp <- dat[dat$Sample == sample, ]
supercell_samp <- runSuperCellCyto(
dt = dat_samp,
markers = markers_col,
sample_colname = sample_col,
cell_id_colname = cell_id_col,
gam = gam
)
return(supercell_samp)
})
Subsequently, to extract and combine the supercell_expression_matrix and
supercell_cell_map, we will need to use rbind:
supercell_expression_matrix <- do.call(
"rbind", lapply(
supercells_diff_gam, function(x) x[["supercell_expression_matrix"]]
)
)
supercell_cell_map <- do.call(
"rbind", lapply(
supercells_diff_gam, function(x) x[["supercell_cell_map"]]
)
)
rbind(
head(supercell_expression_matrix, n = 3),
tail(supercell_expression_matrix, n = 3)
)
#> Marker_1 Marker_2 Marker_3 Marker_4 Marker_5 Marker_6 Marker_7 Marker_8
#> <num> <num> <num> <num> <num> <num> <num> <num>
#> 1: 4.186726 10.811831 5.184253 15.54458 15.468691 15.64087 11.292195 12.51551
#> 2: 6.828708 9.039226 6.506419 15.66295 15.451165 15.26961 12.858617 14.79188
#> 3: 5.810587 9.978808 6.773288 15.63469 16.463423 15.64754 13.279402 12.93802
#> 4: 18.204621 11.758806 4.529428 18.96008 6.523332 12.47470 7.248837 14.97409
#> 5: 18.823289 11.266956 7.107383 16.99650 6.594569 12.21298 7.622260 14.44879
#> 6: 17.373217 11.584687 3.595834 19.05080 7.542699 12.44470 6.101374 14.55203
#> Marker_9 Marker_10 Sample SuperCellId
#> <num> <num> <char> <char>
#> 1: 9.439128 17.738778 Sample_1 SuperCell_1_Sample_Sample_1
#> 2: 10.452472 18.787770 Sample_1 SuperCell_2_Sample_Sample_1
#> 3: 9.681174 18.081401 Sample_1 SuperCell_3_Sample_Sample_1
#> 4: 8.352184 8.556783 Sample_2 SuperCell_498_Sample_Sample_2
#> 5: 8.494558 11.187681 Sample_2 SuperCell_499_Sample_Sample_2
#> 6: 9.739182 11.323073 Sample_2 SuperCell_500_Sample_Sample_2
rbind(head(supercell_cell_map, n = 3), tail(supercell_cell_map, n = 3))
#> SuperCellID CellId Sample
#> <char> <char> <char>
#> 1: SuperCell_777_Sample_Sample_1 Cell_1 Sample_1
#> 2: SuperCell_992_Sample_Sample_1 Cell_2 Sample_1
#> 3: SuperCell_721_Sample_Sample_1 Cell_3 Sample_1
#> 4: SuperCell_170_Sample_Sample_2 Cell_19998 Sample_2
#> 5: SuperCell_120_Sample_Sample_2 Cell_19999 Sample_2
#> 6: SuperCell_270_Sample_Sample_2 Cell_20000 Sample_2
If for whatever reason you don’t mind (or perhaps more to the point want)
each supercell to contain cells from different biological samples,
you still need to have the sample column in your data.table.
However, what you need to do is essentially set the value in the column
to exactly one unique value.
That way, SuperCellCyto will treat all cells as coming from one sample.
Just note, the parallel processing feature in SuperCellCyto won’t work for this as you will essentially only have 1 sample and nothing for SuperCellCyto to parallelise.
Is your dataset so huge that you are constantly running out of RAM when generating supercells? This thing happens and we have a solution for it.
Since supercells are generated for each sample independent of others you can easily break up the process. For example:
supercell_expression_matrix and supercell_cell_map,
and export them out as a csv file using data.table’s fwrite function.Once you have processed all the samples, you can then load all
supercell_expression_matrix and supercell_cell_map csv files and
analyse them.
If you want to regenerate the supercells using different gamma values,
load the relevant output saved using the qs package and the relevant data
(remember to note which output belongs to which sets of samples!), and run
recomputeSupercells function.
sessionInfo()
#> R version 4.5.1 Patched (2025-09-10 r88807)
#> Platform: x86_64-apple-darwin20
#> Running under: macOS Monterey 12.7.6
#>
#> Matrix products: default
#> BLAS: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.5-x86_64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1
#>
#> locale:
#> [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#>
#> time zone: America/New_York
#> tzcode source: internal
#>
#> attached base packages:
#> [1] parallel stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] BiocParallel_1.44.0 SuperCellCyto_1.0.0 BiocStyle_2.38.0
#>
#> loaded via a namespace (and not attached):
#> [1] cli_3.6.5 knitr_1.50 rlang_1.1.6
#> [4] xfun_0.53 jsonlite_2.0.0 data.table_1.17.8
#> [7] plyr_1.8.9 htmltools_0.5.8.1 sass_0.4.10
#> [10] rmarkdown_2.30 grid_4.5.1 evaluate_1.0.5
#> [13] jquerylib_0.1.4 fastmap_1.2.0 yaml_2.3.10
#> [16] lifecycle_1.0.4 bookdown_0.45 BiocManager_1.30.26
#> [19] compiler_4.5.1 igraph_2.2.1 codetools_0.2-20
#> [22] Rcpp_1.1.0 pkgconfig_2.0.3 lattice_0.22-7
#> [25] digest_0.6.37 SuperCell_1.0.1 R6_2.6.1
#> [28] RANN_2.6.2 magrittr_2.0.4 bslib_0.9.0
#> [31] Matrix_1.7-4 tools_4.5.1 cachem_1.1.0