SeSAMe implements inference of sex, age, ethnicity. These are valuable information for checking the integrity of the experiment and detecting sample swaps.

Sex, XCI

Sex is inferred based on our curated X-linked probes and Y chromosome probes excluding pseudo-autosomal regions and XCI escapes.

Human:

sdf = sesameDataGet('EPIC.1.SigDF')
inferSex(sdf)
## [1] "MALE"
inferSexKaryotypes(sdf)
## [1] "XaY"

Mouse:

sdf = sesameDataGet("MM285.1.SigDF")
inferSex(sdf)
## [1] Male
## Levels: Female Male

Ethnicity

Ethnicity is inferred using a random forest model trained based on both the built-in SNPs (rs probes) and channel-switching Type-I probes.

sdf = sesameDataGet('EPIC.1.SigDF')
inferEthnicity(sdf)
## [1] "WHITE"

Age & Epigenetic Clock

SeSAMe provides age regression through multiple previously established models, e.g., the well-known Horvath 353 model (Horvath 2013) which returns the chronological age in the number of years. Here is an example:

betas <- sesameDataGet('HM450.1.TCGA.PAAD')$betas
## download clock file from http://zwdzwd.github.io/InfiniumAnnotation
model <- readRDS("~/Downloads/Clock_Horvath353.rds")
predictAge(betas, model)

And MM285 mouse array data using a set of 347 CpGs (see Zhou et al. 2022) The function returns the age in the number of months. We recommend using SeSAMe preprocessed data as input to the function. Here’s an example:

library(SummarizedExperiment)
betas <- assay(sesameDataGet("MM285.10.SE.tissue"))[,1]
## download clock file from http://zwdzwd.github.io/InfiniumAnnotation
model <- readRDS("~/Downloads/Clock_Zhou347.rds")
predictAge(betas, model)

This indicates that this mouse is approximately 1.41 months old. The function looks for overlapping probes and estimates age using the corresponding clock models.

Copy Number

See Supplemental Vignette

Cell Count Deconvolution

SeSAMe estimates leukocyte fraction using a two-component model.This function works for samples whose targeted cell-of-origin is not related to white blood cells.

betas.tissue <- sesameDataGet('HM450.1.TCGA.PAAD')$betas
estimateLeukocyte(betas.tissue)
## [1] 0.2007592

Genomic Privacy

The goal of data sanitization is to modifiy IDAT files in place, so they can be released to public domain without privacy leak. This will be achieved by deIdentification.

Let’s take DNA methylation data from the HM450 platform for example. First, let’s download test IDATs from https://github.com/zhou-lab/InfiniumAnnotationV1/tree/main/Test

De-identify by Masking

This first method of deIdentification masks SNP probe intensity mean by zero. As a consequence, the allele frequency will be 0.5.

deIdentify("~/Downloads/3999492009_R01C01_Grn.idat",
    "~/Downloads/deidentified_Grn.idat")
deIdentify("~/Downloads/3999492009_R01C01_Red.idat",
    "~/Downloads/deidentified_Red.idat")

betas1 = getBetas(readIDATpair("~/Downloads/3999492009_R01C01"))
betas2 = getBetas(readIDATpair("~/Downloads/deidentified"))

head(betas1[grep('rs',names(betas1))]) 
head(betas2[grep('rs',names(betas2))])

Note that before deIdentify, the rs values will all be different. After deIdentify, the rs values will all be masked at an intensity of 0.5.

De-identify by Scrambling

This second method of deIdentification will scramble the intensities using a secret key to help formalize a random number. Therefore, randomize needs to be set to TRUE.

my_secret <- 13412084
set.seed(my_secret)

deIdentify("~/Downloads/3999492009_R01C01_Grn.idat",
    "~/Downloads/deidentified_Grn.idat", randomize=TRUE)

my_secret <- 13412084
set.seed(my_secret)
deIdentify("~/Downloads/3999492009_R01C01_Red.idat",
    "~/Downloads/deidentified_Red.idat", randomize=TRUE)

betas1 = getBetas(readIDATpair("~/Downloads/3999492009_R01C01"))
betas2 = getBetas(readIDATpair("~/Downloads/deidentified"))

head(betas1[grep('rs',names(betas1))]) 
head(betas2[grep('rs',names(betas2))]) 

Note that the rs values are scrambled after deIdentify.

Re-identify

To restore order of the deIdentified intensities, one can re-identify IDATs. The reIdentify function can thus restore the scrambled SNP intensities.

my_secret <- 13412084
set.seed(my_secret)

reIdentify(sprintf("%s/deidentified_Grn.idat", tmp),
    sprintf("%s/reidentified_Grn.idat", tmp))

