Skip to contents

The iSEEfgseaResults class is used to provide an common interface to pathway analysis results produced by the fgsea package. It provides methods to access the set of features in each pathway.

Value

iSEEfgseaResults() returns an object of class iSEEfgseaResults.

Slot overview

This class inherits all its slots directly from its parent class iSEEpathwaysResults.

Constructor

iSEEfgseaResults(data, pathwayType, pathwaysList = NULL, featuresStats = NULL) creates an instance of a iSEEfgseaResults class, with:

data

A data.frame produced by fgsea::fgsea().

pathwayType

A character scalar specifying the type of pathway (e.g., "GO"). See embedPathwaysResults.

pathwaysList

A named list of pathways and associated feature identifiers.

featuresStats

Feature-level statistics used in the pathway analysis.

Supported methods

In the following code snippets, x is an instance of a iSEEfgseaResults class. Refer to the documentation for each method for more details on the remaining arguments.

  • embedPathwaysResults(x, se, name, pathwayType, ...) embeds x in se under the identifier name. See embedPathwaysResults() for more details.

  • pathwaysList(x) returns the named list of pathway identifiers and associated

Author

Kevin Rue-Albrecht

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_4 0.08684211 0.2857143 0.24891111 -0.3611377 -1.3079994    77
#> 2: pathway_3 0.11428571 0.2857143 0.20429476 -0.4161521 -1.3133662    90
#> 3: pathway_2 0.49532710 0.7552083 0.08809450 -0.2666742 -0.9741345    37
#> 4: pathway_5 0.60416667 0.7552083 0.08312913 -0.2459574 -0.9332602    61
#> 5: pathway_1 0.83769634 0.8376963 0.06658921 -0.2175310 -0.8140116    59
#>     leadingEdge
#>          <list>
#> 1: feature_....
#> 2: feature_....
#> 3: feature_....
#> 4: feature_....
#> 5: feature_....

##
# iSEEfgseaResults ----
##

# Embed the FGSEA results in the SummarizedExperiment object
se <- embedPathwaysResults(fgseaRes, se, name = "fgsea", class = "fgsea",
  pathwayType = "simulated", pathwaysList = pathways_list, featuresStats = features_stats)
se
#> class: SummarizedExperiment 
#> dim: 100 8 
#> metadata(1): iSEEpathways
#> assays(1): counts
#> rownames(100): feature_1 feature_2 ... feature_99 feature_100
#> rowData names(0):
#> colnames(8): sample_1 sample_2 ... sample_7 sample_8
#> colData names(0):

##
# Access ----
##

pathwaysResultsNames(se)
#> [1] "fgsea"
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.8376963  0.837696 0.0665892 -0.217531 -0.814012        59
#> pathway_2 0.4953271  0.755208 0.0880945 -0.266674 -0.974134        37
#> pathway_3 0.1142857  0.285714 0.2042948 -0.416152 -1.313366        90
#> pathway_4 0.0868421  0.285714 0.2489111 -0.361138 -1.307999        77
#> pathway_5 0.6041667  0.755208 0.0831291 -0.245957 -0.933260        61
#>                                    leadingEdge
#>                                         <list>
#> pathway_1  feature_84,feature_10,feature_4,...
#> pathway_2 feature_42,feature_36,feature_66,...
#> pathway_3 feature_42,feature_36,feature_84,...
#> pathway_4 feature_42,feature_36,feature_84,...
#> pathway_5 feature_42,feature_84,feature_10,...
#> 
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.8376963  0.837696 0.0665892 -0.217531 -0.814012        59
#> pathway_2 0.4953271  0.755208 0.0880945 -0.266674 -0.974134        37
#> pathway_3 0.1142857  0.285714 0.2042948 -0.416152 -1.313366        90
#> pathway_4 0.0868421  0.285714 0.2489111 -0.361138 -1.307999        77
#> pathway_5 0.6041667  0.755208 0.0831291 -0.245957 -0.933260        61
#>                                    leadingEdge
#>                                         <list>
#> pathway_1  feature_84,feature_10,feature_4,...
#> pathway_2 feature_42,feature_36,feature_66,...
#> pathway_3 feature_42,feature_36,feature_84,...
#> pathway_4 feature_42,feature_36,feature_84,...
#> pathway_5 feature_42,feature_84,feature_10,...

pathwayType(pathwaysResults(se, "fgsea"))
#> [1] "simulated"
head(lengths(pathwaysList(pathwaysResults(se, "fgsea"))))
#> pathway_1 pathway_2 pathway_3 pathway_4 pathway_5 
#>        59        37        90        77        61 
head(featuresStats(pathwaysResults(se, "fgsea")))
#>   feature_1   feature_2   feature_3   feature_4   feature_5   feature_6 
#>  0.23742535 -1.31281425  0.74702859 -1.56251843  0.07105336 -0.63953477