Package 'cpvSNP'

Title: Gene set analysis methods for SNP association p-values that lie in genes in given gene sets
Description: Gene set analysis methods exist to combine SNP-level association p-values into gene sets, calculating a single association p-value for each gene set. This package implements two such methods that require only the calculated SNP p-values, the gene set(s) of interest, and a correlation matrix (if desired). One method (GLOSSI) requires independent SNPs and the other (VEGAS) can take into account correlation (LD) among the SNPs. Built-in plotting functions are available to help users visualize results.
Authors: Caitlin McHugh, Jessica Larson, and Jason Hackney
Maintainer: Caitlin McHugh <[email protected]>
License: Artistic-2.0
Version: 1.37.0
Built: 2024-07-19 10:45:56 UTC
Source: https://github.com/bioc/cpvSNP

Help Index


Create a Density Plot of P-Values, Highlighting SNPs Within a Gene Set

Description

Plots association test p-values, highlighting p-values within gene sets.

Usage

assocPvalBySetPlot(pvals, set, title = "", geneSetColor = "red",
   xlab = expression(-log[10]~p-value), 
   ylab = "Density", 
   xlim = NULL, ylim = NULL,
  ...)

Arguments

pvals

A numerical vector of p-values with names corresponding to elements listed in the geneIds slot in the set object.

set

An object of class GeneSet where geneIds holds the items of each set corresponding to the pvals.

title

A character string containing the title of the plot. Default is "", a blank title.

geneSetColor

A character vector holding the color the highlighted SNPs in the set object should be. Default is red.

xlab

A character vector with the label for the x-axis on the plot. Default is -log10(PValue).

ylab

A character vector holding the label for the y-axis on the plot. Default is Density.

xlim

A vector of two numeric values with the limits of plotting on the x-axis. By default, the range of the x-axis values will be determined from the data.

ylim

A vector with two numeric values used as limits for plotting on the y-axis. By default, the range of the y-axis values will be determined from the data.

...

Further arguments to be passed to the plotting methods, such as graphical parameters.

Details

Creates a density plot of all -log10(p-value), and overlaying a density plot of -log10(p-value) for SNPs within a gene set of interest.

Value

Creates a plot.

Author(s)

Jason Hackney, Jessica Larson, Caitlin McHugh [email protected]


Create a GRanges Object for a GWAS data.frame

Description

Creates a GRanges object used for SNP set analysis.

Usage

createArrayData(arrayData, positionName = NULL,
  chromosomeName = "chromosome", chromosomeNameConvention = "NCBI",
  verbose = TRUE)

Arguments

arrayData

A data.frame containing array data from which the GRanges object will be created. This object is expected to include the position and chromosome, along with a SNP identifier and corresponding p-value.

positionName

The name of the column in the filepathData object that holds position data for each probe. By default, this value is NULL and there are columns Start and End which hold this information.

chromosomeName

The name of the column in the filepathData object that holds chromosome data for each probe. By default, this value is chromosome.

chromosomeNameConvention

The naming convention used for the chromosomes, either NCBI (default), UCSC, or a user specified mapping. By default, this is NCBI, indicating the autosomes are stored as integers and the mitochondrial and sex chromosomes are stored as MT and X, Y, respectively. The UCSC naming convention appends "chr" to each chromosome, and the mitochondrial coding is "chrM." Alternatively, a user can supply the mapping from each chromosome to the UCSC naming convention in the form of a character vector where the elements are the UCSC names and the corresponding names of the elements are the current chromosome names. An ‘NA’ chromosome will be mapped to unknown. For more information, please see the GenomeInfoDb BioConductor R package.

verbose

A logical argument indicating whether output should be printed. The default is FALSE.

Details

This function takes a data.frame and creates a GRanges object used for SNP set analysis.

Value

A GRanges object.

Author(s)

Jason Hackney, Jessica Larson, Caitlin McHugh [email protected]

Examples

data(geneSetAnalysis)
    head(geneSetAnalysis$arrayData)
    arrayDataGR <- createArrayData(geneSetAnalysis[["arrayData"]], positionName="Position")

~~ Methods for Function degreesOfFreedom ~~

Description

This function returns the degrees of freedom used to calculate p-values from the corresponding GeneSet or a list of degrees of freedom for a GeneSetCollection.

Usage

degreesOfFreedom(object)

Arguments

object