my_secret <- 13412084
set.seed(my_secret)
reIdentify("~/Downloads/deidentified_Red.idat",
    "~/Downloads/reidentified_Red.idat")

betas1 = getBetas(readIDATpair("~/Downloads/3999492009_R01C01"))
betas2 = getBetas(readIDATpair("~/Downloads/reidentified"))

head(betas1[grep('rs',names(betas1))]) 
head(betas2[grep('rs',names(betas2))]) 

Note that reIdentify restored the values. Subsequently, they are the same as betas1.

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] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
##  [1] ggplot2_3.4.4               tibble_3.2.1               
##  [3] SummarizedExperiment_1.32.0 Biobase_2.62.0             
##  [5] MatrixGenerics_1.14.0       matrixStats_1.0.0          
##  [7] knitr_1.44                  S4Vectors_0.40.0           
##  [9] IRanges_2.36.0              GenomicRanges_1.54.0       
## [11] sesame_1.20.0               sesameData_1.19.0          
## [13] ExperimentHub_2.10.0        AnnotationHub_3.10.0       
## [15] BiocFileCache_2.10.0        dbplyr_2.3.4               
## [17] BiocGenerics_0.48.0        
## 
## loaded via a namespace (and not attached):
##  [1] DBI_1.1.3                     bitops_1.0-7                 
##  [3] rlang_1.1.1                   magrittr_2.0.3               
##  [5] e1071_1.7-13                  compiler_4.3.1               
##  [7] RSQLite_2.3.1                 png_0.1-8                    
##  [9] vctrs_0.6.4                   reshape2_1.4.4               
## [11] stringr_1.5.0                 pkgconfig_2.0.3              
## [13] crayon_1.5.2                  fastmap_1.1.1                
## [15] XVector_0.42.0                ellipsis_0.3.2               
## [17] fontawesome_0.5.2             labeling_0.4.3               
## [19] utf8_1.2.4                    promises_1.2.1               
## [21] rmarkdown_2.25                tzdb_0.4.0                   
## [23] preprocessCore_1.64.0         purrr_1.0.2                  
## [25] bit_4.0.5                     xfun_0.40                    
## [27] randomForest_4.7-1.1          zlibbioc_1.48.0              
## [29] cachem_1.0.8                  GenomeInfoDb_1.38.0          
## [31] jsonlite_1.8.7                blob_1.2.4                   
## [33] later_1.3.1                   DelayedArray_0.28.0          
## [35] BiocParallel_1.36.0           interactiveDisplayBase_1.40.0
## [37] parallel_4.3.1                R6_2.5.1                     
## [39] bslib_0.5.1                   stringi_1.7.12               
## [41] RColorBrewer_1.1-3            jquerylib_0.1.4              
## [43] Rcpp_1.0.11                   wheatmap_0.2.0               
## [45] readr_2.1.4                   httpuv_1.6.12                
## [47] Matrix_1.6-1.1                tidyselect_1.2.0             
## [49] abind_1.4-5                   yaml_2.3.7                   
## [51] codetools_0.2-19              curl_5.1.0                   
## [53] lattice_0.22-5                plyr_1.8.9                   
## [55] shiny_1.7.5.1                 withr_2.5.1                  
## [57] KEGGREST_1.42.0               evaluate_0.22                
## [59] proxy_0.4-27                  Biostrings_2.70.0            
## [61] pillar_1.9.0                  BiocManager_1.30.22          
## [63] filelock_1.0.2                stats4_4.3.1                 
## [65] generics_0.1.3                RCurl_1.98-1.12              
## [67] BiocVersion_3.18.0            hms_1.1.3                    
## [69] munsell_0.5.0                 scales_1.2.1                 
## [71] BiocStyle_2.30.0              xtable_1.8-4                 
## [73] class_7.3-22                  glue_1.6.2                   
## [75] tools_4.3.1                   grid_4.3.1                   
## [77] AnnotationDbi_1.64.0          colorspace_2.1-0             
## [79] GenomeInfoDbData_1.2.11       cli_3.6.1                    
## [81] rappdirs_0.3.3                fansi_1.0.5                  
## [83] S4Arrays_1.2.0                dplyr_1.1.3                  
## [85] gtable_0.3.4                  sass_0.4.7                   
## [87] digest_0.6.33                 SparseArray_1.2.0            
## [89] ggrepel_0.9.4                 farver_2.1.1                 
## [91] memoise_2.0.1                 htmltools_0.5.6.1            
## [93] lifecycle_1.0.3               httr_1.4.7                   
## [95] mime_0.12                     MASS_7.3-60                  
## [97] bit64_4.0.5