Package 'stageR'

Title: stageR: stage-wise analysis of high throughput gene expression data in R
Description: The stageR package allows automated stage-wise analysis of high-throughput gene expression data. The method is published in Genome Biology at https://genomebiology.biomedcentral.com/articles/10.1186/s13059-017-1277-0
Authors: Koen Van den Berge and Lieven Clement
Maintainer: Koen Van den Berge <[email protected]>
License: GNU General Public License version 3
Version: 1.27.0
Built: 2024-06-30 04:44:29 UTC
Source: https://github.com/bioc/stageR

Help Index


Get adjusted significance level from the screening stage.

Description

This functions returns the adjusted significance level from the screening stage that should be used to compare confirmation stage FWER adjusted p-values against.

Usage

adjustedAlphaLevel(object, ...)

## S4 method for signature 'stageR'
adjustedAlphaLevel(object)

## S4 method for signature 'stageRTx'
adjustedAlphaLevel(object)

Arguments

object

an object of the stageRClass class.

Details

The adjusted significance level is calculated as the fraction of significant features in the screening stage multiplied the alpha level.

Value

Scalar, the adjusted significance level from the screening stage.

Methods (by class)

  • stageRTx: Get adjusted significance level from the screening stage.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

R. Heller, E. Manduchi, G. R. Grant, and W. J. Ewens, "A flexible two-stage procedure for identifying gene sets that are differentially expressed." Bioinformatics (Oxford, England), vol. 25, pp. 1019-25, 2009.

See Also

stageR, stageRClass

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation, pScreenAdjusted=FALSE)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
adjustedAlphaLevel(stageRObj)
# @method stageR-method

Transcript-level abundance estimates in 14 Chinese prostate cancer patients

Description

A dataset containing 14 matched samples of tumoral prostate cancer and normal tissue, both derived from the same Chinese patient. The dataset has been prefiltered to reduce the computational burden of the vignette.

Usage

esetProstate

Format

An ExpressionSet object.

Source

http://pachterlab.github.io/lair/

References

Ren, Shancheng, Zhiyu Peng, Jian-Hua Mao, Yongwei Yu, Changjun Yin, Xin Gao, Zilian Cui, et al. 2012. "RNA-seq analysis of prostate cancer in the Chinese population identifies recurrent gene fusions, cancer-associated long noncoding RNAs and aberrant alternative splicings." Cell Research 22 (5): 806-21. doi:10.1038/cr.2012.30.


Retrieve the stage-wise adjusted p-values.

Description

This functions returns the stage-wise adjusted p-values for an object from the stageRClass class. Note, that the p-values should have been adjusted with the stageWiseAdjustment,stageR,character,numeric-method function prior to calling this function.

Usage

getAdjustedPValues(object, onlySignificantGenes, order, ...)

## S4 method for signature 'stageR,logical,logical'
getAdjustedPValues(object,
  onlySignificantGenes, order, ...)

## S4 method for signature 'stageRTx,logical,logical'
getAdjustedPValues(object,
  onlySignificantGenes, order, ...)

Arguments

object

an object of the stageRClass class.

onlySignificantGenes

logical. If FALSE (default), all genes are returned. If TRUE, only the genes significant for the screening hypothesis are returned.

order

logical. If TRUE (default), the returned matrix of adjusted p-values are ordered based on the screening hypothesis adjusted p-value.

...

Other arguments passed to .getAdjustedP or .getAdjustedPTx

Details

The function returns FDR adjusted p-values for the screening hypothesis and stage-wise adjusted p-values for the confirmation hypothesis p-values. For features that were not significant in the screening hypothesis, the confirmation stage adjusted p-values are set to NA. The adjusted p-values in the output of getAdjustedPValues can directly be compared to alpha, the OFDR level specified in stageWiseAdjustment, to flag significant features.

Value

For complex DGE experiments (stageR object), a matrix of adjusted p-values where every row corresponds to a gene, and every column corresponds to a contrast. The first column will be the BH FDR adjusted p-value from the screening step. For transcript-level experiments (stageRTx object), a matrix of adjusted p-values where every row corresponds to a transcript.

Methods (by class)

  • object = stageRTx,onlySignificantGenes = logical,order = logical: Retrieve the stage-wise adjusted p-values.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
head(getAdjustedPValues(stageRObj, onlySignificantGenes=TRUE, order=TRUE))

Retrieve the significance level for the stage-wise adjustment.

Description

