Last data update: 2014.03.03

R: Low-level functions to fit dispersion estimates
estimateDispersionsGeneEstR Documentation

Low-level functions to fit dispersion estimates

Description

Normal users should instead use estimateDispersions. These low-level functions are called by estimateDispersions, but are exported and documented for non-standard usage. For instance, it is possible to replace fitted values with a custom fit and continue with the maximum a posteriori dispersion estimation, as demonstrated in the examples below.

Usage

estimateDispersionsGeneEst(object, minDisp = 1e-08, kappa_0 = 1,
  dispTol = 1e-06, maxit = 100, quiet = FALSE, modelMatrix = NULL,
  niter = 1)

estimateDispersionsFit(object, fitType = c("parametric", "local", "mean"),
  minDisp = 1e-08, quiet = FALSE)

estimateDispersionsMAP(object, outlierSD = 2, dispPriorVar, minDisp = 1e-08,
  kappa_0 = 1, dispTol = 1e-06, maxit = 100, modelMatrix = NULL,
  quiet = FALSE)

estimateDispersionsPriorVar(object, minDisp = 1e-08, modelMatrix = NULL)

Arguments

object

a DESeqDataSet

minDisp

small value for the minimum dispersion, to allow for calculations in log scale, one order of magnitude above this value is used as a test for inclusion in mean-dispersion fitting

kappa_0

control parameter used in setting the initial proposal in backtracking search, higher kappa_0 results in larger steps

dispTol

control parameter to test for convergence of log dispersion, stop when increase in log posterior is less than dispTol

maxit

control parameter: maximum number of iterations to allow for convergence

quiet

whether to print messages at each step

modelMatrix

for advanced use only, a substitute model matrix for gene-wise and MAP dispersion estimation

niter

number of times to iterate between estimation of means and estimation of dispersion

fitType

either "parametric", "local", or "mean" for the type of fitting of dispersions to the mean intensity. See estimateDispersions for description.

outlierSD

the number of standard deviations of log gene-wise estimates above the prior mean (fitted value), above which dispersion estimates will be labelled outliers. Outliers will keep their original value and not be shrunk using the prior.

dispPriorVar

the variance of the normal prior on the log dispersions. If not supplied, this is calculated as the difference between the mean squared residuals of gene-wise estimates to the fitted dispersion and the expected sampling variance of the log dispersion

Value

a DESeqDataSet with gene-wise, fitted, or final MAP dispersion estimates in the metadata columns of the object.

estimateDispersionsPriorVar is called inside of estimateDispersionsMAP and stores the dispersion prior variance as an attribute of dispersionFunction(dds), which can be manually provided to estimateDispersionsMAP for parallel execution.

See Also

estimateDispersions

Examples


dds <- makeExampleDESeqDataSet()
dds <- estimateSizeFactors(dds)
dds <- estimateDispersionsGeneEst(dds)
dds <- estimateDispersionsFit(dds)
dds <- estimateDispersionsMAP(dds)
plotDispEsts(dds) 

# after having run estimateDispersionsFit()
# the dispersion prior variance over all genes
# can be obtained like so:

dispPriorVar <- estimateDispersionsPriorVar(dds)

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(DESeq2)
Loading required package: S4Vectors
Loading required package: stats4
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit


Attaching package: 'S4Vectors'

The following objects are masked from 'package:base':

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: GenomicRanges
Loading required package: GenomeInfoDb
Loading required package: SummarizedExperiment
Loading required package: Biobase
Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/DESeq2/estimateDispersionsGeneEst.Rd_%03d_medium.png", width=480, height=480)
> ### Name: estimateDispersionsGeneEst
> ### Title: Low-level functions to fit dispersion estimates
> ### Aliases: estimateDispersionsFit estimateDispersionsGeneEst
> ###   estimateDispersionsMAP estimateDispersionsPriorVar
> 
> ### ** Examples
> 
> 
> dds <- makeExampleDESeqDataSet()
> dds <- estimateSizeFactors(dds)
> dds <- estimateDispersionsGeneEst(dds)
> dds <- estimateDispersionsFit(dds)
> dds <- estimateDispersionsMAP(dds)
> plotDispEsts(dds) 
> 
> # after having run estimateDispersionsFit()
> # the dispersion prior variance over all genes
> # can be obtained like so:
> 
> dispPriorVar <- estimateDispersionsPriorVar(dds)
> 
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>