An object of type GeneSetResult, VEGASResult, GLOSSIResult, GeneSetResultCollection, VEGASResultCollection, or GLOSSIResultCollection

Methods

Defined methods include:

signature(object = "GeneSetResult")

Returns the degrees of freedom for a GeneSetResult object

signature(object = "VEGASResult")

Returns the degrees of freedom for a VEGASResult object

signature(object = "GLOSSIResult")

Returns the degrees of freedom for a GLOSSIResult object

signature(object = "GeneSetResultCollection")

Returns a list of degrees of freedom for a GeneSetResultCollection objects (1 for each set)

signature(object = "VEGASResultCollection")

Returns a list of degrees of freedom for a VEGASResultCollection objects (1 for each set)

signature(object = "GLOSSIResultCollection")

Returns a list of degrees of freedom for a GLOSSIResultCollection objects (1 for each set)

Value

Returns an integer or a list of integers indicating the degrees of freedom in the corresponding object.

Author(s)

Caitlin McHugh

See Also

VEGASResult-class, pValue

Examples

## Not run: 
	degreesOfFreedom( vegas(geneSet, assoc_table, ldMatrix) ) 
	
## End(Not run)

SNP Array Data to Use for Gene Set Analysis

Description

This dataset provides simulated association p-values for an example set of SNP rsIDs, including their chromosome and position along the chromosome. The data in this file are all SNPs that lie on chromosome 1.

Usage

exampleArrayData

Format

A data.frame containing 119,487 rows and the following 4 columns:

P simulated p-values
SNP SNP rsIDs
Position SNP position
chromosome SNP chromosome

Data to run gene set analysis methods

Description

A list including the following components:

arrayData – a data.frame of dimension 119487 rows and 4 columns including P-values, SNP ids, position and chromosome for a GWAS experiment;

geneSets – a GeneSetCollection containing gene Entrez ids;

ldMat – a matrix of dimension 871 x 871 holding linkage disequilibrium values;

indepSNPs – a list containing 302 independent SNP ids;

Usage

geneSetAnalysis

Examples

data(geneSetAnalysis)
    head(geneSetAnalysis$arrayData)

~~ Methods for Function geneSetName ~~

Description

This function returns the gene set name for a specified GeneSet or a list of gene set names for a GeneSetCollection.

Usage

geneSetName(object)

Arguments

object

An object of type GeneSetResult, VEGASResult, GLOSSIResult, GeneSetResultCollection, VEGASResultCollection, or GLOSSIResultCollection

Methods

Defined methods include:

signature(object = "GeneSetResult")

Returns the gene set name for a specified GeneSetResult object

signature(object = "VEGASResult")

Returns the gene set name for a specified VEGASResult object

signature(object = "GLOSSIResult")

Returns the gene set name for a specified GLOSSIResult object

signature(object = "GeneSetResultCollection")

Returns a list of gene set names for a GeneSetResultCollection object (1 for each set)

signature(object = "VEGASResultCollection")

Returns a list of gene set names for a VEGASResultCollection object (1 for each set)

signature(object = "GLOSSIResultCollection")

Returns a list of gene set names for a GLOSSIResultCollection object (1 for each set)

Value

Returns the name or a list of names of the gene sets in the corresponding object.

Author(s)

Caitlin McHugh

See Also

VEGASResult-class, pValue

Examples

## Not run: 
	geneSetName( vegas(geneSet, assoc_table, ldMatrix) ) 
	
## End(Not run)

Class "GeneSetResult"

Description

Objects of this class store results from running GeneSet methods.

Objects from the Class

Objects can be created by calls of glossi or vegas.

Slots

geneSetName:

Object of class "character", the name of the geneSet

pValue:

Object of class "numeric", the p-value

degreesOfFreedom:

Object of class "integer", the degrees of freedom used to calculate the p-value

statistic:

Object of class "numeric", the test statistic for the set

Methods

No methods defined with class "GeneSetResult" in the signature.

Accessors

In the code snippets below, object is a GeneSetResult object.

geneSetName(object):

The name of the gene set.

pValue(object):

A numeric p-value.

degreesOfFreedom(object):

Integer value of the degrees of freedom used in the statistical test.

statistic(object):

A numeric value corresponding to the statistic calculated from the test performed.

Author(s)

Caitlin McHugh

See Also

GeneSetResultCollection, GLOSSIResult

Examples

showClass("GeneSetResult")

