As sesame and sesameData are under active development, this documentation is specific to the following version of R, sesame, sesameData and ExperimentHub:
## SeSAMe requires matched versions of R, sesame, sesameData and ExperimentHub.
## Here is the current versions installed:
## R: 4.3.0
## Bioconductor: 3.17
## sesame: 1.18.4
## sesameData: 1.18.0
## ExperimentHub: 2.8.0
We recommend updating your R, ExperimentHub, sesameData and sesame to use this documentation consistently. If you have installed directly from github, please make sure the compatible ExperimentHub is installed.
If you use a previous version, please checkout the vignette that corresponds to the right version here https://zhou-lab.github.io/sesame/dev/supplemental.html#Versions
CRITICAL: After a new installation, one must cache the associated annotation data using the following command. This needs to be done only once per SeSAMe installation/update. Caching data to local guarantees proper data retrieval and saves internet traffic.
This function caches the needed SeSAMe annotations. SeSAMe annotation data is managed by the sesameData package which uses the ExperimentHub infrastructure. You can find the location of the cached annotation data on your local computer using:
## [1] "/home/biocbuild/.cache/R/ExperimentHub"
The openSesame
function provides end-to-end processing that converts IDATs to DNA methylation level (aka β value) matrices in R. The function can take either one of the following input:
The following code uses a directory that contains built-in two HM27 IDAT pairs to demonstrates the use of openSesame
:
## The following use the idat_dir in system.file,
## Replace it with your actual IDAT path
idat_dir = system.file("extdata/", package = "sesameData")
betas = openSesame(idat_dir, BPPARAM = BiocParallel::MulticoreParam(2))
The BPPARAM=
option is from the BiocParallel
package and controls parallel processing (in this case, we are using two cores). Under the hood, the function performs a series of tasks including: searching IDAT files from the directory (the searchIDATprefixes
function), reading IDAT data in as SigDF
objects (the readIDATpair
function), preprocessing the signals (the prepSesame
function), and finally converting them to DNA methylation levels (β values, the getBetas
function). Alternatively, one can run the following command to get the same results, while gaining more refined control:
## The above openSesame call is equivalent to:
betas = do.call(cbind, BiocParallel::bplapply(
searchIDATprefixes(idat_dir), function(pfx) {
getBetas(prepSesame(readIDATpair(pfx), "QCDPB"))
}, BPPARAM = BiocParallel::MulticoreParam(2)))
## or even more explicitly (if one needs to control argument passed
## to a specific preprocessing function)
betas = do.call(cbind, BiocParallel::bplapply(
searchIDATprefixes(idat_dir), function(pfx) {
getBetas(noob(pOOBAH(dyeBiasNL(inferInfiniumIChannel(qualityMask(
readIDATpair(pfx)))))))
}, BPPARAM = BiocParallel::MulticoreParam(2)))
The openSesame
function is highly customizable. The prep=
argument is the same argument one gives to the prepSesame
function (see Data Preprocessing for detail) which openSesame
calls internally. The argument uniquely specifies a preprocessing procedure. The func=
option specifies the signal extraction function. It can be either be getBetas
(DNA methylation) or getAFs
(allele frequencies of SNP probes) or NULL (returns SigDF
). The manifest=
option allows one to provide an array manifest when handling data from platform not supported natively. Finally, the BPPARAM=
argument is the same argument taken by BiocParallel::bplapply
to allow parallel processing. See Supplemental Vignette for details of these component functions of openSesame.
The output of openSesame
can also be customized. It can either be beta values, which are the end DNA methylation readings, as shown above. It can also be a list of SigDF
s which stores the signal intensities and can be further put back to openSesame for more processing. The openSesame(func=)
argument specifies whether the output is a SigDF list or beta values. The following shows some usage:
The prep=
argument instructs the openSesame
function to call the prepSesame
function to preprocess signal intensity under the hood. This can be skipped by using prep=""
. The prepSesame
function takes a single SigDF
as input and returns a processed SigDF
. When prep=
is non-empty, it selects the preprocessing functions (see Preprocessing Function Code) and specifies the order of their execution. For example,
performs dye bias correction (D
) followed by background subtraction (B
). In other words, prepSesame(sdf, "DB")
is equivalent to noob(dyeBiasNL(sdf))
. All the preprocessing functions take a SigDF
as input and return an updated SigDF
. Therefore, these functions can be chained together. The choice of preprocessing functions and the order of their chaining is important (see Supplemental Vignette) for detailed discussions of these functions). The following table lists the best preprocessing strategy based on our experience.
Platform | Sample Organism | Prep Code |
---|---|---|
EPICv2/EPIC/HM450 | human | QCDPB |
EPICv2/EPIC/HM450 | non-human organism | SQCDPB |
MM285 | mouse | TQCDPB |
MM285 | non-mouse organism | SQCDPB |
Mammal40 | human | HCDPB |
Mammal40 | non-human organism | SHCDPB |
The optimal strategy of preprocessing depends on:
The array platform. For example, certain array platforms (e.g., the Mammal40) do not have enough Infinium-I probes for background estimation and dye bias correction, therefore background subtraction (where the out-of-band signals are from) might not work most optimally;
The expected sample property. For example, some samples have the signature bimodal distribution of methylation of most mammalian cells. Others may undergo global loss of methylation (germ cells, tumors etc). Other important factors include high-input vs low-input, tumor vs normal, somatic vs germ cells, human vs model organisms, mouse strains etc. Some platforms (e.g., Mammal40 and MM285) are designed for multiple species and strains. Therefore S
and T
would be important when those arrays are used on non-reference organisms (see Working with Nonhuman Arrays).
The prepSesameList
function lists all the available codes and the associated preprocessing functions.
Here are some consideration when determining the preprocessing order. Species (S)
and strain (T)
inference resets the mask and color channels based on probe alignment and presence of genetic variants. Therefore when they are used, they need to be called first. Q
masks non-uniquely mapped probes which may inflate the out-of-band signal for background estimation. Therefore Q
should be used before detection p-value calculation (P
) and background subtraction (B
) when necessary. Channel inference (C
) and dye bias correction (D
) should take place early since dye bias effect is global. C
should be placed before D
because dye bias correction uses in-band signal the identification of which relies on correct channel designation. Detection p-value (P
) should happen before background subtraction (B
) since background subtraction modifies signal and may affect out-of-band signal assumption used in P
. Lastly, functions that explicitly normalizes β value distribution (M
) should happen last if they even need to be used.
See Supplemental Vignette for details of preprocessing functions.
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.17-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] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] wheatmap_0.2.0 tidyr_1.3.0
## [3] dplyr_1.1.2 ggplot2_3.4.2
## [5] tibble_3.2.1 SummarizedExperiment_1.30.1
## [7] Biobase_2.60.0 GenomicRanges_1.52.0
## [9] GenomeInfoDb_1.36.0 IRanges_2.34.0
## [11] S4Vectors_0.38.1 MatrixGenerics_1.12.0
## [13] matrixStats_0.63.0 knitr_1.42
## [15] sesame_1.18.4 sesameData_1.18.0
## [17] ExperimentHub_2.8.0 AnnotationHub_3.8.0
## [19] BiocFileCache_2.8.0 dbplyr_2.3.2
## [21] BiocGenerics_0.46.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.0
## [7] RSQLite_2.3.1 mgcv_1.8-42
## [9] png_0.1-8 vctrs_0.6.2
## [11] reshape2_1.4.4 stringr_1.5.0
## [13] pkgconfig_2.0.3 crayon_1.5.2
## [15] fastmap_1.1.1 XVector_0.40.0
## [17] ellipsis_0.3.2 fontawesome_0.5.1
## [19] labeling_0.4.2 utf8_1.2.3
## [21] promises_1.2.0.1 rmarkdown_2.21
## [23] tzdb_0.3.0 preprocessCore_1.62.1
## [25] purrr_1.0.1 bit_4.0.5
## [27] xfun_0.39 randomForest_4.7-1.1
## [29] zlibbioc_1.46.0 cachem_1.0.8
## [31] jsonlite_1.8.4 blob_1.2.4
## [33] highr_0.10 later_1.3.1
## [35] DelayedArray_0.26.2 BiocParallel_1.34.1
## [37] interactiveDisplayBase_1.38.0 parallel_4.3.0
## [39] R6_2.5.1 bslib_0.4.2
## [41] stringi_1.7.12 RColorBrewer_1.1-3
## [43] jquerylib_0.1.4 Rcpp_1.0.10
## [45] readr_2.1.4 splines_4.3.0
## [47] httpuv_1.6.9 Matrix_1.5-4
## [49] tidyselect_1.2.0 yaml_2.3.7
## [51] codetools_0.2-19 curl_5.0.0
## [53] lattice_0.21-8 plyr_1.8.8
## [55] shiny_1.7.4 withr_2.5.0
## [57] KEGGREST_1.40.0 evaluate_0.20
## [59] proxy_0.4-27 Biostrings_2.68.0
## [61] pillar_1.9.0 BiocManager_1.30.20
## [63] filelock_1.0.2 generics_0.1.3
## [65] RCurl_1.98-1.12 BiocVersion_3.17.1
## [67] hms_1.1.3 munsell_0.5.0
## [69] scales_1.2.1 BiocStyle_2.28.0
## [71] xtable_1.8-4 class_7.3-22
## [73] glue_1.6.2 tools_4.3.0
## [75] grid_4.3.0 AnnotationDbi_1.62.1
## [77] colorspace_2.1-0 nlme_3.1-162
## [79] GenomeInfoDbData_1.2.10 cli_3.6.1
## [81] rappdirs_0.3.3 fansi_1.0.4
## [83] S4Arrays_1.0.1 gtable_0.3.3
## [85] sass_0.4.6 digest_0.6.31
## [87] ggrepel_0.9.3 farver_2.1.1
## [89] memoise_2.0.1 htmltools_0.5.5
## [91] lifecycle_1.0.3 httr_1.4.5
## [93] mime_0.12 MASS_7.3-60
## [95] bit64_4.0.5