Skip to contents

This method provides a concise summary of an anglemania_object, including dataset and batch information, the number of intersected genes, and other relevant details.

Retrieves the list of gene expression matrices stored in the anglemania_object object.

Assigns a new list of gene expression matrices to the anglemania_object.

Retrieves the dataset key used in the anglemania_object.

Retrieves the batch key used in the anglemania_object.

Retrieves the data frame summarizing the selected anglemania gene pairs anglemania_object.

Retrieves the weights assigned to each dataset or batch in the anglemania_object

Assigns new weights to the datasets or batches in the anglemania_object.

Retrieves the list of statistical measures computed across datasets in the anglemania_object.

Assigns a new list of statistical measures to the anglemania_object.

Retrieves the vector of genes that are expressed in at least the specified number of cells across all batches.

Assigns a new vector of intersected genes to the anglemania_object.

Retrieves the list of genes selected for integration from the anglemania_object

This function adds a unique batch identifier to the metadata of a Seurat object by combining specified dataset and batch keys. This is useful for distinguishing samples during integration or analysis.

Usage

# S4 method for class 'anglemania_object'
show(object)

matrix_list(object)

# S4 method for class 'anglemania_object'
matrix_list(object)

matrix_list(object) <- value

# S4 method for class 'anglemania_object'
matrix_list(object) <- value

dataset_key(object)

# S4 method for class 'anglemania_object'
dataset_key(object)

batch_key(object)

# S4 method for class 'anglemania_object'
batch_key(object)

data_info(object)

# S4 method for class 'anglemania_object'
data_info(object)

angl_weights(object)

# S4 method for class 'anglemania_object'
angl_weights(object)

angl_weights(object) <- value

# S4 method for class 'anglemania_object'
angl_weights(object) <- value

list_stats(object)

# S4 method for class 'anglemania_object'
list_stats(object)

list_stats(object) <- value

# S4 method for class 'anglemania_object'
list_stats(object) <- value

intersect_genes(object)

# S4 method for class 'anglemania_object'
intersect_genes(object)

intersect_genes(object) <- value

# S4 method for class 'anglemania_object'
intersect_genes(object) <- value

get_anglemania_genes(object)

# S4 method for class 'anglemania_object'
get_anglemania_genes(object)

add_unique_batch_key(
  seurat_object,
  dataset_key = NA_character_,
  batch_key,
  new_unique_batch_key = "batch"
)

Arguments

object

An anglemania_object.

value

A character vector of gene names.

seurat_object

A Seurat object.

dataset_key

A character string specifying the column name in the metadata that identifies the dataset. If NA, only the batch_key is used.

batch_key

A character string specifying the column name in the metadata that identifies the batch.

new_unique_batch_key

A character string for the new unique batch key to be added to the metadata. Default is "batch".

Value

Prints a summary to the console.

A list of FBM objects containing gene expression matrices.

The updated anglemania_object.

A character string representing the dataset key.

A character string representing the batch key.

A data frame containing dataset and batch information.

A named numeric vector of weights.

The updated anglemania_object.

A list containing statistical matrices such as mean z-scores and SNR z-scores

The updated anglemania_object.

A character vector of intersected gene names from multiple Seurat objects.

The updated anglemania_object object.

A character vector of integration gene names.

A Seurat object with an additional metadata column containing the unique batch key.

Functions

  • show(anglemania_object): show anglemania_object info

  • matrix_list(): Access matrix list

  • matrix_list(object) <- value: set matrix list in anglemania_object

  • dataset_key(): Access dataset key of anglemania_object

  • batch_key(): Access batch key of anglemania_object

  • data_info(): Access info of selected gene pairs

  • angl_weights(): Access weights

  • angl_weights(object) <- value: Set weights

  • list_stats(): Access statistics of the gene-gene matrices

  • list_stats(object) <- value: Set statistics of the gene-gene matrices

  • intersect_genes(): Access the intersection of genes of all batches

  • intersect_genes(object) <- value: Set the intersection of genes of all batches

  • get_anglemania_genes(): Access the genes extracted by anglemania

  • add_unique_batch_key(): Temporarily add a unique batch key to the dataset

Examples

se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(se, batch_key = "groups")
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl <- anglemania(angl)
#> Computing angles and transforming to z-scores...
#> 
  |                                                  | 0 % elapsed=00s   
  |=========================                         | 50% elapsed=00s, remaining~00s
  |==================================================| 100% elapsed=00s, remaining~00s
