cpg.annotate | R Documentation |
Annotate a matrix/GenomicRatioSet representing EPICv2, EPICv1 or 450K data with probe weights and chromosomal position. Provides replicate filtering and remapping functions for EPICv2 probes.
cpg.annotate(datatype = c("array", "sequencing"), object,
what = c("Beta", "M"), arraytype = c("EPICv2", "EPICv1", "EPIC",
"450K"), epicv2Remap = TRUE, epicv2Filter = c("mean",
"sensitivity", "precision", "random"), analysis.type =
c("differential", "variability", "ANOVA", "diffVar"),
design, contrasts = FALSE, cont.matrix = NULL, fdr = 0.05, coef,
varFitcoef = NULL, topVarcoef = NULL, ...)
datatype |
Character string representing the type of data being analysed. |
object |
Either: - A matrix of M-values, with unique Illumina probe IDs as rownames and unique sample IDs as column names or, - A GenomicRatioSet, appropriately annotated. |
what |
Does the data matrix contain Beta or M-values? Not needed if object is a GenomicRatioSet. |
arraytype |
Is the data matrix sourced from EPIC or 450K data? Not needed if object is a GenomicRatioSet. |
epicv2Remap |
Logical indicating whether to remap 11,878 cross-hybridising EPICv2 probes to their more likely CpG target (see Peters et al. 2024). |
epicv2Filter |
Strategy for filtering probe replicates that map to the same CpG site.
|
analysis.type |
|
design |
Study design matrix. Identical context to differential analysis
pipeline in |
contrasts |
Logical denoting whether a |
cont.matrix |
|
fdr |
FDR cutoff (Benjamini-Hochberg) for which CpG sites are individually called
as significant. Used to index default thresholding in dmrcate(). Highly
recommended as the primary thresholding parameter for calling DMRs.
Not used when |
coef |
The column index in |
varFitcoef |
The columns of the design matrix containing the comparisons to test for
differential variability. If left |
topVarcoef |
Column number or column name specifying which coefficient of the linear
model fit is of interest. It should be the same coefficient that
the differential variability testing was performed on. Default is last
column of fit object. Identical context to |
... |
Extra arguments passed to the |
A CpGannotated-class
.
Tim J. Peters <t.peters@garvan.org.au>
Ritchie, M. E., Phipson, B., Wu, D., Hu, Y., Law, C. W., Shi, W., & Smyth, G. K. (2015). limma powers differential expression analyses for RNA-sequencing and microarray studies. Nucleic Acids Research, 43(7), e47.
Phipson, B., & Oshlack, A. (2014). DiffVar: a new method for detecting differential variability with application to methylation in cancer and aging. Genome Biol, 15(9), 465.
Peters T.J., Buckley M.J., Statham, A., Pidsley R., Samaras K., Lord R.V., Clark S.J. and Molloy P.L. De novo identification of differentially methylated regions in the human genome. Epigenetics & Chromatin 2015, 8:6, doi:10.1186/1756-8935-8-6.
Peters, T.J., Meyer, B., Ryan, L., Achinger-Kawecka, J., Song, J., Campbell, E.M., Qu, W., Nair, S., Loi-Luu, P., Stricker, P., Lim, E., Stirzaker, C., Clark, S.J. and Pidsley, R. (2024). Characterisation and reproducibility of the HumanMethylationEPIC v2.0 BeadChip for DNA methylation profiling. BMC Genomics, 25, 251. doi:10.1186/s12864-024-10027-5.
library(AnnotationHub)
ah <- AnnotationHub()
EPICv2manifest <- ah[["AH116484"]]
object <- minfi::logit2(matrix(rbeta(10000, 3, 1), 1000, 10))
rownames(object) <- sample(rownames(EPICv2manifest), 1000)
type <- rep(c("Ctrl", "Treat"), each=5)
design <- model.matrix(~type)
myannotation <- cpg.annotate("array", object, what = "M", arraytype = "EPICv2",
analysis.type="differential", design=design, coef=2)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.