This functions returns the significance level on which the stage-wise adjustment is based.

Usage

getAlpha(object, ...)

## S4 method for signature 'stageR'
getAlpha(object, ...)

## S4 method for signature 'stageRTx'
getAlpha(object, ...)

Arguments

object

an object of the stageRClass or stageRTxClass class.

...

Additional arguments

Value

Returns a calar vector with the OFDR alpha level that was specified by the user.

Methods (by class)

  • stageRTx: Retrieve the significance level for the stage-wise adjustment.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
getAlpha(stageRObj)

Retrieve FWER correction method.

Description

This functions retrieves the method used for FWER multiple testing correction in the confirmation stage of a stage-wise analysis.

Usage

getMethod(object, ...)

## S4 method for signature 'stageR'
getMethod(object, ...)

## S4 method for signature 'stageRTx'
getMethod(object, ...)

Arguments

object

an object of the stageRClass or stageRTxClass class.

...

Additional arguments

Value

Returns a character vector of length 1 specifying the FWER correction method that is used in the confirmation stage of the stage-wise analysis.

Methods (by class)

  • stageRTx: Retrieve FWER correction method.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
getMethod(stageRObj)

Return unadjusted confirmation hypothesis p-values from a stageRClass object.

Description

Return unadjusted confirmation hypothesis p-values from a stageRClass object.

Usage

getPConfirmation(object, ...)

## S4 method for signature 'stageR'
getPConfirmation(object)

## S4 method for signature 'stageRTx'
getPConfirmation(object)

Arguments

object

an object of the stageRClass class.

Value

A matrix of the unadjusted p-values to be used in the confirmation stage.

Methods (by class)

  • stageRTx: Return unadjusted confirmation hypothesis p-values from a stageRClass object.

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
getPConfirmation(stageRObj)

Return screening hypothesis p-values from a stageRClass object.

Description

Return screening hypothesis p-values from a stageRClass object.

Usage

getPScreen(object, ...)

## S4 method for signature 'stageR'
getPScreen(object)

## S4 method for signature 'stageRTx'
getPScreen(object)

Arguments

object

an object of the stageRClass class.

...

Additional arguments

Value

A vector of screening stage (aggregated) p-values.

Methods (by class)

  • stageRTx: Return screening hypothesis p-values from a stageRClass object.

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
getPScreen(stageRObj)

Get significance results according to a stage-wise analysis.

Description

This functions returns a matrix that indicates whether a specific feature is significant for a specific hypothesis of interest according to a stage-wise analysis. The function is not applicable to transcript-level analysis.

Usage

getResults(object, ...)

## S4 method for signature 'stageR'
getResults(object)

Arguments

object

an object of the stageRClass class.

Details

The FDR adjusted screening hypothesis p-values are compared to the alpha level specified. The FWER adjusted confirmation stage p-values are compared to the adjusted significance level from the screening stage.

Value

A logical matrix with rows corresponding to genes and columns corresponding to contrasts, where the first column represents the screening stage on the aggregated p-values. A 0 represents a non-significant test, a 1 represents a significant test according to the stage-wise analysis.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
head(getResults(stageRObj))

Return significant genes when performing transcript-level analysis.

Description

This functions returns a matrix with significant genes by aggregated testing of its respective transcripts.

Usage

getSignificantGenes(object, ...)

## S4 method for signature 'stageRTx'
getSignificantGenes(object)

Arguments

object

an object of the stageRClass class.

Value

A matrix with significant genes and their corresponding FDR-adjusted screening stage (aggregated) p-value.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

#make identifiers linking transcripts to genes
set.seed(1)
genes=paste0("gene",sample(1:200,1000,replace=TRUE))
nGenes=length(table(genes))
transcripts=paste0("tx",1:1000)
tx2gene=data.frame(transcripts,genes)
#gene-wise q-values
pScreen=c(seq(1e-10,1e-2,length.out=nGenes-100),seq(1e-2,.2,length.out=50),seq(50))
names(pScreen)=names(table(genes)) #discards genes that are not simulated
pConfirmation=matrix(runif(1000),nrow=1000,ncol=1)
rownames(pConfirmation)=transcripts
stageRObj <- stageRTx(pScreen=pScreen, pConfirmation=pConfirmation ,pScreenAdjusted=TRUE, tx2gene=tx2gene)
stageRObj <- stageWiseAdjustment(stageRObj, method="dte", alpha=0.05)
head(getSignificantGenes(stageRObj))

