Last data update: 2014.03.03

R: Binomial and poisson goodness of fit statistics for BSSeq...
GoodnessOfFitR Documentation

Binomial and poisson goodness of fit statistics for BSSeq objects

Description

Binomial and poisson goodness of fit statistics for BSSeq objects, including plotting capability.

Usage

poissonGoodnessOfFit(BSseq, nQuantiles = 10^5)
binomialGoodnessOfFit(BSseq, method = c("MLE"), nQuantiles = 10^5)
## S3 method for class 'chisqGoodnessOfFit'
print(x, ...)
## S3 method for class 'chisqGoodnessOfFit'
plot(x, type = c("chisq", "pvalue"), plotCol = TRUE, qqline = TRUE,
  pch = 16, cex = 0.75, ...)

Arguments

BSseq

An object of class BSseq.

x

A chisqGoodnessOfFit object (as produced by poissonGoodnessOfFit or binomialGoodnessOfFit).

nQuantiles

The number of (evenly-spaced) quantiles stored in the return object.

method

How is the parameter estimated.

type

Are the chisq or the p-values being plotted.

plotCol

Should the extreme quantiles be colored.

qqline

Add a qqline.

pch, cex

Plotting symbols and size.

...

Additional arguments being passed to qqplot (for plot) or ignored (for print).

Details

These functions compute and plot goodness of fit statistics for BSseq objects. For each methylation loci, the Poisson goodness of fit statistic tests whether the coverage (at that loci) is independent and identically Poisson distributed across the samples. In a similar fashion, the Binomial goodness of fit statistic tests whether the number of reads supporting methylation are independent and identically binomial distributed across samples (with different size parameters given by the coverage vector).

These functions do not handle NA values.

Value

The plotting method is invoked for its side effect. Both poissonGoodnessOfFit and binomialGoodnessOfFit returns an object of class chisqGoodnessOfFit which is a list with components

chisq

a vector of Chisq values.

quantiles

a vector of quantiles (of the chisq values).

df

degress of freedom

Author(s)

Kasper Daniel Hansen khansen@jhsph.edu

See Also

For the plotting method, see qqplot.

Examples

if(require(bsseqData)) {
  data(BS.cancer.ex)
  BS.cancer.ex <- updateObject(BS.cancer.ex)
  gof <- poissonGoodnessOfFit(BS.cancer.ex)
  plot(gof)
}

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(bsseq)
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

Loading required package: GenomicRanges
Loading required package: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

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

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
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")'.

Loading required package: limma

Attaching package: 'limma'

The following object is masked from 'package:BiocGenerics':

    plotMA

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/bsseq/goodnessOfFit.Rd_%03d_medium.png", width=480, height=480)
> ### Name: GoodnessOfFit
> ### Title: Binomial and poisson goodness of fit statistics for BSSeq
> ###   objects
> ### Aliases: poissonGoodnessOfFit poissonGoodnessOfFit chisqGoodnessOfFit
> ###   binomialGoodnessOfFit print.chisqGoodnessOfFit
> ###   plot.chisqGoodnessOfFit
> 
> ### ** Examples
> 
> if(require(bsseqData)) {
+   data(BS.cancer.ex)
+   BS.cancer.ex <- updateObject(BS.cancer.ex)
+   gof <- poissonGoodnessOfFit(BS.cancer.ex)
+   plot(gof)
+ }
Loading required package: bsseqData
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>