Adjust Methylation Data for Principal Components
Source:R/Examine_Region_Methylation.R
adjustRegionMeth.Rd
adjustRegionMeth()
adjusts region methylation data for the top
principal components, transposes it, and then saves it as a .rds file.
Usage
adjustRegionMeth(
meth,
PCs,
save = TRUE,
file = "Adjusted_Region_Methylation.rds",
verbose = TRUE
)
Arguments
- meth
A
numeric matrix
, where each row is a region and each column is a sample. This is typically obtained fromgetRegionMeth()
.- PCs
A
numeric matrix
, where each row is a sample and each column is a principal component. This is typically obtained fromgetPCs()
.- save
A
logical(1)
indicating whether to save thematrix
.- file
A
character(1)
giving the file name (.rds) for the savedmatrix
.- verbose
A
logical(1)
indicating whether messages should be printed.
Details
adjustRegionMeth()
regresses out the top principal components
generated by getPCs()
. This is the same approach as taken by
sva::sva_network()
. More information on the function and approach is given
in the documentation and publications related to the sva package.
See also
getRegionMeth()
to extract region methylation values.getPCs()
to calculate top principal components for region methylation values.getMEtraitCor()
to compare these top PCs to sample traits.getDendro()
andplotDendro()
to generate and visualize dendrograms.getSoftPower()
andplotSoftPower()
to estimate the best soft-thresholding power and visualize scale-free topology fit and connectivity.getModules()
to build a comethylation network and identify modules of comethylated regions.
Examples
if (FALSE) {
# Get Methylation Data
meth <- getRegionMeth(regions, bs = bs, file = "Region_Methylation.rds")
# Adjust Methylation Data for PCs
mod <- model.matrix(~1, data = pData(bs))
PCs <- getPCs(meth, mod = mod, file = "Top_Principal_Components.rds")
methAdj <- adjustRegionMeth(meth, PCs = PCs,
file = "Adjusted_Region_Methylation.rds")
# Compare Top PCs to Sample Traits
MEtraitCor <- getMEtraitCor(PCs, colData = colData, corType = "bicor",
file = "PC_Trait_Correlation_Stats.txt")
PCdendro <- getDendro(PCs, distance = "bicor")
PCtraitDendro <- getCor(PCs, y = colData, corType = "bicor", robustY = FALSE) %>%
getDendro(transpose = TRUE)
plotMEtraitCor(PCtraitCor, moduleOrder = PCdendro$order,
traitOrder = PCtraitDendro$order,
file = "PC_Trait_Correlation_Heatmap.pdf")
# Assess Sample Similarity
getDendro(methAdj, distance = "euclidean") %>%
plotDendro(file = "Sample_Dendrogram.pdf", expandY = c(0.25,0.08))
# Select Soft Power Threshold
sft <- getSoftPower(methAdj, corType = "pearson", file = "Soft_Power.rds")
plotSoftPower(sft, file = "Soft_Power_Plots.pdf")
# Get Comethylation Modules
modules <- getModules(methAdj, power = sft$powerEstimate, regions = regions,
corType = "pearson", file = "Modules.rds")
}