Apply dimensionality reduction using Weighted Gene Correlation Network Analysis
Source:R/modules.R
identifyModules.Rd
Performs Weighted gene correlation network analysis (WGCNA) and packages both the input data and subsequent results into a ModularExperiment. Calls runWGCNA to perform the analysis; see its documentation page for more information on the ICA method, parameters and outputs.
Arguments
- X
Either a SummarizedExperiment object or a matrix containing data to be subject to WGCNA.
X
should have rows as features and columns as samples.- power
An integer representing the soft-thresholding power to be used to define modules. See the assessSoftThreshold function for aid in determining this parameter.
- center_X
If
TRUE
,X
is centered (i.e., features / rows are transformed to have a mean of 0) prior to WGCNA.- scale_X
If
TRUE
,X
is scaled (i.e., features / rows are transformed to have a standard deviation of 1) before WGCNA.- assay_name
If
X
is a SummarizedExperiment, then this should be the name of the assay to be subject to WGCNA.- ...
Additional arguments to be passed to runWGCNA.
Value
A ModularExperiment is returned
containing the input data (i.e., the original data matrix in addition to
other slots if a SummarizedExperiment was used
as input). Additionally contains the results of module analysis, stored in
the reduced
and assignments
slots. The center_X
, scale_X
,
loadings
, threshold
and dendrogram
slots may also be filled depending
on the arguments given to identifyModules
.
Examples
# Get the airway data as a SummarizedExperiment (with a subset of features)
set.seed(2)
airway_se <- ReducedExperiment:::.getAirwayData(n_features = 500)
# Select soft-thresholding power to use (use capture.output to hide WGCNA's prints)
WGCNA::disableWGCNAThreads()
invisible(capture.output(fit_indices <- assessSoftThreshold(airway_se)))
estimated_power <- fit_indices$Power[fit_indices$estimated_power]
# Identify modules using WGCNA
airway_me <- identifyModules(airway_se, verbose = 0, power = estimated_power)
airway_me
#> class: ModularExperiment
#> dim: 500 8 6
#> metadata(1): ''
#> assays(3): counts normal transformed
#> rownames(500): ENSG00000186523 ENSG00000203871 ... ENSG00000147439
#> ENSG00000116661
#> rowData names(10): gene_id gene_name ... seq_coord_system symbol
#> colnames(8): SRR1039508 SRR1039509 ... SRR1039520 SRR1039521
#> colData names(9): SampleName cell ... Sample BioSample
#> 6 components