getModulePreservation() examines module replication between a
reference and a test data set by estimating various preservation statistics,
which are then saved as a .txt file. Correlations are performed using either
pearson or bicor methods.
Usage
getModulePreservation(
  meth_disc,
  regions_disc,
  meth_rep,
  regions_rep,
  corType = c("pearson", "bicor"),
  maxPOutliers = 0.1,
  nPermutations = 100,
  save = TRUE,
  file = "Module_Preservation_Stats.txt",
  verbose = TRUE
)Arguments
- meth_disc
 A
numeric matrix, where each row is a sample and each column is a region. This is typically obtained fromadjustRegionMeth()and is related to the discovery (reference) data set.- regions_disc
 A
data.frameof regions with module assignments, which is typically obtained fromgetModules()and is related to the discovery data set.- meth_rep
 A
numeric matrix, where each row is a sample and each column is a region. This is typically obtained fromadjustRegionMeth()and is related to the replication (test) data set.- regions_rep
 A
data.frameof regions with module assignments, which is typically obtained fromgetModules()and is related to the replication data set.- corType
 A
character(1)indicating which correlation statistic to use in the adjacency calculation.- maxPOutliers
 A
numeric(1)specifying the maximum percentile that can be considered outliers on each side of the median for thebicorstatistic.- nPermutations
 A
numeric(1)indicating the number of permutations to perform in the permutation test.- save
 A
logical(1)indicating whether to save thedata.frame.- file
 A
character(1)giving the file name (.txt) for the saveddata.frame.- verbose
 A
logical(1)indicating whether messages should be printed.
Details
Identical sets of regions should be assessed and assigned modules within
discovery (reference) and replication (test) data sets, though the replication
regions may be a subset of the discovery regions due to low coverage. It's
also recommended to filter CpGs so identical loci are also assessed within
regions. Network parameters should be as similar as possible, although
modules should be identified independently between the discovery and
replication datasets. Preservation statistics are calculated by
WGCNA::modulePreservation(), with corFnc set to either cor or
bicor. More information is given in the documentation for
WGCNA::modulePreservation().
See also
getModules()to build a comethylation network and identify modules of comethylated regions.plotModulePreservation()to visualize module preservation statistics.
Examples
if (FALSE) {
# Calculate Module Preservation
regions_disc <- modules_disc$regions
regions_rep <- modules_rep$regions
preservation <- getModulePreservation(methAdj_disc,
                                      regions_disc = regions_disc,
                                      meth_rep = methAdj_rep,
                                      regions_rep = regions_rep,
                                      corType = "pearson",
                                      file = "Module_Preservation_Stats.txt")
# Visualize Module Preservation
plotModulePreservation(preservation, file = "Module_Preservation_Plots.pdf")
}