Class "GeneSetResultCollection"

Description

Objects of this class store results from running gene set methods such as GLOSSI or VEGAS. GeneSetResultCollection objects contain a list of GeneSetResult objects (one result for each GeneSet)

Objects from the Class

Objects can be created by calls of glossi or vegas.

Methods

No methods defined with class "GeneSetResultCollection" in the signature.

Accessors

In the code snippets below, object is a GeneSetResultCollection object.

geneSetName(object):

A list of the names of the gene sets.

pValue(object):

A list of numeric p-values.

degreesOfFreedom(object):

List of integer values of the degrees of freedom used in the statistical tests performed on each gene set.

statistic(object):

A list of numeric values corresponding to the statistics calculated from the tests performed on each gene set.

Author(s)

Caitlin McHugh

See Also

GeneSetResult, GLOSSIResultCollection

Examples

showClass("GeneSetResultCollection")

Gene Set Collection object

Description

This gene set collection object provides four gene sets of varying sizes to be used as examples of how to run the cpvSNP package functions.

Usage

geneSets

Format

A GeneSetCollection object of length 4. The values stored in the geneIds slot are Entrez ids. Gene sets 1-3 were randomly selected with sizes 100, 200 and 50, respectively. The fourth GeneSet holds genes in which the simulated p-values for corresponding SNPs are significant.


Translate a List of Gene Ids to Their Corresponding SNP Ids

Description

Translates a set of gene ids to their corresponding SNPs.

Usage

geneToSNPList(geneList, arrayData, genes, maxgap = 20000, quiet = TRUE)

Arguments

geneList

An object of class GeneSetCollection or GeneSet where geneIds holds the gene Entrez ids.

arrayData

A GRanges object that corresponds to the probes for which data are present.

genes

A GRanges object holding the gene transcripts for which the mapping will occur. The ids of the genes should match the ids that are listed in the geneList object.

maxgap

An integer value indicating the flanking region around the gene for which SNPs are considered part of the gene. Default is 20000 (20kb).

quiet

A logical variable indicating whether warnings should be output from the SNP to gene mapping step.

Details

Translates a GeneSet or GeneSetCollection of gene Entrez ids to the corresponding SNPs that lie within a prespecified region of the gene.

Value

A GeneSetCollection object containing all SNP ids that lie within the genes listed in geneList.

Author(s)

Jason Hackney, Jessica Larson, Caitlin McHugh [email protected]

Examples

data(geneSetAnalysis)
    head(geneSetAnalysis$arrayData)
    arrayDataGR <- createArrayData(geneSetAnalysis[["arrayData"]], positionName="Position")
    
    library(TxDb.Hsapiens.UCSC.hg19.knownGene)
    genesHg19 <- genes(TxDb.Hsapiens.UCSC.hg19.knownGene)
	geneSets <- geneSetAnalysis[["geneSets"]]
	
	snpsGSC <- geneToSNPList(geneSets, arrayDataGR, genesHg19)

Calculate a Chi-squared Statistic and P-Value for Independent Items of a Set Using the GLOSSI Method

Description

Calculates the p-value representing the association of the set with the phenotype of interest, assuming all items in the set are independent, using the GLOSSI method.

Usage

glossi(pvals, snp.gsc)

Arguments

pvals

A numerical vector of p-values with names corresponding to elements listed in the geneIds slot in the snp.gsc object.

snp.gsc

An object of class GeneSetCollection where geneIds holds the items of each set corresponding to the pvals.

Details

This function calculates a p-value for sets of SNPs that reside within a gene set collection. We calculate the chi- square p-values and statistic by applying Fisher's transformation to the observed p-values.

Value

An object with the corresponding GLOSSI results. If snp.gsc is a GeneSetCollection (i.e., multiple sets of interest), then the corresponding GLOSSIResultCollection is returned. If snp.gsc is a GeneSet, a GLOSSIResult object will be returned.

Author(s)

Jason Hackney, Jessica Larson, Caitlin McHugh [email protected]

References

Chai, High-Seng and Sicotte, Hughes et al. GLOSSI: a method to assess the association of genetic loci-sets with complex diseases. BMC Bioinformatics, 2009.

Examples

set.seed(30)
    pvals <- runif(100)
    names(pvals) <- paste0("rs", 1:100)
    snpGS5 <- GeneSet(geneIds=names(pvals)[1:5], setName="set5")
    res <- glossi (pvals, snpGS5)

