Rqc 1.20.0
Rqc is an optimized tool designed for quality control and assessment of high-throughput sequencing data. It performs parallel processing of entire files and produces a report which contains a set of high-resolution graphics that can be used for quality assessment.
This version of Rqc produces high-quality images for the following statistics:
The main goal of Rqc is to provide graphical tools for quality
assessment of reads contained in FASTQ files. This package is
designed focusing on simplicity of use. Therefore, the Rqc package
allows the user to call one single function called rqc
. The rqc
method processes a set of input files and generates an HTML report
containing several plots that can be used for quality assessment.
To access this functionality, the user needs to load Rqc package.
library(Rqc)
The next step is to determine the location of the FASTQ files that should be analyzed. The example below, uses sample files provided by the ShortRead package, but the user must modify this location accordingly, in order to reflect the actual location of the files that need QA.
folder <- system.file(package="ShortRead", "extdata/E-MTAB-1147")
The basic usage of the rqc
function requires the definition of 2
arguments. One, path
, is the location where the files of interest
are saved at (this was defined on the step above). The other argument,
pattern
, is a regular expression that identifies all files of
interest. Below, we use .fastq.gz
to specify that all files
containing that string are to be processed.
rqc(path = folder, pattern = ".fastq.gz")
At this point, the user’s default Internet browser will open an HTML file. This file is the report generated by Rqc, which, by default, is stored in a temporary directory. A sample report is shown below:
This table describes input files. reads
column can be total number of reads (sample=FALSE
) or sample size.
knitr::kable(perFileInformation(qa))
filename | pair | format | group | reads | total.reads | path |
---|---|---|---|---|---|---|
ERR127302_1_subset.fastq.gz | 1 | FASTQ | None | 20000 | 20000 | /home/biocbuild/bbs-3.10-bioc/R/library/ShortRead/extdata/E-MTAB-1147 |
ERR127302_2_subset.fastq.gz | 2 | FASTQ | None | 20000 | 20000 | /home/biocbuild/bbs-3.10-bioc/R/library/ShortRead/extdata/E-MTAB-1147 |
This plot describe an overview of per read mean quality distribution of all files
rqcReadQualityBoxPlot(qa)
This plot describes the average quality pattern by showing on the X-axis quality thresholds and on the Y-axis the percentage of reads that exceed that quality level.
rqcReadQualityPlot(qa)
This describes the average quality scores for each cycle of sequencing.
rqcCycleAverageQualityPlot(qa)
This plot shows the proportion of reads that appeared many times.
rqcReadFrequencyPlot(qa)
This heatmap plot shows dstance matrix between top represented reads. This functon only works with one result file (and not a list).
rqcFileHeatmap(qa[[1]])
Barplot that presents the distribuition of the lengths of the reads available in the FASTQ file.
rqcReadWidthPlot(qa)
Line plot showing the average GC content for every cycle of sequencing.
rqcCycleGCPlot(qa)
Bar plot showing the proportion of quality calls per cycle. Colors are presented in a gradient Red-Blue, where red identifies calls of lower quality. This visualization is preferred as it is cleaner than the boxplots described below.
rqcCycleQualityPlot(qa)
Biplot from Principal Component Analysis (PCA) of cycle-specific read average quality.
rqcCycleAverageQualityPcaPlot(qa)
## Warning in geom_hline(aes(0), size = 0.2, yintercept = 0): Using both `yintercept` and `mapping` may not have the desired result as mapping is overwritten if `yintercept` is specified
## Warning in geom_vline(aes(0), size = 0.2, xintercept = 0): Using both `xintercept` and `mapping` may not have the desired result as mapping is overwritten if `xintercept` is specified
Boxplots describing empirical patterns of quality distribution on each cycle of sequencing.
rqcCycleQualityBoxPlot(qa)
This bar plot describes the proportion of each nucleotide called for every cycle of sequencing.
rqcCycleBaseCallsPlot(qa)
The line plot shows a more detailed view.
rqcCycleBaseCallsLinePlot(qa)
It is important to note that the rqc
function samples 1 million
records from the FASTQ files. This can be set by adjusting the n
argument for this function. If the user desires to have the file
processed as a whole (rather than sampling records from it), (s)he
must set the argument sample
to FALSE
.
The rqc
function wraps a set of functions to generate a quick report
that can be used for quality assessment. However, users can perform a
step-by-step analysis by using the information described below.
If one wants to process a set of files that are not located at the
same directory, the users needs to create a vector containing the
absolute path of files. The list.files
function can be useful.
fastqDir <- system.file(package="ShortRead", "extdata/E-MTAB-1147")
files <- list.files(fastqDir, "fastq.gz", full.names=TRUE)
The example input files are samples from a public data set. These samples are available through the ShortRead package. More information regarding these data can be found on the vignette of that package.
To process the files without generating an HTML report, the user
should use rqcQA
function instead rqc
. This function receives a
vector containing the paths of the input files.
qa <- rqcQA(files, workers=1)
The rqcQA
function returns a list that contains the required
information to create the plots present on the standard HTML
report. Actually, both rqc
and rqcQA
return a named list of
RqcResultSet
objects. This output can be used directly as input to
other Rqc package functions. Examples of functions that can use these
objects are rqcReport
and plot
. RqcResultSet
is a class that
extends .QA
class of ShortRead package. An RqcResultSet
object
contains information in two different perspectives: read-specific
data and cycle-specific data. They can be accessed using [[
brackets.
To create a final HTML report, the user can apply the rqcReport
function to the qa
object.
reportFile <- rqcReport(qa)
browseURL(reportFile)
The report generated by rqcReport
contains the all plots described
on the beginning of this document. By default, it is written to a
temporary directory. This behavior can by modified by setting the
outdir
argument.
The Rqc package performs parallel processing of the samples by
interfacing with the BiocParallel package. By default,
Rqc calls multicoreWorkers
function that returns the maximum number
of workers available. You can change the number of cores setting workers
parameter on rqc
and rqcQA
functions. It is possible to process the
input files serially by setting workers=1
.
For each plot generated by Rqc, there is a function that shapes the data appropriately. The shaped information is then used to produce the final plot. The example below shows how the user can access these data to generate plots using other tools.
df <- rqcCycleAverageQualityCalc(qa)
cycle <- as.numeric(levels(df$cycle))[df$cycle]
plot(cycle, df$quality, col = df$filename, xlab='Cycle', ylab='Quality Score')
One can also process a subset of result data using subsets of a list.
sublist <- qa[1]
rqcCycleQualityPlot(sublist)
Rqc package accepts R Markdown file format as template file for generating custom reports. Markdown is a markup language for web development. R Markdown files are regular Markdown files with R codes. Every chunk of code is executed during compilation done by knitr package. Knitr takes R Markdown file and generates Markdown file merged with R codes and their outputs such as text, tables and figures. Result Markdown file is used by Rqc for generating final report in HTML or PDF formats. The source file of Rqc’s default report is a good reference for writing new template reports. Run code below returns the system file path to this source file.
system.file(package = "Rqc", "templates", "rqc_report.Rmd")
Basic concepts for writing Markdown documents with R code are available at http://rmarkdown.rstudio.com/. Advanced concepts about the compilation process of R Markdown files are available at http://yihui.name/knitr/. The Bioconductor project provides style guidelines for HTML documents at http://www.bioconductor.org/packages/release/bioc/vignettes/BiocStyle/inst/doc/HtmlStyle.html.
Rqc result data is available inside template files through rqcResultSet
object.
This object is a list of summarized statistics about input files and it is used by all accessor functions and plots provided by Rqc package.
The rqcReport
function takes template file path as argument and generates personalized reports.
rqcReport(qa, templateFile = "custom_report.Rmd")
The Rqc package provides a simple interface to generate plots often used for quality assessment of high-throughput sequencing data. It uses the standard Bioconductor parallelization framework to add efficiency to data processing. The images produced by the package are high-quality figures that can be directly used on publications.
## R version 3.6.1 (2019-07-05)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.10-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 parallel stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] Rqc_1.20.0 ggplot2_3.2.1
## [3] ShortRead_1.44.0 GenomicAlignments_1.22.0
## [5] SummarizedExperiment_1.16.0 DelayedArray_0.12.0
## [7] matrixStats_0.55.0 Biobase_2.46.0
## [9] Rsamtools_2.2.0 GenomicRanges_1.38.0
## [11] GenomeInfoDb_1.22.0 Biostrings_2.54.0
## [13] XVector_0.26.0 IRanges_2.20.0
## [15] S4Vectors_0.24.0 BiocGenerics_0.32.0
## [17] BiocParallel_1.20.0 BiocStyle_2.14.0
##
## loaded via a namespace (and not attached):
## [1] colorspace_1.4-1 hwriter_1.3.2
## [3] biovizBase_1.34.0 htmlTable_1.13.2
## [5] markdown_1.1 base64enc_0.1-3
## [7] dichromat_2.0-0 rstudioapi_0.10
## [9] bit64_0.9-7 AnnotationDbi_1.48.0
## [11] splines_3.6.1 knitr_1.25
## [13] zeallot_0.1.0 Formula_1.2-3
## [15] cluster_2.1.0 dbplyr_1.4.2
## [17] shiny_1.4.0 BiocManager_1.30.9
## [19] compiler_3.6.1 httr_1.4.1
## [21] backports_1.1.5 fastmap_1.0.1
## [23] assertthat_0.2.1 Matrix_1.2-17
## [25] lazyeval_0.2.2 later_1.0.0
## [27] acepack_1.4.1 htmltools_0.4.0
## [29] prettyunits_1.0.2 tools_3.6.1
## [31] gtable_0.3.0 glue_1.3.1
## [33] GenomeInfoDbData_1.2.2 reshape2_1.4.3
## [35] dplyr_0.8.3 rappdirs_0.3.1
## [37] Rcpp_1.0.2 vctrs_0.2.0
## [39] rtracklayer_1.46.0 xfun_0.10
## [41] stringr_1.4.0 mime_0.7
## [43] ensembldb_2.10.0 XML_3.98-1.20
## [45] zlibbioc_1.32.0 scales_1.0.0
## [47] BSgenome_1.54.0 VariantAnnotation_1.32.0
## [49] hms_0.5.1 promises_1.1.0
## [51] ProtGenerics_1.18.0 AnnotationFilter_1.10.0
## [53] RColorBrewer_1.1-2 yaml_2.2.0
## [55] curl_4.2 memoise_1.1.0
## [57] gridExtra_2.3 biomaRt_2.42.0
## [59] rpart_4.1-15 latticeExtra_0.6-28
## [61] stringi_1.4.3 RSQLite_2.1.2
## [63] highr_0.8 checkmate_1.9.4
## [65] GenomicFeatures_1.38.0 rlang_0.4.1
## [67] pkgconfig_2.0.3 GenomicFiles_1.22.0
## [69] bitops_1.0-6 evaluate_0.14
## [71] lattice_0.20-38 purrr_0.3.3
## [73] labeling_0.3 htmlwidgets_1.5.1
## [75] bit_1.1-14 tidyselect_0.2.5
## [77] plyr_1.8.4 magrittr_1.5
## [79] bookdown_0.14 R6_2.4.0
## [81] Hmisc_4.2-0 DBI_1.0.0
## [83] pillar_1.4.2 foreign_0.8-72
## [85] withr_2.1.2 survival_2.44-1.1
## [87] RCurl_1.95-4.12 nnet_7.3-12
## [89] tibble_2.1.3 crayon_1.3.4
## [91] BiocFileCache_1.10.0 rmarkdown_1.16
## [93] progress_1.2.2 grid_3.6.1
## [95] data.table_1.12.6 blob_1.2.0
## [97] digest_0.6.22 xtable_1.8-4
## [99] httpuv_1.5.2 openssl_1.4.1
## [101] munsell_0.5.0 askpass_1.1