R: Get a Data Frame of Quality Data from a 'FASTQSummary' object
getQual-methods
R Documentation
Get a Data Frame of Quality Data from a FASTQSummary
object
Description
An object of class FASTQSummary contains quality data (binned by readSeqFile). getQual is
an accessor function that reshapes the data into a data frame.
This accessor function is useful if you want to map variables to
custom ggplot2 aesthetics.
Usage
getQual(x)
Arguments
x
an S4 object of class FASTQSummary from
readSeqFile.
Value
getQual returns a data.frame with columns:
position
the position in the read.
ymin
the minimum quality found per a position in the read.
alt.lower
the 10% quantile found per a position in the read.
lower
the 25% quartile found per a position in the read.
middle
the median found per a position in the read.
upper
the 75% quartile found per a position in the read.
alt.upper
the 90% quantile found per a position in the read.
ymax
the maximum quality found per a position in the read.
mean
the mean quality (calculated from the binned data by using a weighted mean function)
per the position in the read.
Methods
signature(x="FASTQSummary")
getQual is an accessor function that only works if there is
quality data, thus it only works with objects of class
FASTQSummary.
## Load a FASTQ file, with sequence hashing.
s.fastq <- readSeqFile(system.file('extdata', 'test.fastq', package='qrqc'))
## Mean quality by position
p <- ggplot(getQual(s.fastq)) + geom_line(aes(x=position, y=mean), color="blue")
p <- p + scale_y_continuous(limits=c(0, 42))
p
## A different type of quality plot
p <- ggplot(getQual(s.fastq)) + geom_linerange(aes(x=position,
ymin=lower, ymax=upper, color=mean))
p <- p + scale_color_gradient("mean quality", low="red", high="green")
p + scale_y_continuous("quality")
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/getQual-methods.Rd_%03d_medium.png", width=480, height=480)
> ### Name: getQual-methods
> ### Title: Get a Data Frame of Quality Data from a 'FASTQSummary' object
> ### Aliases: getQual getQual-methods getQual,FASTQSummary-method
> ### Keywords: methods accessor
>
> ### ** Examples
>
> ## Load a FASTQ file, with sequence hashing.
> s.fastq <- readSeqFile(system.file('extdata', 'test.fastq', package='qrqc'))
>
> ## Mean quality by position
> p <- ggplot(getQual(s.fastq)) + geom_line(aes(x=position, y=mean), color="blue")
> p <- p + scale_y_continuous(limits=c(0, 42))
> p
>
> ## A different type of quality plot
> p <- ggplot(getQual(s.fastq)) + geom_linerange(aes(x=position,
+ ymin=lower, ymax=upper, color=mean))
> p <- p + scale_color_gradient("mean quality", low="red", high="green")
> p + scale_y_continuous("quality")
>
>
>
>
>
> dev.off()
null device
1
>