Calculate a Chi-squared Statistic and P-Value for Items in a Set

Description

A helper function to calculate the chi-squared statistic corresponding to an observed set of p-values.

Usage

glossiMarginal(pval, set)

Arguments

pval

A vector of p-values of length equal to the number of geneIds in the GeneSetCollection object.

set

An element of a GeneSetCollection object.

Details

This function calculates a chi-squared statistic from a set of pvalues.

Value

An object of type GLOSSIResult.

Author(s)

Jason Hackney, Jessica Larson

Examples

set.seed(30)
    pvals <- runif(100)
    names(pvals) <- paste0("rs", 1:100)
    snpGS5 <- GeneSet(geneIds=names(pvals)[1:5], setName="set5")
    res <- glossiMarginal (pvals, snpGS5)

Class "GLOSSIResult"

Description

Objects of this class store results from running GLOSSI.

Objects from the Class

Objects can be created by calls of glossi.

Slots

geneSetName:

Object of class "character", the name of the geneSet

pValue:

Object of class "numeric", the p-value

degreesOfFreedom:

Object of class "integer", the degrees of freedom used to calculate the p-value

statistic:

Object of class "numeric", the test statistic for the set

Methods

No methods defined with class "GLOSSIResult" in the signature.

Accessors

In the code snippets below, object is a GLOSSIResult object.

geneSetName(object):

The name of the gene set.

pValue(object):

A numeric p-value.

degreesOfFreedom(object):

Integer value of the degrees of freedom used in the statistical test.

statistic(object):

A numeric value corresponding to the statistic calculated from the test performed.

Author(s)

Caitlin McHugh

See Also

GLOSSIResultCollection, VEGASResult

Examples

showClass("GLOSSIResult")

Class "GLOSSIResultCollection"

Description

Objects of this class store results from running GLOSSI. GLOSSIResultCollection objects contain a list of GLOSSIResult objects (one result for each GeneSet)

Objects from the Class

Objects can be created by calls of glossi.

Methods

No methods defined with class "GLOSSIResult" in the signature.

Accessors

In the code snippets below, object is a GLOSSIResultCollection object.

geneSetName(object):

A list of the names of the gene sets.

pValue(object):

A list of numeric p-values.

degreesOfFreedom(object):

List of integer values of the degrees of freedom used in the statistical tests performed on each gene set.

statistic(object):

A list of numeric values corresponding to the statistics calculated from the tests performed on each gene set.

Author(s)

Caitlin McHugh

See Also

GLOSSIResult, VEGASResultCollection

Examples

showClass("GLOSSIResultCollection")

Create a Plot of P-Values Against Number of SNPs per Set for All Sets in a Collection

Description

Plots association test p-values against number of SNPs per gene set.

Usage

plotPvals(set, title = "", xlab = "SNPs per Gene Set", ylab = "PValue",
  ...)

Arguments

set

An object of class GeneSetResultCollection.

title

A character string containing the title of the plot. Default is "", a blank title.

xlab

A character vector with the label for the x-axis on the plot. Default is SNPs per Gene Set.

ylab

A character vector holding the label for the y-axis on the plot. Default is PValue.

...

Further arguments to be passed to the plotting methods, such as graphical parameters.

Details

Creates a plot of all p-values against the number of SNPs in the gene set, for all sets within a collection.

Value

Creates a plot.

Author(s)

Caitlin McHugh [email protected]


~~ Methods for Function pValue ~~

Description

This function returns the calculated p-value for a specified GeneSet or a list of p-values for a GeneSetCollection.

Usage

pValue(object)

Arguments

object

An object of type GeneSetResult, VEGASResult, GLOSSIResult, GeneSetResultCollection, VEGASResultCollection, or GLOSSIResultCollection

Methods

Defined methods include:

signature(object = "VEGASResult")

Returns the calculated p-value for a specified VEGASResult object

signature(object = "GeneSetResult")

Returns the calculated p-value for a specified GeneSetResult object

signature(object = "GLOSSIResult")

Returns the calculated p-value for a specified GLOSSIResult object

signature(object = "VEGASResultCollection")

Returns a list of calculated p-values for a VEGASResultCollection object (1 for each set)

signature(object = "GeneSetResultCollection")

Returns a list of calculated p-values for a GeneSetResultCollection object (1 for each set)

signature(object = "GLOSSIResultCollection")