Return significant transcripts when performing transcript-level analysis.

Description

This functions returns a matrix with significant transctripts according to a stage-wise analysis.

Usage

getSignificantTx(object, ...)

## S4 method for signature 'stageRTx'
getSignificantTx(object)

Arguments

object

an object of the stageRClass class.

Value

A matrix of significant transcripts with their corresponding stage-wise adjusted p-value (i.e. FDR and FWER correction).

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

#make identifiers linking transcripts to genes
set.seed(1)
genes=paste0("gene",sample(1:200,1000,replace=TRUE))
nGenes=length(table(genes))
transcripts=paste0("tx",1:1000)
tx2gene=data.frame(transcripts,genes)
#gene-wise q-values
pScreen=c(seq(1e-10,1e-2,length.out=nGenes-100),seq(1e-2,.2,length.out=50),seq(50))
names(pScreen)=names(table(genes)) #discards genes that are not simulated
pConfirmation=matrix(runif(1000),nrow=1000,ncol=1)
rownames(pConfirmation)=transcripts
stageRObj <- stageRTx(pScreen=pScreen, pConfirmation=pConfirmation ,pScreenAdjusted=TRUE, tx2gene=tx2gene)
stageRObj <- stageWiseAdjustment(stageRObj, method="dte", alpha=0.05)
head(getSignificantTx(stageRObj))

Retrieve the data frame linking genes to transcripts.

Description

This functions returns a data frame that links the genes with the transcripts being analysed.

Usage

getTx2gene(object, ...)

## S4 method for signature 'stageRTx'
getTx2gene(object, ...)

Arguments

object

an object of the stageRTxClass class.

...

Additional arguments

Value

A matrix linking gene to transcript identifiers.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

#make identifiers linking transcripts to genes
set.seed(1)
genes=paste0("gene",sample(1:200,1000,replace=TRUE))
nGenes=length(table(genes))
transcripts=paste0("tx",1:1000)
tx2gene=data.frame(transcripts,genes)
#gene-wise q-values
pScreen=c(seq(1e-10,1e-2,length.out=nGenes-100),seq(1e-2,.2,length.out=50),seq(50))
names(pScreen)=names(table(genes)) #discards genes that are not simulated
pConfirmation=matrix(runif(1000),nrow=1000,ncol=1)
rownames(pConfirmation)=transcripts
stageRObj <- stageRTx(pScreen=pScreen, pConfirmation=pConfirmation ,pScreenAdjusted=TRUE, tx2gene=tx2gene)
getTx2gene(stageRObj)

Hammer dataset

Description

A gene expression dataset from an experiment on spinal nerve ligation in rats, comparing this treatment to control samples in two timepoints, i.e. two weeks and two months post treatment. 2 Biological replicates available in every treatment x time combination.

Usage

hammer.eset

Format

An ExpressionSet object.

Source

http://bowtie-bio.sourceforge.net/recount/

References

Hammer P, Banck MS, Amberg R, et al. mRNA-seq with agnostic splice site discovery for nervous system transcriptomics tested in chronic pain. Genome Research. 2010;20(6):847-860. doi:10.1101/gr.101204.109.


Has stage-wise adjustment already been performed on the object?

Description

This functions returns a logical stating whether the p-values have already been adjusted according to the stage-wise method.

Usage

isAdjusted(object, ...)

## S4 method for signature 'stageR'
isAdjusted(object, ...)

## S4 method for signature 'stageRTx'
isAdjusted(object, ...)

Arguments

object

an object of the stageRClass or stageRTxClass class.

...

Additional arguments

Value

A logical stating whether the p-values have already been adjusted according to the stage-wise method

Methods (by class)

  • stageRTx: Has stage-wise adjustment already been performed on the object?

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
isAdjusted(stageRObj)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
isAdjusted(stageRObj)

Are the screening p-values adjusted for multiplicity?

Description

This functions returns a logical stating whether the screening hypothesis p-values are already adjusted for multiple testing according to the BH FDR criterion.

Usage

isPScreenAdjusted(object, ...)

## S4 method for signature 'stageR'
isPScreenAdjusted(object, ...)

## S4 method for signature 'stageRTx'
isPScreenAdjusted(object, ...)

Arguments

object

an object of the stageRClass or stageRTxClass class.

...

Additional arguments

Value

A logical stating whether the screening hypothesis p-values are already adjusted for multiple testing according to the BH FDR criterion.

