nbinomLRT: Likelihood ratio test (chi-squared test) for GLMs

Description Usage Arguments Details Value See Also Examples

View source: R/core.R

Description

This function tests for significance of change in deviance between a full and reduced model which are provided as formula. Fitting uses previously calculated sizeFactors (or normalizationFactors) and dispersion estimates.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
nbinomLRT(
  object,
  full = design(object),
  reduced,
  betaTol = 1e-08,
  maxit = 100,
  useOptim = TRUE,
  quiet = FALSE,
  useQR = TRUE,
  minmu = if (type == "glmGamPoi") 1e-06 else 0.5,
  type = c("DESeq2", "glmGamPoi")
)

Arguments

object

a DESeqDataSet

full

the full model formula, this should be the formula in design(object). alternatively, can be a matrix

reduced

a reduced formula to compare against, e.g. the full model with a term or terms of interest removed. alternatively, can be a matrix

betaTol

control parameter defining convergence

maxit

the maximum number of iterations to allow for convergence of the coefficient vector

useOptim

whether to use the native optim function on rows which do not converge within maxit

quiet

whether to print messages at each step

useQR

whether to use the QR decomposition on the design matrix X while fitting the GLM

minmu

lower bound on the estimated count while fitting the GLM

type

either "DESeq2" or "glmGamPoi". If type = "DESeq2" a classical likelihood ratio test based on the Chi-squared distribution is conducted. If type = "glmGamPoi" and previously the dispersion has been estimated with glmGamPoi as well, a quasi-likelihood ratio test based on the F-distribution is conducted. It is supposed to be more accurate, because it takes the uncertainty of dispersion estimate into account in the same way that a t-test improves upon a Z-test.

Details

The difference in deviance is compared to a chi-squared distribution with df = (reduced residual degrees of freedom - full residual degrees of freedom). This function is comparable to the nbinomGLMTest of the previous version of DESeq and an alternative to the default nbinomWaldTest.

Value

a DESeqDataSet with new results columns accessible with the results function. The coefficients and standard errors are reported on a log2 scale.

See Also

DESeq, nbinomWaldTest

Examples

1
2
3
4
5
dds <- makeExampleDESeqDataSet()
dds <- estimateSizeFactors(dds)
dds <- estimateDispersions(dds)
dds <- nbinomLRT(dds, reduced = ~ 1)
res <- results(dds)

DESeq2 documentation built on Feb. 22, 2021, 10 a.m.