View source: R/Build_Comethylation_Network.R
| getSoftPower | R Documentation |
getSoftPower() analyzes scale-free topology to estimate the best
soft-thresholding power from a vector of powers, calculate fit indices, and
then saves this as a .rds file. Possible correlation statistics include
pearson and bicor.
getSoftPower(
meth,
powerVector = 1:20,
corType = c("pearson", "bicor"),
maxPOutliers = 0.1,
RsquaredCut = 0.8,
blockSize = 40000,
gcInterval = blockSize - 1,
save = TRUE,
file = "Soft_Power.rds",
verbose = TRUE
)
meth |
A |
powerVector |
A |
corType |
A |
maxPOutliers |
A |
RsquaredCut |
A |
blockSize |
A |
gcInterval |
A |
save |
A |
file |
A |
verbose |
A |
Soft power is estimated by WGCNA::pickSoftThreshold(), with corFnc
set to either cor or bicor. Calculations are performed for a
signed network in blocks of regions of size blockSize (default = 40000).
The best soft power threshold is chosen as the lowest power where fit
(R-squared) is greater than RsquaredCut (default = 0.8). More
information is given in the documentation for WGCNA::pickSoftThreshold().
A list with two elements: powerEstimate, which gives the
estimated best soft-thresholding power, and fitIndices, which
is a data.frame with statistics on scale-free topology,
including fit and connectivity, along with network density,
centralization, and heterogeneity.
getRegionMeth(), getPCs(), and adjustRegionMeth() to
extract methylation data and then adjust it for the top
principal components.
plotSoftPower() to visualize fit and connectivity for soft
power estimation.
getModules() to build a comethylation network and identify
modules of comethylated regions.
## Not run:
# 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")
# 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")
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.