Methods (by class)

  • stageRTx: Are the screening p-values adjusted for multiplicity?

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
isPScreenAdjusted(stageRObj)

Create stageR object

Description

Constructor function for stageRClass. A stageR class is a class used for stage-wise analysis in high throughput settings. In its most basic form, it consists of a vector of p-values for the screening hypothesis and a matrix of p-values for the confirmation hypotheses.

Usage

stageR(pScreen, pConfirmation, pScreenAdjusted = FALSE)

Arguments

pScreen

A vector of screening hypothesis p-values.

pConfirmation

A matrix of confirmation hypothesis p-values. When constructing a stageRClass object, the number of rows should be equal to the length of pScreen. For a stageRTxClass object, the dimensions can be different.

pScreenAdjusted

logical, indicating whether the supplied p-values for the screening hypothesis have already been adjusted for multiplicity according to the FDR.

...

Additional arguments.

Value

An instance of an object of the stageRClass

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

# create a \code{\link{stageRClass}} object
pScreen <- runif(10)
names(pScreen) <- paste0("gene",1:10)
pConfirmation <- matrix(runif(30),nrow=10,ncol=3)
rownames(pConfirmation) <-  paste0("gene",1:10)
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
pConfirmationTx <- matrix(runif(10),ncol=1)
names(pScreen) <- paste0("gene",rep(1:2,each=5))
stageRObj <- stageRTx(pScreen=pScreen, pConfirmation=pConfirmationTx, tx2gene=data.frame(transcripts=paste0("transcript",1:10),genes=paste0("gene",rep(1:2,each=5))))

The stageR class

Description

This class is used for adjusting p-values with stage-wise testing for high-throughput studies.

Slots

pScreen

A vector of p-values for the screening hypothesis.

pConfirmation

A matrix of p-values for the confirmation hypotheses.

adjustedP

A matrix of adjusted p-values. This slot should be accessed through getAdjustedPValues,stageR,logical,logical-method. Alternatively, significance results can be accessed through getResults,stageR-method.

method

Character string indicating the method used for FWER correction in the confirmation stage of the stage-wise analysis. Can be any of "none", "holm", "dte", "dtu", "user". "none" will not adjust the p-values in the confirmation stage. "holm" is an adapted Holm procedure for a stage-wise analysis, where the method takes into account the fact that genes in the confirmation stage have already passed the screening stage, hence the procedure will be more powerful for the most significant p-value as compared to the standard Holm procedure. "dte" is the adjusted Holm-Shaffer procedure for differential transcript expression analysis. "dtu" is the adjusted Holm-Shaffer procedure for differential transcript usage. "user" indicates a user-defined adjustment that should be specified with the adjustment argument.

alpha

the OFDR level on which the stage-wise analysis should be controlled.

alphaAdjusted

the adjusted significance level to compare against FWER-adjusted p-values of the confirmation stage to decide on significance of the hypothesis test.

pScreenAdjusted

logical, indicating whether the supplied p-values for the screening hypothesis have already been adjusted for multiplicity according to the FDR.

tx2gene

matrix with transcript IDs in the first column and gene IDs in the second column to be used for DTE and DTU analysis. All rownames from pConfirmation should match with a transcript ID and all names from pScreen should match with a gene ID.

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0 R. Heller, E. Manduchi, G. R. Grant, and W. J. Ewens, "A flexible two-stage procedure for identifying gene sets that are differentially expressed." Bioinformatics (Oxford, England), vol. 25, pp. 1019-25, 2009. S. Holm, "A Simple Sequentially Rejective Multiple Test Procedure," Scandinavian Journal of Statistics, vol. 6, no. 2, pp. 65-70, 1979. J. P. Shaffer, "Modified Sequentially Rejective Multiple Test Procedures," Journal of the American Statistical Association, vol. 81, p. 826, 1986.


Create stageRTx object.

Description

Constructor function for stageRTxClass. A stageR class is a class used for stage-wise analysis in high throughput settings. In its most basic form, it consists of a vector of p-values for the screening hypothesis, a matrix of p-values for the confirmation hypotheses and a tx2gene object for linking genes to transcripts.

Usage

stageRTx(pScreen, pConfirmation, pScreenAdjusted = FALSE, tx2gene)

Arguments

pScreen

A vector of screening hypothesis p-values.

pConfirmation

A matrix of confirmation hypothesis p-values. The number of rows should be equal to the length of pScreen.

pScreenAdjusted

