Last data update: 2014.03.03

R: Get a Data Frame of Simulated Qualitied from a 'FASTQSummary'...
getMCQual-methodsR Documentation

Get a Data Frame of Simulated Qualitied from a FASTQSummary object

Description

An object that inherits from class FASTQSummary contains base quality data by position gathered by readSeqFile. getMCQual generates simulated quality data for each base from this binned quality data that can be used for adding smoothed lines via lowess.

This accessor function is useful if you want to map variables to custom ggplot2 aesthetics.

Usage

  getMCQual(x, n=100)

Arguments

x

an S4 object that inherits from FASTQSummary from readSeqFile.

n

a numeric value indicating the number of quality values to draw per base.

Value

getMCQual returns a data.frame with columns:

position

the position in the read.

quality

simulated quality scores.

Methods

signature(x = "FASTQSummary")

getMCQual is a function that works on any object with class FASTQSummary read in with readSeqFile.

Author(s)

Vince Buffalo <vsbuffalo@ucdavis.edu>

See Also

getGC, getSeqlen, getBase, getBaseProp, getQual, qualPlot

Examples

  ## Load a FASTQ file, with sequence hashing.
  s.fastq <- readSeqFile(system.file('extdata', 'test.fastq',
    package='qrqc'))

  # A custom quality plot
  ggplot(getQual(s.fastq)) + geom_linerange(aes(x=position, ymin=lower,
    ymax=upper), color="grey") + geom_smooth(aes(x=position, y=quality),
    data=getMCQual(s.fastq), color="blue", se=FALSE)

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(qrqc)
Loading required package: reshape
Loading required package: ggplot2
Loading required package: Biostrings
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: S4Vectors
Loading required package: stats4

Attaching package: 'S4Vectors'

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

    expand, rename

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

    colMeans, colSums, expand.grid, rowMeans, rowSums

Loading required package: IRanges
Loading required package: XVector
Loading required package: biovizBase
Loading required package: brew
Loading required package: xtable
Loading required package: Rsamtools
Loading required package: GenomeInfoDb
Loading required package: GenomicRanges
Loading required package: testthat

Attaching package: 'testthat'

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

    compare

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/qrqc/getMCQual-methods.Rd_%03d_medium.png", width=480, height=480)
> ### Name: getMCQual-methods
> ### Title: Get a Data Frame of Simulated Qualitied from a 'FASTQSummary'
> ###   object
> ### Aliases: getMCQual getMCQual-methods getMCQual,FASTQSummary-method
> ### Keywords: methods accessor
> 
> ### ** Examples
> 
>   ## Load a FASTQ file, with sequence hashing.
>   s.fastq <- readSeqFile(system.file('extdata', 'test.fastq',
+     package='qrqc'))
> 
>   # A custom quality plot
>   ggplot(getQual(s.fastq)) + geom_linerange(aes(x=position, ymin=lower,
+     ymax=upper), color="grey") + geom_smooth(aes(x=position, y=quality),
+     data=getMCQual(s.fastq), color="blue", se=FALSE)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>