Title: | Visualize Intensities Produced by Illumina's Human Methylation 450k BeadChip |
---|---|
Description: | The skewr package is a tool for visualizing the output of the Illumina Human Methylation 450k BeadChip to aid in quality control. It creates a panel of nine plots. Six of the plots represent the density of either the methylated intensity or the unmethylated intensity given by one of three subsets of the 485,577 total probes. These subsets include Type I-red, Type I-green, and Type II.The remaining three distributions give the density of the Beta-values for these same three subsets. Each of the nine plots optionally displays the distributions of the "rs" SNP probes and the probes associated with imprinted genes as series of 'tick' marks located above the x-axis. |
Authors: | Ryan Putney [cre, aut], Steven Eschrich [aut], Anders Berglund [aut] |
Maintainer: | Ryan Putney <[email protected]> |
License: | GPL-2 |
Version: | 1.39.0 |
Built: | 2024-10-31 05:34:40 UTC |
Source: | https://github.com/bioc/skewr |
A convenience function for retrieving simple barcodes from idat file names.
getBarcodes(path = getwd(), recurse = FALSE)
getBarcodes(path = getwd(), recurse = FALSE)
path |
The path or a character vector to the directory or directories in which to find the idat files. |
recurse |
logical; should the function check subdirectories to derive barcodes from any found idat files. The default is |
Barcodes will be generated by all found idats in path(s). The default path is the current working directory.
A character vector of barcodes.
Ryan Putney [email protected]
if(require(minfiData)){ path <- system.file("extdata/5723646052", package="minfiData") barcodes <- getBarcodes(path = path) }
if(require(minfiData)){ path <- system.file("extdata/5723646052", package="minfiData") barcodes <- getBarcodes(path = path) }
MethyLumiSet
object
This a wrapper function for methylumIDAT
that does not require a vector of barcodes to be provided.
getMethyLumiSet(path = getwd(), barcodes = NULL, norm = c("none", "illumina", "SWAN", "dasen"), bg.corr = TRUE)
getMethyLumiSet(path = getwd(), barcodes = NULL, norm = c("none", "illumina", "SWAN", "dasen"), bg.corr = TRUE)
path |
The path to the directory containing the idat files. |
barcodes |
A vector of barcodes specifying which idat's to read. |
norm |
Should normalization be done on the resulting |
bg.corr |
logical; if |
If only path
is provided, all idat's found in the given folder will be pulled. If only barcodes
is given, corresponding idat's will be pulled from the current working directory. Both path
and barcodes
may be passed for finer control. The default is to pull all idat's found in the current working directory.
A MethyLumiSet
object
One would probably not normally want to use the preprocess option at this stage. It is more likely that a MethyLumiSet
of the raw data will be desired. Then the preprocess
method may be used to normalize the raw data or use background subtraction only on the raw data. See the vignette for example workflow.
Ryan Putney [email protected]
Davis S, Du P, Bilke S, Triche T, Jr. and Bootwalla M (2014). methylumi: Handle Illumina methylation data. R package version 2.12.0.
if(require('minfiData')) { path <- system.file("extdata/5723646052", package="minfiData") barcodes <- getBarcodes(path = path) methylumiset.raw <- getMethyLumiSet(path = path, barcodes = barcodes[1]) }
if(require('minfiData')) { path <- system.file("extdata/5723646052", package="minfiData") barcodes <- getBarcodes(path = path) methylumiset.raw <- getMethyLumiSet(path = path, barcodes = barcodes[1]) }
Utilizes smsn.mix
from the mixsmsn
package to find the parameters for a finite mixture of skew normal distributions to model the overall distribution of signal intensities for a subset of probes on the Illumina Infinium HumanMethylation450. The probes may be subset by type and methylated or unmethylated. It can also be specified whether the SNP(rs), imprinted(idmr), or ch probes should be included or filtered out prior to parameter estimation.
getSNparams(MethyLumiSet, allele = c('M', 'U'), type = c('I-red', 'I-green', 'II'), snps = TRUE, idmr = TRUE, ch = FALSE)
getSNparams(MethyLumiSet, allele = c('M', 'U'), type = c('I-red', 'I-green', 'II'), snps = TRUE, idmr = TRUE, ch = FALSE)
MethyLumiSet |
A |
allele |
Should parameter estimation be done on the methylated or unmethylated signal intensities |
type |
Use the signal intensities for which probe type |
snps |
logical; should the rs probes be included in the dataset. The default is |
idmr |
logical; should the probes of imprinted gene loci be included in the dataset. The default is |
ch |
logical; should the ch probes be included in the dataset. The default is |
A Skew.normal
object as returned by smsn.mix
from the mixsmsn
package with the means and modes of the components added.
Ryan Putney [email protected]
Pidsley R, Wong CCY, Volta M, Lunnon K, Mill J, Schalwyk LC(2013). A data-driven approach to preprocessing Illumina 450k methylation array data. BMC Genomics, 14:293.
Prates MO, Cabral CRB, Lachos VH (2013).mixsmsn: Fitting Finite Mixture of Scale Mixture of Skew-Normal Distributions. Journal of Statistical Software, 54(12), 1-20. http://www.jstatsoft.org/v54/i12/
if(require('wateRmelon')) { data(melon) mixes.raw.meth.II <- getSNparams(melon[,1], 'M', 'II') }
if(require('wateRmelon')) { data(melon) mixes.raw.meth.II <- getSNparams(melon[,1], 'M', 'II') }
Creates a panel of nine plots. Six of the plots represent the density of either the methylated intensity or the unmethylated intensity given by one of three subsets of the 485,577 total probes. These subsets include Type I-red, Type I-green, and Type II.The remaining three distributions give the density of the beta-values for these same three subsets. Each of the nine plots optionally displays the distributions of the "rs" SNP probes and the probes associated with imprinted genes(Pidsley,2013) as a series of 'tick' marks located above the x-axis.
panelPlots(MethyLumiSet, typeIRedModels, typeIGreenModels, typeIIModels, plot = c("panel", "frames"), samp.num = NULL, frame.nums = 1:9, norm = "", idmr = TRUE, snps = TRUE)
panelPlots(MethyLumiSet, typeIRedModels, typeIGreenModels, typeIIModels, plot = c("panel", "frames"), samp.num = NULL, frame.nums = 1:9, norm = "", idmr = TRUE, snps = TRUE)
MethyLumiSet |
The |
typeIRedModels |
A |
typeIGreenModels |
A |
typeIIModels |
A |
plot |
Should the output consist of panel plots–one panel per sample or a single panel if |
samp.num |
If plotting for a single sample is desired, for which sample. The number given simply refers to the |
frame.nums |
if |
norm |
A character string which will be displayed as part of the main title for each plot. Useful in indicated which normalization method was used for the modeled and plotted data |
idmr |
logical; should the intensities of the idmr probes be plotted as a series of tick-marks above the x-axis. The default is |
snps |
logical; should the intensities of the rs probes be plotted as a series of tick-marks above the x-axis. The default is |
No return value. Only plots are generated.
Please refer to the vignette for an example workflow.
Ryan Putney [email protected]
Prates MO, Cabral CRB, Lachos VH (2013).mixsmsn: Fitting Finite Mixture of Scale Mixture of Skew-Normal Distributions. Journal of Statistical Software, 54(12), 1-20. http://www.jstatsoft.org/v54/i12/
if(require('minfiData')) { path <- system.file("extdata/5723646052", package="minfiData") methylumiset.raw <- getMethyLumiSet(path = path) mixes.raw.meth.I.red <- getSNparams(methylumiset.raw, 'M', 'I-red') mixes.raw.meth.I.green <- getSNparams(methylumiset.raw, 'M', 'I-green') mixes.raw.meth.II <- getSNparams(methylumiset.raw, 'M', 'II') mixes.raw.unmeth.I.red <- getSNparams(methylumiset.raw, 'U', 'I-red') mixes.raw.unmeth.I.green <- getSNparams(methylumiset.raw, 'U', 'I-green') mixes.raw.unmeth.II <- getSNparams(methylumiset.raw, 'U', 'II') mixes.I.red <- list(mixes.raw.meth.I.red, mixes.raw.unmeth.I.red) mixes.I.green <- list(mixes.raw.meth.I.green, mixes.raw.unmeth.I.green) mixes.II <- list(mixes.raw.meth.II, mixes.raw.unmeth.II) panelPlots(methylumiset.raw, mixes.I.red, mixes.I.green, mixes.II) }
if(require('minfiData')) { path <- system.file("extdata/5723646052", package="minfiData") methylumiset.raw <- getMethyLumiSet(path = path) mixes.raw.meth.I.red <- getSNparams(methylumiset.raw, 'M', 'I-red') mixes.raw.meth.I.green <- getSNparams(methylumiset.raw, 'M', 'I-green') mixes.raw.meth.II <- getSNparams(methylumiset.raw, 'M', 'II') mixes.raw.unmeth.I.red <- getSNparams(methylumiset.raw, 'U', 'I-red') mixes.raw.unmeth.I.green <- getSNparams(methylumiset.raw, 'U', 'I-green') mixes.raw.unmeth.II <- getSNparams(methylumiset.raw, 'U', 'II') mixes.I.red <- list(mixes.raw.meth.I.red, mixes.raw.unmeth.I.red) mixes.I.green <- list(mixes.raw.meth.I.green, mixes.raw.unmeth.I.green) mixes.II <- list(mixes.raw.meth.II, mixes.raw.unmeth.II) panelPlots(methylumiset.raw, mixes.I.red, mixes.I.green, mixes.II) }
This is a wrapper function that allows normalizing of a MethyLumiSet using either a BeadStudio approximation, SWAN, or dasen. If desired, background correction only may be performed on the raw data.
preprocess(MethyLumiSet, norm = c("none", "illumina", "SWAN", "dasen"), bg.corr = TRUE)
preprocess(MethyLumiSet, norm = c("none", "illumina", "SWAN", "dasen"), bg.corr = TRUE)
MethyLumiSet |
A |
norm |
The normalization method to be used |
bg.corr |
If TRUE, background subtarction using negative controls is performed. Ignored unless |
Both Illumina style normalization via controls and the background correct method are handled by methylumi
. The SWAN
and dasen
normalization methods are both performed by wateRmelon
A MethyLumiSet
Ryan Putney [email protected]
Davis S, Du P, Bilke S, Triche T, Jr. and Bootwalla M (2014). methylumi: Handle Illumina methylation data. R package version 2.12.0.
Maksimovic J, Gordon L, Oshlack A (2012). SWAN: Subset Quantile Within-Array Normalization for Illumina Infinium HumanMethylation450 BeadChips. Genome Biology, 13:R44.
Pidsley R, Wong CCY, Volta M, Lunnon K, Mill J, Schalwyk LC(2013). A data-driven approach to preprocessing Illumina 450k methylation array data. BMC Genomics, 14:293.
Schalkwyk LC, Pidsley R, Wong CC, Touleimat N, Defrance M, Teschendorff A and Maksimovic J (2013). wateRmelon: Illumina 450 methylation array normalization and metrics. R package version 1.5.1.
if(require('wateRmelon')) { data(melon) melon.dasen <- preprocess(melon, norm = 'dasen') }
if(require('wateRmelon')) { data(melon) melon.dasen <- preprocess(melon, norm = 'dasen') }
Thus function accepts a MethyLumiSet
object generated by methylumi
or a MethylSet
object generated by minfi
. It will subset the probes by type–"I-red"
, "I-green"
, or "II"
–and return a matrix of the methylated, "M"
, or unmethylated, "U"
signal intensities. It is also possible to include or filter out probes according to whether they are CpG sites(cg), SNPs(rs), imprinted(idmr) gene sites, or non-CpG loci(ch).
subsetProbes(object, allele = c("M", "U"), type = c("I-red", "I-green", "II"), cg = TRUE, snps = TRUE, idmr = TRUE, ch = FALSE)
subsetProbes(object, allele = c("M", "U"), type = c("I-red", "I-green", "II"), cg = TRUE, snps = TRUE, idmr = TRUE, ch = FALSE)
object |
A |
allele |
Should methylated or unmethylated data for the probes be returned. |
type |
May be |
cg |
Logical; Should the returned dataset contain the CpG probes. The default is |
snps |
Logical; Should the returned dataset conain the rs probes. The default is |
idmr |
Logical; should the returned dataset include probes that interrogate imprinted gene sites as given by Pidsley et al.(2013). The default is |
ch |
Logical; should the returned dataset include the non-CpG (ch) probes. The default if |
A matrix
Ryan Putney [email protected]
Pidsley R, Wong CCY, Volta M, Lunnon K, Mill J, Schalwyk LC(2013). A data-driven approach to preprocessing Illumina 450k methylation array data. BMC Genomics, 14:293.
if(require('wateRmelon')) { data(melon) melon.meth.II <- subsetProbes(melon, 'M', 'II') }
if(require('wateRmelon')) { data(melon) melon.meth.II <- subsetProbes(melon, 'M', 'II') }