Returns a list of calculated p-values for a GLOSSIResultCollection object (1 for each set)

Value

Returns a decimal p-value corresponding to the gene set or a list of p-values corresponding to the gene set collection.

Author(s)

Caitlin McHugh

See Also

VEGASResult-class, geneSetName

Examples

## Not run: 
	pValue( vegas(geneSet, assoc_table, ldMatrix) ) 
	
## End(Not run)

Simulate a Chi-squared Statistic from a Distribution with a Specified Correlation Structure

Description

Simulates a value from a chi-squared distribution with a specified correlation structure.

Usage

simulate_chisq(corr_matrix, num_sims)

Arguments

corr_matrix

A square matrix of correlations.

num_sims

An integer value for the number of simulations to be performed.

Details

This function simulates a value from a chi-squared distribution with a specified correlation structure and degrees of freedom corresponding to the number of rows in the correlation matrix.

Value

A numeric value corresponding to a chi-squared statistic simulated from a chi-squared-like distribution with specified correlation structure and degrees of freedom corresponding to the number of rows in the correlation matrix.

Author(s)

Caitlin McHugh [email protected]


~~ Methods for Function simulatedStats ~~

Description

This function returns the simulated statistics for a specified GeneSet or a list of simulated statistics for a GeneSetCollection.

Usage

simulatedStats(object)

Arguments

object

An object of type VEGASResult or VEGASResultCollection

Methods

Defined methods include:

signature(object = "VEGASResult")

Returns the simulated statistics for a specified VEGASResult object

signature(object = "VEGASResultCollection")

Returns a list of simulated statistics for a VEGASResultCollection object (1 for each set)

Value

Returns decimal values of simulated statistics from a VEGAS function call or a list of simulated statistics from a VEGAS function call on a gene set collection.

Author(s)

Caitlin McHugh

See Also

VEGASResult-class, observedStat

Examples

## Not run: 
	simulatedStats( vegas(geneSet, assoc_table, ldMatrix) ) 
	
## End(Not run)

~~ Methods for Function statistic ~~

Description

This function returns the calculated statistic for a specified GeneSet or a list of calculated statistics for a GeneSetCollection.

Usage

statistic(object)

Arguments

object

An object of type GeneSetResult, VEGASResult, GLOSSIResult, GeneSetResultCollection, VEGASResultCollection or GLOSSIResultCollection

Methods

Defined methods include:

signature(object = "GeneSetResult")

Returns the calculated statistic for a specified GeneSetResult object

signature(object = "GeneSetResultCollection")

Returns a list of calculated statistics for a GeneSetResultCollection object (1 for each set)

signature(object = "VEGASResult")

Returns the calculated statistic for a specified VEGASResult object

signature(object = "VEGASResultCollection")

Returns a list of calculated statistics for a VEGASResultCollection object (1 for each set)

signature(object = "GLOSSIResult")

Returns the calculated statistic for a specified GLOSSIResult object

signature(object = "GLOSSIResultCollection")

Returns a list of calculated statistics for a GLOSSIResultCollection object (1 for each set)

Value

Returns decimal values or a list of decimal values corresponding to calculated statistics from a GLOSSI function call.

Author(s)

Caitlin McHugh

See Also

GLOSSIResult-class, pValue

Examples

## Not run: 
	statistic( glossi(geneSet, assoc_table, ldMatrix) ) 
	
## End(Not run)

~~ Methods for Function show in Package base ~~

Description

~~ Methods for function show in package base ~~

Methods

Defined methods include:

signature(x = "VEGASResult"), signature(x = "GLOSSIResult")

These methods display the corresponding statistics and p-values for the corresponding VEGAS or GLOSSI analysis.

Author(s)

Caitlin McHugh


Calculate the P-Value for a Set Using the VEGAS Method

Description

Calculates the p-value representing the association of the set with the phenotype of interest.

Usage

vegas(set, assoc_table, ldMatrix, num_sims = 1000, correction = TRUE,
  seed = NULL, verbose = FALSE)

Arguments

set

A GeneSet object containing a set of genes of interest or a GeneSetCollection object containing a collection of GeneSets.

assoc_table

An object of class GRanges. This object should at least contain columns SNP and P which hold SNP rsIDs and their corresponding association test p-values, respectively.

ldMatrix

A square, symmetric matrix of LD values, with each row and column corresponding to each of the items in the set. The diagonal entries should be 1, indicating the LD between an item in the set and itself is 1.