logical, indicating whether the supplied p-values for the screening hypothesis have already been adjusted for multiplicity according to the FDR.

tx2gene

Only applicable for transcript-level analysis. A data.frame with transcript IDs in the first columns and gene IDs in the second column. The rownames from pConfirmation must be contained in the transcript IDs from tx2gene, and the names from pScreen must be contained in the gene IDs.

...

Additional arguments.

Value

An instance of an object of the stageRTxClass

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0

Examples

# create a \code{\link{stageRClass}} object
pScreen <- runif(10)
names(pScreen) <- paste0("gene",1:10)
pConfirmation <- matrix(runif(30),nrow=10,ncol=3)
rownames(pConfirmation) <-  paste0("gene",1:10)
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
pConfirmationTx <- matrix(runif(10),ncol=1)
names(pScreen) <- paste0("gene",rep(1:2,each=5))
stageRObj <- stageRTx(pScreen=pScreen, pConfirmation=pConfirmationTx, tx2gene=data.frame(transcripts=paste0("transcript",1:10),genes=paste0("gene",rep(1:2,each=5))))

adjust p-values in a two-stage analysis

Description

This function will adjust p-values according to a hierarchical two-stage testing paradigm.

Usage

stageWiseAdjustment(object, method, alpha, ...)

## S4 method for signature 'stageR,character,numeric'
stageWiseAdjustment(object, method, alpha,
  adjustment = NULL, ...)

## S4 method for signature 'stageRTx,character,numeric'
stageWiseAdjustment(object, method,
  alpha, tx2gene, ...)

Arguments

object

an object of the stageRClass class.

method

Character string indicating the method used for FWER correction in the confirmation stage of the stage-wise analysis. Can be any of "none", "holm", "dte", "dtu", "user". "none" will not adjust the p-values in the confirmation stage. "holm" is an adapted Holm procedure for a stage-wise analysis, where the method takes into account the fact that genes in the confirmation stage have already passed the screening stage, hence the procedure will be more powerful for the most significant p-value as compared to the standard Holm procedure. "dte" is the adjusted Holm-Shaffer procedure for differential transcript expression analysis. "dtu" is the adjusted Holm-Shaffer procedure for differential transcript usage. "user" indicates a user-defined adjustment that should be specified with the adjustment argument.

alpha

the OFDR on which to control the two-stage analysis.

...

Additional arguments passed to .stageWiseTest

adjustment

a user-defined adjustment of the confirmation stage p-values. Only applicable when method is "user" and ignored otherwise.

tx2gene

Only applicable when method is "dte" or "dtu". A data.frame with transcript IDs in the first columns and gene IDs in the second column. The rownames from pConfirmation must be contained in the transcript IDs from tx2gene, and the names from pScreen must be contained in the gene IDs.

Value

A stageR/stageRTx object with stage-wise adjusted p-values.

Methods (by class)

  • object = stageRTx,method = character,alpha = numeric: Adjust p-values in a two-stage analysis

References

Van den Berge K., Soneson C., Robinson M.D., Clement L. (2017). stageR: a general stage-wise method for controlling the gene-level false discovery rate in differential expression and differential transcript usage. Genome Biology 18:151. https://doi.org/10.1186/s13059-017-1277-0 R. Heller, E. Manduchi, G. R. Grant, and W. J. Ewens, "A flexible two-stage procedure for identifying gene sets that are differentially expressed." Bioinformatics (Oxford, England), vol. 25, pp. 1019-25, 2009.

S. Holm, "A Simple Sequentially Rejective Multiple Test Procedure," Scandinavian Journal of Statistics, vol. 6, no. 2, pp. 65-70, 1979. J. P. Shaffer, "Modified Sequentially Rejective Multiple Test Procedures," Journal of the American Statistical Association, vol. 81, p. 826, 1986.

Examples

pScreen=c(seq(1e-10,1e-2,length.out=100),seq(1e-2,.2,length.out=100),seq(.2,1,length.out=100))
names(pScreen)=paste0("gene",1:300)
pConfirmation=matrix(runif(900),nrow=300,ncol=3)
dimnames(pConfirmation)=list(paste0("gene",1:300),c("H1","H2","H3"))
stageRObj <- stageR(pScreen=pScreen, pConfirmation=pConfirmation)
stageRObj <- stageWiseAdjustment(stageRObj, method="holm", alpha=0.05)
getAdjustedPValues(stageRObj, onlySignificantGenes=TRUE, order=TRUE)