Extract contrast results embedded in a SummarizedExperiment object
Source:R/utils-SummarizedExperiment.R
pathwaysResults.Rd
pathwaysResults
returns either all pathway analysis results stored in object
or a single pathway analysis result by name.
pathwaysResultsNames
returns the names of pathway analysis results embedded in object
.
Arguments
- object
A SummarizedExperiment object.
- name
(Optional) Name of a single pathway analysis result name to extract. Use
pathwaysResultsNames(object)
to list available names.
Value
For pathwaysResultsNames
: the names of embedded contrast results available.
For pathwaysResults
: a list
of differential expression statistics.
If name
is missing, pathwaysResults
returns a list
in which each item contains the results of a single pathway analysis.
If name
is given, pathwaysResults
returns a DataFrame
that contains the results of a single pathway analysis.
Examples
library("fgsea")
##
# Simulate example data
##
simulated_data <- simulateExampleData(n_pathways = 5, n_features = 100, pathway_sizes = 15:100)
pathways_list <- simulated_data$pathwaysList
features_stats <- simulated_data$featuresStat
se <- simulated_data$summarizedexperiment
##
# Run pathway analysis ----
##
set.seed(42)
fgseaRes <- fgsea(pathways = pathways_list,
stats = features_stats,
minSize = 15,
maxSize = 500)
head(fgseaRes[order(pval), ])
#> pathway pval padj log2err ES NES size
#> <char> <num> <num> <num> <num> <num> <int>
#> 1: pathway_2 0.1629213 0.8146067 0.18470647 0.3232076 1.2078961 39
#> 2: pathway_1 0.5502249 0.9067055 0.05998925 -0.2922706 -0.9574770 79
#> 3: pathway_4 0.6109375 0.9067055 0.05700595 -0.2849186 -0.9191748 32
#> 4: pathway_3 0.7906977 0.9067055 0.04477489 -0.2704209 -0.7902678 88
#> 5: pathway_5 0.9067055 0.9067055 0.06799226 0.1875875 0.7315551 63
#> leadingEdge
#> <list>
#> 1: feature_....
#> 2: feature_....
#> 3: feature_....
#> 4: feature_....
#> 5: feature_....
##
# iSEEfgseaResults ----
##
se <- embedPathwaysResults(fgseaRes, se, name = "fgsea", class = "fgsea", pathwayType = "GO",
pathwaysList = pathways_list, featuresStats = features_stats)
##
# List result names ---
##
pathwaysResultsNames(se)
#> [1] "fgsea"
##
# Extract results ---
##
pathwaysResults(se)
#> $fgsea
#> iSEEfgseaResults with 5 rows and 7 columns
#> pval padj log2err ES NES size
#> <numeric> <numeric> <numeric> <numeric> <numeric> <integer>
#> pathway_1 0.550225 0.906706 0.0599892 -0.292271 -0.957477 79
#> pathway_2 0.162921 0.814607 0.1847065 0.323208 1.207896 39
#> pathway_3 0.790698 0.906706 0.0447749 -0.270421 -0.790268 88
#> pathway_4 0.610938 0.906706 0.0570060 -0.284919 -0.919175 32
#> pathway_5 0.906706 0.906706 0.0679923 0.187588 0.731555 63
#> leadingEdge
#> <list>
#> pathway_1 feature_73,feature_22,feature_4,...
#> pathway_2 feature_83,feature_80,feature_57,...
#> pathway_3 feature_73,feature_75,feature_36,...
#> pathway_4 feature_22,feature_11,feature_18,...
#> pathway_5 feature_83,feature_80,feature_23,...
#>
pathwaysResults(se, "fgsea")
#> iSEEfgseaResults with 5 rows and 7 columns
#> pval padj log2err ES NES size
#> <numeric> <numeric> <numeric> <numeric> <numeric> <integer>
#> pathway_1 0.550225 0.906706 0.0599892 -0.292271 -0.957477 79
#> pathway_2 0.162921 0.814607 0.1847065 0.323208 1.207896 39
#> pathway_3 0.790698 0.906706 0.0447749 -0.270421 -0.790268 88
#> pathway_4 0.610938 0.906706 0.0570060 -0.284919 -0.919175 32
#> pathway_5 0.906706 0.906706 0.0679923 0.187588 0.731555 63
#> leadingEdge
#> <list>
#> pathway_1 feature_73,feature_22,feature_4,...
#> pathway_2 feature_83,feature_80,feature_57,...
#> pathway_3 feature_73,feature_75,feature_36,...
#> pathway_4 feature_22,feature_11,feature_18,...
#> pathway_5 feature_83,feature_80,feature_23,...