num_sims

A positive integer value for the number of simulations to be performed. Default is 1000.

correction

A logical argument indicating whether a value of one should be added to the numerator and denominator when calculating the p-value based upon the simulated statistics. By default, the correction is added and this argument is TRUE.

seed

An integer argument indicating what the random seed should be set to. This allows for replication of results. The default is NULL, and a random seed will be set internally.

verbose

A logical argument indicating whether output should be printed. The default is FALSE.

Details

This function calculates a p-value for sets of SNPs that reside within a gene set collection. We calculate the null distribution by taking into account the observed correlation among the SNPs and simulating a specified number of statistics from which the resulting p-value is calculated.

Value

An object with the corresponding VEGAS results. If set is a GeneSetCollection (i.e., multiple sets of interest), then the corresponding VEGASResultCollection is returned. If set is a GeneSet, a VEGASResult object will be returned.

Author(s)

Caitlin McHugh [email protected]

References

Liu, Jimmy Z. and Mcrae, Allan F. et al. A Versatile Gene-Based Test for Genome-Wide Association Studies. The American Journal of Human Genetics, 2010.


Calculate the P-Value for a Set Using the VEGAS Method

Description

Calculates the p-value representing the association of the set with the phenotype of interest.

Usage

vegasMarginal(pvals, ld_matrix, num_sims, correction = TRUE, seed = NULL,
  verbose = FALSE)

Arguments

pvals

A vector of p-values corresponding to items in the set.

ld_matrix

A square, symmetric matrix of LD values, with each row and column corresponding to each of the items in the set. The diagonal entries should be 1, indicating the LD between an item in the set and itself is 1.

num_sims

An integer value for the number of simulations to be performed.

correction

A logical argument indicating whether a value of one should be added to the numerator when calculating the p-value based upon the simulated statistics. By default, the correction is added. An argument of FALSE will not add one to the numerator.

seed

An integer argument indicating what the random seed should be set to. This allows for replication of results. The default is NULL, and a random seed will be set internally.

verbose

A logical argument indicating whether periodic output should be printed. Defaults to FALSE, indicating no output will be printed.

Details

This is a helper function to calculate the p-value for a set of SNPs that reside within a gene set collection. The correlation among the SNPs is taken into account by the LD matrix. The resulting p-value is calculated from a null distribution that is simulated num_sims times based upon the specified correlation structure.

Value

A VEGASResult object with the corresponding VEGAS results.

Author(s)

Caitlin McHugh [email protected]


Class "VEGASResult"

Description

Objects of this class store results from running VEGAS.

Objects from the Class

Objects can be created by calls of vegas.

Slots

geneSetName:

Object of class "character", the name of the geneSet

pValue:

Object of class "numeric", the p-value

degreesOfFreedom:

Object of class "integer", the degrees of freedom used to calculate the p-value

statistic:

Object of class "numeric", the test statistic for the set

simulatedStats:

Object of class "vector", the simulated statistics for the set used to calculate the p-value

Methods

No methods defined with class "VEGASResult" in the signature.

Accessors

In the code snippets below, object is a VEGASResult object.

geneSetName(object):

The name of the gene set.

pValue(object):

A numeric p-value.

degreesOfFreedom(object):

Integer value of the degrees of freedom used in the statistical test.

statistic(object):

A numeric value corresponding to the statistic calculated from the test performed.

Author(s)

Caitlin McHugh

See Also

VEGASResultCollection, GLOSSIResult

Examples

showClass("VEGASResult")

Class "VEGASResultCollection"

Description

Objects of this class store results from running VEGAS. VEGASResultCollection objects contain a list of VEGASResult objects (one result for each GeneSet)

Objects from the Class

Objects can be created by calls of vegas.

Methods

No methods defined with class "VEGASResult" in the signature.

Accessors

In the code snippets below, object is a VEGASResultCollection object.

geneSetName(object):

A list of the names of the gene sets.

pValue(object):

A list of numeric p-values.

degreesOfFreedom(object):

List of integer values of the degrees of freedom used in the statistical tests performed on each gene set.

statistic(object):

A list of numeric values corresponding to the statistics calculated from the tests performed on each gene set.

Author(s)

Caitlin McHugh

See Also

VEGASResult, GLOSSIResultCollection

Examples

showClass("VEGASResultCollection")