Constructs an anglemania_object-class from a given Seurat or SingleCellExperiment object. This includes extracting and processing count matrices, filtering genes based on expression in a minimum number of cells, and storing results along with dataset and batch information. It also calculates weights for each dataset or batch based on the number of samples.

create_anglemania_object(object, ...)

# S4 method for class 'Seurat'
create_anglemania_object(
  object,
  dataset_key = NA,
  batch_key,
  min_cells_per_gene = 1
)

# S4 method for class 'SingleCellExperiment'
create_anglemania_object(
  object,
  dataset_key = NA,
  batch_key,
  min_cells_per_gene = 1
)

# S4 method for class 'list'
create_anglemania_object(object, min_cells_per_gene = 1)

Arguments

object

A Seurat or SingleCellExperiment object containing single-cell RNA-seq data.

...

Additional arguments

dataset_key

A character string indicating the column name in the object metadata that identifies the dataset to which each cell belongs. If NA, all cells are assumed to belong to the same dataset.

batch_key

A character string indicating the column name(s) in the object metadata that identify the batch to which each cell belongs.

min_cells_per_gene

A numeric value indicating the minimum number of cells in which a gene must be expressed to be included in the analysis. Default is 1.

Value

An anglemania_object-class containing:

matrix_list

A list of filtered count matrices for each unique batch.

dataset_key

The dataset key used for splitting the object.

batch_key

The batch key used for splitting the object.

data_info

A data frame summarizing the number of samples per dataset and their weights.

weights

A numeric vector of weights for each dataset or batch based on the number of samples.

intersect_genes

A character vector of genes expressed in at least the specified number of cells across all batches.

min_cells_per_gene

The minimum number of cells per gene threshold used for filtering.

Details

The function performs the following steps:

  1. Adds a unique batch key to the metadata using add_unique_batch_key.

  2. Extracts count matrices for each batch.

  3. Filters genes based on the min_cells_per_gene threshold.

  4. Identifies intersected genes present across all batches.

  5. Converts count matrices to FBM objects.

  6. Computes weights for each batch or dataset.

Examples

sce <- sce_example()
angl <- create_anglemania_object(sce, batch_key = "batch")
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: batch
#> Extracting count matrices...
#> Filtering each batch to at least 1 cells per gene...
#> Using the intersection of filtered genes from all batches...
#> Number of genes in intersected set: 300
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: batch 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> batch1, batch2 
#> Number of intersected genes: 300 
#> Intersected genes (showing first 10):
#> gene1, gene2, gene3, gene4, gene5, gene6, gene7, gene8, gene9, gene10 , ...
#> Min cells per gene: 1 
sce <- sce_example()
se <- Seurat::as.Seurat(sce, data = "counts")
se <- SeuratObject::RenameAssays(se, "originalexp", "RNA")
#> Renaming default assay from originalexp to RNA
#> Warning: Key ‘originalexp_’ taken, using ‘rna_’ instead
angl <- create_anglemania_object(se, batch_key = "batch")
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: batch
#> Extracting count matrices...
#> Filtering each batch to at least 1 cells per gene...
#> Using the intersection of filtered genes from all batches...
#> Number of genes in intersected set: 300
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: batch 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> batch1, batch2 
#> Number of intersected genes: 300 
#> Intersected genes (showing first 10):
#> gene1, gene2, gene3, gene4, gene5, gene6, gene7, gene8, gene9, gene10 , ...
#> Min cells per gene: 1 
sce <- sce_example()
angl <- create_anglemania_object(sce, batch_key = "batch")
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: batch
#> Extracting count matrices...
#> Filtering each batch to at least 1 cells per gene...
#> Using the intersection of filtered genes from all batches...
#> Number of genes in intersected set: 300
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: batch 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> batch1, batch2 
#> Number of intersected genes: 300 
#> Intersected genes (showing first 10):
#> gene1, gene2, gene3, gene4, gene5, gene6, gene7, gene8, gene9, gene10 , ...
#> Min cells per gene: 1 
sce <- sce_example()
sce_list <- list(sce1 = sce, sce2 = sce)
angl <- create_anglemania_object(sce_list)
#> Extracting count matrices...
#> Filtering each batch to at least 1 cells per gene...
#> Using the intersection of filtered genes from all batches...
#> Number of genes in intersected set: 300
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: batch 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> sce1, sce2 
#> Number of intersected genes: 300 
#> Intersected genes (showing first 10):
#> gene1, gene2, gene3, gene4, gene5, gene6, gene7, gene8, gene9, gene10 , ...
#> Min cells per gene: 1 
se <- Seurat::as.Seurat(sce, data = "counts")
se <- SeuratObject::RenameAssays(se, "originalexp", "RNA")
#> Renaming default assay from originalexp to RNA
#> Warning: Key ‘originalexp_’ taken, using ‘rna_’ instead
se_list <- list(se1 = se, se2 = se)
angl <- create_anglemania_object(se_list)
#> Extracting count matrices...
#> Filtering each batch to at least 1 cells per gene...
#> Using the intersection of filtered genes from all batches...
#> Number of genes in intersected set: 300
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: batch 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> se1, se2 
#> Number of intersected genes: 300 
#> Intersected genes (showing first 10):
#> gene1, gene2, gene3, gene4, gene5, gene6, gene7, gene8, gene9, gene10 , ...
#> Min cells per gene: 1