View source: R/sequencing.annotate.R
sequencing.annotate | R Documentation |
Either:
- Annotate a BSseq object with chromosome position and test statistic, or
- Parse output from DSS::DMLtest()
or DSS::DMLtest.multiFactor()
into a CpGannotated object.
sequencing.annotate(obj, methdesign, all.cov=FALSE, contrasts = FALSE,
cont.matrix = NULL, fdr = 0.05, coef, ...)
obj |
A BSseq object or data.frame output from |
methdesign |
Methylation study design matrix describing samples and groups. Use of
edgeR::modelMatrixMeth() to make this matrix is highly recommended, since it
transforms a regular model.matrix (as one would construct for a microarray or
RNA-Seq experiment) into a “two-channel” matrix representing methylated and
unmethylated reads for each sample. Only applicable when |
all.cov |
If |
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.
Only applicable when |
coef |
The column index in |
... |
Extra arguments passed to the |
A CpGannotated-class
.
Tim J. Peters <t.peters@garvan.org.au>
Peters, T. J., Buckley, M.J., Chen, Y., Smyth, G.K., Goodnow, C. C. and Clark, S. J. (2021). Calling differentially methylated regions from whole genome bisulphite sequencing with DMRcate. Nucleic Acids Research, 49(19), e109.
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.
library(ExperimentHub)
library(SummarizedExperiment)
library(bsseq)
library(GenomeInfoDb)
eh = ExperimentHub()
bis_1072 <- eh[["EH1072"]]
pData(bis_1072) <- data.frame(replicate=gsub(".*-", "", colnames(bis_1072)),
tissue=substr(colnames(bis_1072), 1, nchar(colnames(bis_1072))-3),
row.names=colnames(bis_1072))
colData(bis_1072)$tissue <- gsub("-", "_", colData(bis_1072)$tissue)
bis_1072 <- renameSeqlevels(bis_1072, mapSeqlevels(seqlevels(bis_1072), "UCSC"))
bis_1072 <- bis_1072[seqnames(bis_1072)=="chr19",]
bis_1072 <- bis_1072[240201:240300,]
tissue <- factor(pData(bis_1072)$tissue)
tissue <- relevel(tissue, "Liver_Treg")
design <- model.matrix(~tissue)
colnames(design) <- gsub("tissue", "", colnames(design))
colnames(design)[1] <- "Intercept"
rownames(design) <- colnames(bis_1072)
methdesign <- edgeR::modelMatrixMeth(design)
cont.mat <- limma::makeContrasts(treg_vs_tcon=Lymph_N_Treg-Lymph_N_Tcon,
fat_vs_ln=Fat_Treg-Lymph_N_Treg,
skin_vs_ln=Skin_Treg-Lymph_N_Treg,
fat_vs_skin=Fat_Treg-Skin_Treg,
levels=methdesign)
seq_annot <- sequencing.annotate(bis_1072, methdesign, all.cov = TRUE,
contrasts = TRUE, cont.matrix = cont.mat,
coef = "treg_vs_tcon", fdr=0.05)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.