#> Computing statistics...
#> Weighting matrix_list...
#> Calculating mean...
#> Calculating sds...
#> Filtering features...
show(angl)
#> anglemania_object
#> --------------
#> Dataset key: NA 
#> Batch key: groups 
#> Number of datasets: 1 
#> Total number of batches: 2 
#> Batches (showing first 5):
#> g2, g1 
#> Number of intersected genes: 228 
#> Intersected genes (showing first 10):
#> MS4A1, CD79B, CD79A, HLA-DRA, TCL1A, HLA-DQB1, HVCN1, HLA-DMB, LTB, LINC00926 , ...
#> Min cells per gene: 1 
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(se, batch_key = "groups")
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl <- anglemania(angl)
#> Computing angles and transforming to z-scores...
#> 
  |                                                  | 0 % elapsed=00s   
  |=========================                         | 50% elapsed=00s, remaining~00s
  |==================================================| 100% elapsed=00s, remaining~00s
#> Computing statistics...
#> Weighting matrix_list...
#> Calculating mean...
#> Calculating sds...
#> Filtering features...
str(matrix_list(angl))
#> List of 2
#>  $ g1:Reference class 'FBM' [package "bigstatsr"] with 15 fields
#>   ..$ extptr      :<externalptr> 
#>   ..$ extptr_rw   :<externalptr> 
#>   ..$ nrow        : int 228
#>   ..$ ncol        : int 228
#>   ..$ type        : Named int 8
#>   .. ..- attr(*, "names")= chr "double"
#>   ..$ backingfile : chr "/tmp/Rtmpb1qOQN/file170928707d7a.bk"
#>   ..$ is_read_only: logi FALSE
#>   ..$ address     :<externalptr> 
#>   ..$ address_rw  :<externalptr> 
#>   ..$ bk          : chr "/tmp/Rtmpb1qOQN/file170928707d7a.bk"
#>   ..$ rds         : chr "/tmp/Rtmpb1qOQN/file170928707d7a.rds"
#>   ..$ is_saved    : logi FALSE
#>   ..$ type_chr    : chr "double"
#>   ..$ type_size   : int 8
#>   ..$ file_size   : num 415872
#>   ..and 22 methods, of which 8 are  possibly relevant:
#>   ..  add_columns, bm, bm.desc, check_dimensions, check_write_permissions,
#>   ..  initialize, save, show#envRefClass
#>  $ g2:Reference class 'FBM' [package "bigstatsr"] with 15 fields
#>   ..$ extptr      :<externalptr> 
#>   ..$ extptr_rw   :<externalptr> 
#>   ..$ nrow        : int 228
#>   ..$ ncol        : int 228
#>   ..$ type        : Named int 8
#>   .. ..- attr(*, "names")= chr "double"
#>   ..$ backingfile : chr "/tmp/Rtmpb1qOQN/file1709dda5d42.bk"
#>   ..$ is_read_only: logi FALSE
#>   ..$ address     :<externalptr> 
#>   ..$ address_rw  :<externalptr> 
#>   ..$ bk          : chr "/tmp/Rtmpb1qOQN/file1709dda5d42.bk"
#>   ..$ rds         : chr "/tmp/Rtmpb1qOQN/file1709dda5d42.rds"
#>   ..$ is_saved    : logi FALSE
#>   ..$ type_chr    : chr "double"
#>   ..$ type_size   : int 8
#>   ..$ file_size   : num 415872
#>   ..and 22 methods, of which 8 are  possibly relevant:
#>   ..  add_columns, bm, bm.desc, check_dimensions, check_write_permissions,
#>   ..  initialize, save, show#envRefClass
se <- SeuratObject::pbmc_small
se[[]]$Dataset <- rep(c("A", "B"), each = ncol(se) / 2)
angl <- create_anglemania_object(
  se,
  dataset_key = "Dataset",
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> Using dataset_key: Dataset
#> 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: 156
#> 
  |                                                  | 0 % elapsed=00s   
  |=========================                         | 50% elapsed=00s, remaining~00s
  |==================================================| 100% elapsed=00s, remaining~00s
dataset_key(angl)
#> [1] "Dataset"
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
batch_key(angl)
#> [1] "groups"
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
batch_key(angl)
#> [1] "groups"
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
batch_key(angl)
#> [1] "groups"
angl_weights(angl)
#>  g2  g1 
#> 0.5 0.5 
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl <- anglemania(angl)
#> Computing angles and transforming to z-scores...
#> 
  |                                                  | 0 % elapsed=00s   
  |=========================                         | 50% elapsed=00s, remaining~00s
  |==================================================| 100% elapsed=01s, remaining~00s
#> Computing statistics...
#> Weighting matrix_list...
#> Calculating mean...
#> Calculating sds...
#> Filtering features...
stats <- list_stats(angl)
str(stats)
#> List of 3
#>  $ mean_zscore: num [1:228, 1:228] NA 4.02 4.88 2.85 4.48 ...
#>  $ sds_zscore : num [1:228, 1:228] NA 0.601 0.21 0.461 0.86 ...
#>  $ sn_zscore  : num [1:228, 1:228] NA 6.69 23.19 6.17 5.21 ...

se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
intersect_genes(angl)
#>   [1] "MS4A1"         "CD79B"         "CD79A"         "HLA-DRA"      
#>   [5] "TCL1A"         "HLA-DQB1"      "HVCN1"         "HLA-DMB"      
#>   [9] "LTB"           "LINC00926"     "FCER2"         "SP100"        
#>  [13] "NCF1"          "PPP3CC"        "EAF2"          "PPAPDC1B"     
#>  [17] "CD19"          "KIAA0125"      "CYB561A3"      "CD180"        
#>  [21] "RP11-693J15.5" "FAM96A"        "CXCR4"         "STX10"        
#>  [25] "SNHG7"         "NT5C"          "BANK1"         "IGLL5"        
#>  [29] "CD200"         "FCRLA"         "CD3D"          "NOSIP"        
#>  [33] "SAFB2"         "CD2"           "IL7R"          "PIK3IP1"      
#>  [37] "MPHOSPH6"      "KHDRBS1"       "MAL"           "CCR7"         
#>  [41] "THYN1"         "TAF7"          "LDHB"          "TMEM123"      
#>  [45] "CCDC104"       "EPC1"          "EIF4A2"        "CD3E"         
#>  [49] "TMUB1"         "BLOC1S4"       "ACSM3"         "TMEM204"      
#>  [53] "SRSF7"         "ACAP1"         "TNFAIP8"       "CD7"          
#>  [57] "TAGAP"         "DNAJB1"        "ASNSD1"        "S1PR4"        
#>  [61] "CTSW"          "GZMK"          "NKG7"          "IL32"         
#>  [65] "DNAJC2"        "LYAR"          "CST7"          "LCK"          
#>  [69] "CCL5"          "HNRNPH1"       "SSR2"          "DLGAP1-AS1"   
#>  [73] "GIMAP1"        "MMADHC"        "ZNF76"         "CD8A"         
#>  [77] "PTPN22"        "GYPC"          "HNRNPF"        "RPL7L1"       
#>  [81] "KLRG1"         "CRBN"          "SATB1"         "SIT1"         
#>  [85] "PMPCB"         "NRBP1"         "TCF7"          "HNRNPA3"      
#>  [89] "S100A8"        "S100A9"        "LYZ"           "CD14"         
#>  [93] "FCN1"          "TYROBP"        "ASGR1"         "NFKBIA"       
#>  [97] "TYMP"          "CTSS"          "TSPO"          "RBP7"         
#> [101] "CTSB"          "LGALS1"        "FPR1"          "VSTM1"        
#> [105] "BLVRA"         "MPEG1"         "BID"           "SMCO4"        
#> [109] "CFD"           "LINC00936"     "LGALS2"        "MS4A6A"       
#> [113] "FCGRT"         "LGALS3"        "NUP214"        "SCO2"         
#> [117] "IL17RA"        "IFI6"          "HLA-DPA1"      "FCER1A"       
#> [121] "CLEC10A"       "HLA-DMA"       "RGS1"          "HLA-DPB1"     
#> [125] "HLA-DQA1"      "RNF130"        "HLA-DRB5"      "HLA-DRB1"     
#> [129] "CST3"          "IL1B"          "POP7"          "HLA-DQA2"     
#> [133] "CD1C"          "GSTP1"         "EIF3G"         "VPS28"        
#> [137] "LY86"          "ZFP36L1"       "ANXA2"         "GRN"          
#> [141] "CFP"           "HSP90AA1"      "LST1"          "AIF1"         
#> [145] "PSAP"          "YWHAB"         "MYO1G"         "SAT1"         
#> [149] "RGS2"          "SERPINA1"      "IFITM3"        "FCGR3A"       
#> [153] "LILRA3"        "S100A11"       "FCER1G"        "TNFRSF1B"     
#> [157] "IFITM2"        "WARS"          "IFI30"         "MS4A7"        
#> [161] "C5AR1"         "HCK"           "COTL1"         "LGALS9"       
#> [165] "CD68"          "RP11-290F20.3" "RHOC"          "CARD16"       
#> [169] "LRRC25"        "COPS6"         "ADAR"          "PPBP"         
#> [173] "GPX1"          "TPM4"          "PF4"           "SDPR"         
#> [177] "NRGN"          "SPARC"         "GNG11"         "CLU"          
#> [181] "HIST1H2AC"     "NCOA4"         "GP9"           "FERMT3"       
#> [185] "ODC1"          "CD9"           "RUFY1"         "TUBB1"        
#> [189] "TALDO1"        "TREML1"        "NGFRAP1"       "PGRMC1"       
#> [193] "CA2"           "ITGA2B"        "MYL9"          "TMEM40"       
#> [197] "PARVB"         "PTCRA"         "ACRBP"         "TSC22D1"      
#> [201] "VDAC3"         "GZMB"          "GZMA"          "GNLY"         
#> [205] "FGFBP2"        "AKR1C3"        "CCL4"          "PRF1"         
#> [209] "GZMH"          "XBP1"          "GZMM"          "PTGDR"        
#> [213] "IGFBP7"        "TTC38"         "KLRD1"         "ARHGDIA"      
#> [217] "IL2RB"         "CLIC3"         "PPP1R18"       "CD247"        
#> [221] "ALOX5AP"       "XCL2"          "C12orf75"      "RARRES3"      
#> [225] "PCMT1"         "LAMP1"         "SPON2"         "S100B"        
se <- SeuratObject::pbmc_small
angl <- create_anglemania_object(
  se,
  batch_key = "groups",
  min_cells_per_gene = 1
)
#> No dataset_key specified.
#> Assuming that all samples belong to the same dataset and are separated by batch_key: groups
#> 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: 228
#> 
  |                                                  | 0 % elapsed=00s   
  |==================================================| 100% elapsed=00s, remaining~00s
angl <- anglemania(angl)
#> Computing angles and transforming to z-scores...
#> 
  |                                                  | 0 % elapsed=00s   
  |=========================                         | 50% elapsed=00s, remaining~00s
  |==================================================| 100% elapsed=00s, remaining~00s
#> Computing statistics...
#> Weighting matrix_list...
#> Calculating mean...
#> Calculating sds...
#> Filtering features...
# extract the genes identified by anglemania()
anglemania_genes <- get_anglemania_genes(angl)

se <- SeuratObject::pbmc_small
se[[]]$Dataset <- rep(c("A", "B"), each = ncol(se)/2)
se <- add_unique_batch_key(
  seurat_object = se,
  dataset_key = "Dataset",
  batch_key = "groups",
  new_unique_batch_key = "batch" 
  )
head(se[[]])
#>                   orig.ident nCount_RNA nFeature_RNA RNA_snn_res.0.8
#> ATGCCAGAACGACT SeuratProject         70           47               0
#> CATGGCCTGTGCAT SeuratProject         85           52               0
#> GAACCTGATGAACC SeuratProject         87           50               1
#> TGACTGGATTCTCA SeuratProject        127           56               0
#> AGTCAGACTGCACA SeuratProject        173           53               0
#> TCTGATACACGTGT SeuratProject         70           48               0
#>                letter.idents groups RNA_snn_res.1 Dataset batch
#> ATGCCAGAACGACT             A     g2             0       A  A:g2
#> CATGGCCTGTGCAT             A     g1             0       A  A:g1
#> GAACCTGATGAACC             B     g2             0       A  A:g2
#> TGACTGGATTCTCA             A     g2             0       A  A:g2
#> AGTCAGACTGCACA             A     g2             0       A  A:g2
#> TCTGATACACGTGT             A     g1             0       A  A:g1