View source: R/Visualize_Comethylation_Network.R
getCor | R Documentation |
getCor()
calculates correlation coefficients using either
pearson
or bicor
methods. Calculations can be done between
columns of a single matrix or between two vectors or matrices.
getCor(
x,
y = NULL,
transpose = FALSE,
corType = c("bicor", "pearson"),
maxPOutliers = 0.1,
robustY = TRUE,
verbose = TRUE
)
x |
A |
y |
A |
transpose |
A |
corType |
A |
maxPOutliers |
A |
robustY |
A |
verbose |
A |
The first input argument can be optionally transposed. The correlation
calculations are performed by WGCNA::cor()
and WGCNA::bicor()
.
A numeric matrix
.
getModules()
to build a comethylation network and identify
modules of comethylated regions.
getDendro()
and plotDendro()
to generate and visualize
dendrograms.
plotHeatmap()
to visualize correlations between samples and
modules.
getMEtraitCor()
to calculate pairwise correlation
coefficients and p-values between module eigennode values.
## Not run:
# Get Comethylation Modules
modules <- getModules(methAdj, power = sft$powerEstimate, regions = regions,
corType = "pearson", file = "Modules.rds")
# Examine Correlations between Modules
MEs <- modules$MEs
moduleDendro <- getDendro(MEs, distance = "bicor")
plotDendro(moduleDendro, labelSize = 4, nBreaks = 5,
file = "Module_ME_Dendrogram.pdf")
moduleCor <- getCor(MEs, corType = "bicor")
plotHeatmap(moduleCor, rowDendro = moduleDendro, colDendro = moduleDendro,
file = "Module_Correlation_Heatmap.pdf")
moduleCorStats <- getMEtraitCor(MEs, colData = MEs, corType = "bicor",
robustY = TRUE,
file = "Module_Correlation_Stats.txt")
# Examine Correlations between Samples
sampleDendro <- getDendro(MEs, transpose = TRUE, distance = "bicor")
plotDendro(sampleDendro, labelSize = 3, nBreaks = 5,
file = "Sample_ME_Dendrogram.pdf")
sampleCor <- getCor(MEs, transpose = TRUE, corType = "bicor")
plotHeatmap(sampleCor, rowDendro = sampleDendro, colDendro = sampleDendro,
file = "Sample_Correlation_Heatmap.pdf")
# Visualize Module Eigennode Values
plotHeatmap(MEs, rowDendro = sampleDendro, colDendro = moduleDendro,
legend.title = "Module\nEigennode",
legend.position = c(0.37,0.89),
file = "Sample_ME_Heatmap.pdf")
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.