Package ‘basecallQC’
May 29, 2024
Type Package
Title Working with Illumina Basecalling and Demultiplexing input and
output files
Version 1.28.0
Author Thomas Carroll and Marian Dore
Maintainer Thomas Carroll <[email protected]>
Description
The basecallQC package provides tools to work with Illumina bcl2Fastq (versions >= 2.1.7) soft-
ware.Prior to basecalling and demultiplexing using the bcl2Fastq software, basecallQC func-
tions allow the user to update Illumina sample sheets from versions <= 1.8.9 to >= 2.1.7 stan-
dards, clean sample sheets of common problems such as invalid sample names and IDs, cre-
ate read and index basemasks and the bcl2Fastq command. Following the generation of base-
called and demultiplexed data, the basecallQC packages allows the user to generate HTML ta-
bles, plots and a self contained report of summary metrics from Illumina XML output files.
LazyData TRUE
biocViews Sequencing, Infrastructure, DataImport, QualityControl
License GPL (>= 3)
Depends R (>= 3.4), stats, utils, methods, rmarkdown, knitr,
prettydoc, yaml
Imports ggplot2, stringr, XML, raster, dplyr, data.table, tidyr,
magrittr, DT, lazyeval, ShortRead
Suggests testthat, BiocStyle
VignetteBuilder knitr
SystemRequirements bcl2Fastq (versions >= 2.1.7)
Collate processIlluminaSamplesheets_Functions.R
processXMLs_Functions.R allClasses.R allMethods.R
processExternalFormats_Functions.R plots.R reporting.R tables.R
FastQCShortRead_Functions.R processInterOps_Functions.R
RoxygenNote 7.1.0
git_url https://git.bioconductor.org/packages/basecallQC
git_branch RELEASE_3_19
1
2 baseCallMetrics
git_last_commit 536aad0
git_last_commit_date 2024-04-30
Repository Bioconductor 3.19
Date/Publication 2024-05-29
Contents
baseCallMetrics . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 2
basecallQC-class . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 3
BCL2FastQparams-class . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 4
createBasemasks . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 5
createBCLcommand . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 6
demultiplexMetrics . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 7
demuxBarplot . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 8
demuxBoxplot . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
indexlengths . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
interOpsReport . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 10
makeFQTable . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 12
passFilterBar . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 13
passFilterBoxplot . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 14
passFilterTilePlot . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 15
readlengths . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 16
reportBCL . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 16
summaryConvStatsTable . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 17
summaryDemuxTable . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 18
validateBCLSheet . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 19
Index 21
baseCallMetrics Gather basecalling metrics from a Run (using Run’s Conversion-
Stats.xml file).
Description
Gather basecalling metrics from a Run (using Run’s ConversionStats.xml file).
Usage
baseCallMetrics(bcl2fastqparams)
Arguments
bcl2fastqparams
A BCL2FastQparams object as created by BCL2FastQparams() constructor.
basecallQC-class 3
Value
A list of length two containing the full basecalling metrics from a Run (using Run’s Conversion-
Stats.xml file). Contains an unsummarised data.frame and basecalling metrics summarised to Sam-
ple, Lane, Sample by lane, and Sample by Lane and Tile.
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
convMetrics <- baseCallMetrics(bcl2fastqparams)
basecallQC-class The basecallQC object and constructor.
Description
Object and method to handle Illumina basecalling/demultiplexing inputs and output files. Provides
sample sheet cleanup, basecall command and summary QC statistics for basecalling/demultiplexing.
The basecallQC object and constructor.
Usage
basecallQC(bcl2fastqparams, RunMetaData = NULL, sampleSheet = NULL,
doFQMetric = FALSE)
Arguments
bcl2fastqparams
A BCL2FastQparams object as created by BCL2FastQparams() constructor.
RunMetaData Any run metadata to attach (data.frame)
sampleSheet A sample sheet for Illumina basecalling using bcl2Fastq (See vignette for more
details).
doFQMetric TRUE or FALSE. Perform ShortRead FastQ quality assessment using Short-
Read’s qa and report function
Details
The basecallQC object contains slots BCL2FastQparams, cleanedSampleSheet, baseMasks, BCLCom-
mand, baseCallMetrics, demultiplexMetrics and fqQCmetrics.
"BCL2FastQparams" A BCL2FastQparams object
"cleanedSampleSheet" A data.frame containing the cleaned sample sheet for Illumina base-
calling using bcl2Fastq versions >= 2.1.7
4 BCL2FastQparams-class
"baseMasks" A data.frame containing basecall masks per lane for use with bcl2Fastq versions
>= 2.1.7. Basemasks in data.frame for reads and indexes as well as the total basemasks for
each lane.
"BCLCommand" A character string containing the command to be used for basecalling using
bcl2Fastq (versions >= 2.1.7).
"baseCallMetrics" A list containing the full basecalling metrics from ConversionStats.xml.
Contains an unsummarised data.frame and basecalling metrics summarised to Sample, Lane,
Sample by lane, and Sample by Lane and Tile
"demultiplexMetrics" A list containing the full demultiplexing metrics from Demultiplex-
ingStats.xml. Contains an unsummarised data.frame and demultiplexing metrics filtered to
per Sample metrics
"fqQCmetrics" A list containing a data.frame of read counts and links to ShortRead QA reports
and a ShortRead QA object containing quality information for generated fastQs.
Value
basecallQC a basecallQC object (See details for more information)
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
BCL2FastQparams-class The Parameters for BCL2FastQparams object.
Description
Parameter class and accessors for use with basecallQC
Parameter class and accessors
Usage
BCL2FastQparams(runXML = NULL, config = NULL, runDir = NULL,
outDir = NULL, verbose = TRUE)
Arguments
runXML file path to runParameters.xml ,if not specified looks in top level of run directory.
config file path to config.ini ,if not specified looks in top level of run directory.
runDir file path to run directory.
outDir file path to out directory.
verbose TRUE or FALSE. Messages on or off. Warnings/errors persist
createBasemasks 5
Details
The BCL2FastQparams object contains slots RunDir, OutDir and RunParameters
"RunDir" Character string specifying the top level Run directory
"OutDir" Character string specifying the output directory
"RunParameters" A data.frame containing the information from runParameters.xml (See vi-
gnette for more details).
Value
A BCL2FastQparams object (See details).
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
createBasemasks Function to create basemasks for basecalling from Illumina sam-
plesheet (for bcl2Fastq versions >= 2.1.7).
Description
Parses the Illumina sample sheet for versions >= 2.1.7 and creates basemasks.
Usage
createBasemasks(cleanedSampleSheet, param)
Arguments
cleanedSampleSheet
Data.frame of cleaned samplesheet for Illumina basecalling using bcl2Fastq ver-
sions >= 2.1.7 (see vignette for more details)
param A BCL2FastQparams object
Value
A data.frame containing basecall masks per lane for reads and indexes as well as per lane complete
basemasks.
Author(s)
Thomas Carroll and Marian Dore
6 createBCLcommand
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
cleanedSampleSheet <- validateBCLSheet(sampleSheet,param=bcl2fastqparams)
basemasks <- createBasemasks(cleanedSampleSheet,param=bcl2fastqparams)
createBCLcommand Function to create command for Illumina basecalling/demultiplexing
using bcl2fastq versions >= 2.1.7.
Description
Creates the command to be used for basecalling/demultiplexing with bcl2fastq versions >= 2.1.7
Usage
createBCLcommand(bcl2fastqparams, cleanedSampleSheet, baseMasks)
Arguments
bcl2fastqparams
A BCL2FastQparams object.
cleanedSampleSheet
Data.frame of cleaned samplesheet for Illumina basecalling/demultiplexing us-
ing bcl2fastq versions >= 2.1.7 (see vignette for more details)
baseMasks A data.frame of basemasks as created by createBasemasks() function
Value
A character vector containing the command for Illumina basecalling using bcl2fastq versions >=
2.1.7
Author(s)
Thomas Carroll and Marian Dore
demultiplexMetrics 7
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
cleanedSampleSheet <- validateBCLSheet(sampleSheet,param=bcl2fastqparams)
baseMasks <- createBasemasks(cleanedSampleSheet,param=bcl2fastqparams)
toSubmit <- createBCLcommand(bcl2fastqparams,cleanedSampleSheet,baseMasks)
demultiplexMetrics Gather demultiplexing metrics from a Run (using Run’s Demultiplex-
ingStats.xml file).
Description
Gather demultiplexing metrics from a Run (using Run’s DemultiplexingStats.xml file).
Usage
demultiplexMetrics(bcl2fastqparams)
Arguments
bcl2fastqparams
A BCL2FastQparams object as created by BCL2FastQparams() constructor.
Value
A list of length two containing the full demultiplexing metrics from a Run (using Run’s Demulti-
plexingStats.xml file). Contains an unsummarised data.frame and demultiplexing metrics filtered to
per Sample metrics
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
demuxMetrics <- demultiplexMetrics(bcl2fastqparams)
8 demuxBarplot
demuxBarplot Barplot of Illumina demultiplexing statistics.
Description
Produces a barplot of demultiplexing statistics of reads with perfect/mismatched barcode.
Usage
## S4 method for signature 'baseCallQC'
demuxBarplot(object,groupBy)
## S4 method for signature 'basecallQC'
demuxBarplot(object = "basecallQC",
groupBy = c("Lane"))
## S4 method for signature 'list'
demuxBarplot(object = "basecallQC", groupBy = c("Lane"))
Arguments
object A basecallQC object or list from call to demultiplexMetrics()
groupBy Character vector of how data is grouped for plotting. Should be either "Project","Sample"
or "Lane".
Value
A ggplot2 object.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
plot <- demuxBarplot(bclQC)
demuxBoxplot 9
demuxBoxplot Boxplot of Illumina demultiplexing statistics.
Description
Produces a boxplot of demultiplexing statistics of reads with perfect/mismatched barcode.
Usage
## S4 method for signature 'baseCallQC'
demuxBoxplot(object,groupBy)
## S4 method for signature 'basecallQC'
demuxBoxplot(object = "basecallQC",
groupBy = c("Lane"))
## S4 method for signature 'list'
demuxBoxplot(object = "basecallQC", groupBy = c("Lane"))
Arguments
object A basecallQC object or list from call to demultiplexMetrics()
groupBy Character vector of how data is grouped for plotting. Should be either "Project","Sample"
or "Lane".
Value
A ggplot2 object.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
plot <- demuxBoxplot(bclQC)
10 interOpsReport
indexlengths Index lengths
Description
Index lengths as defined by runParameters.xml
Usage
## S4 method for signature 'BCL2FastQparams'
indexlengths(object)
## S4 method for signature 'BCL2FastQparams'
indexlengths(object = "BCL2FastQparams")
Arguments
object A BCL2FastQparams object
Value
Index lengths as defined runParameters.xml.
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
indexlength <- indexlengths(bcl2fastqparams)
interOpsReport Function to parse InterOps files and generate summary reports
Description
Parses the InterOps binary files produced by Illumina’s Real Time Analysis sofware and used by
Illumina’s SAV sofware. InterOp binary files contain information on phasing/prephsing, yield,read
numbers and basecalling quality score distributions per cycle. This interOpsReport functions parses
and summarises the InterOps files, TileMetrics.bin and QMetrics.bin, and the Stats directory XML
files, ConversionStats.xml and DemultiplexingStats.xml.
interOpsReport 11
Usage
interOpsReport(bcl2fastqparams, verbose = TRUE)
Arguments
bcl2fastqparams
A BCL2FastQparams object.
verbose TRUE or FALSE . TRUE reports progress through file parsing.
Details
The interOpsReport function returns a list of machine and run information, basecalling quality
information and demultiplexing information. The three named elements are descibed below.
"machineReport" A data.frame containing information machine and software parameters
"sequencingReport" A data.frame of mean cluster density, percentage clusters passing filter,
phasing and prephasing percentages, number of reads total/passing filter and percent of reads
with mean quality score > Q30 grouped by lane and read
"demuxReport" A data.frame of demultiplexing results containing yield, number of reads,
percentage of reads with quality scores greater than >Q30 and the percent of total reads per
lane. Results are summarised per lane for samples, underdetermined indexes and all indexes
(identifed and unidentified).
Value
A named list of length 3 containing machine and run information, basecalling quality information
and demultiplexing information.
Author(s)
Thomas Carroll.
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
# myRes_BCAGJ8ANXX <- interOpsReport(bcl2fastqparams,verbose=TRUE)
12 makeFQTable
makeFQTable Generate an HTML table linking to per sample summary fastq QC
statistics from ShortRead
Description
Creates an HTML table linking to per sample summary fastq QC statistics from ShortRead
Usage
makeFQTable(object, output = "static")
Arguments
object A basecall QC object as returned from basecallQC function
output Whether the report contains frozen or sortable tables. Options are "static" and
"html"
Value
A HTML table for reporting fastq QC results from ShortRead. Table contains read counts and links
to ShortRead QA reports per sample.
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
#makeFQTable(bclQC,output="static")
passFilterBar 13
passFilterBar Barplot of Illumina basecalling statistics for reads passing filter.
Description
Produces a barplot of Illumina basecalling statistics for reads passing filter.
Usage
## S4 method for signature 'baseCallQC'
passFilterBar(object,groupBy,metricToPlot)
## S4 method for signature 'basecallQC'
passFilterBar(object = "basecallQC",
groupBy = c("Lane"), metricToPlot = "Yield")
## S4 method for signature 'list'
passFilterBar(object = "basecallQC", groupBy = c("Lane"),
metricToPlot = "Yield")
Arguments
object A basecallQC object or list from call to baseCallMetrics()
groupBy Character vector of how data is grouped for plotting. Should be either "Project","Sample","Lane","Tile","ReadNumber".
metricToPlot Character vector defining which metric will be displayed in plot. Should be
either "Yield","Yield30","QualityScoreSum" or "ClusterCount".
Value
A ggplot2 object.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
plot <- passFilterBar(bclQC)
14 passFilterBoxplot
passFilterBoxplot Boxplot of Illumina basecalling statistics for reads passing filter.
Description
Produces a boxplot of basecalling statistics for reads passing filter.
Usage
## S4 method for signature 'baseCallQC'
passFilterBoxplot(object,groupBy,metricToPlot)
## S4 method for signature 'basecallQC'
passFilterBoxplot(object = "basecallQC",
groupBy = c("Lane"), metricToPlot = "Yield")
## S4 method for signature 'list'
passFilterBoxplot(object = "basecallQC",
groupBy = c("Lane"), metricToPlot = "Yield")
Arguments
object A basecallQC object or list from call to baseCallMetrics()
groupBy Character vector of how data is grouped for plotting. Should be either "Project","Sample","Lane","Tile","ReadNumber".
metricToPlot Character vector defining which metric will be displayed in plot. Should be
either "Yield","Yield30","QualityScoreSum" or "ClusterCount".
Value
A ggplot2 object.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
plot <- passFilterBoxplot(bclQC,groupBy = "Sample")
passFilterTilePlot 15
passFilterTilePlot Tile plot of Illumina basecalling statistics for reads passing filter.
Description
Produces a plot of metric per Tile for basecalling statistics of reads passing/failing filter.
Usage
## S4 method for signature 'baseCallQC'
passFilterTilePlot(object,metricToPlot)
## S4 method for signature 'basecallQC'
passFilterTilePlot(object = "basecallQC",
metricToPlot = "Yield")
## S4 method for signature 'list'
passFilterTilePlot(object = "basecallQC",
metricToPlot = "Yield")
Arguments
object A basecallQC object or list from call to baseCallMetrics()
metricToPlot Character vector defining which metric will be displayed in plot. Should be
either "Yield","Yield30","QualityScoreSum" or "ClusterCount".
Value
A ggplot2 object.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
plot <- passFilterTilePlot(bclQC,metricToPlot="Yield")
16 reportBCL
readlengths Read lengths
Description
Read lengths as defined by runParameters.xml
Usage
## S4 method for signature 'BCL2FastQparams'
readlengths(object)
## S4 method for signature 'BCL2FastQparams'
readlengths(object = "BCL2FastQparams")
Arguments
object A BCL2FastQparams object
Value
Read lengths as defined runParamaeters.xml
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
readlength <- readlengths(bcl2fastqparams)
reportBCL Generate basecallQC report
Description
Creates a summary report from basecalling and demultiplexing metrics.
summaryConvStatsTable 17
Usage
## S4 method for signature 'basecallQC'
reportBCL(object,reportOut,reportOutDir,output,reportRMDfile,FQQC)
## S4 method for signature 'basecallQC'
reportBCL(object = "basecallQC",
reportOut = "report.html", reportOutDir = getwd(), output = "static",
reportRMDfile = NULL, FQQC = FALSE)
Arguments
object A basecall QC object as returned from basecallQC() function
reportOut Name of report file
reportOutDir Directory for the report file
output Whether the report contains frozen or sortable tables. Options are "static" and
"html"
reportRMDfile RMD to be used for reporting. (Default uses standard report template)
FQQC TRUE or FALSE, whether to run ShortRead fastq QC on any fastQ in output
directory.
Value
An HTML report is written to file.
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
reportBCL(bclQC,"TestReport.html",output="html")
summaryConvStatsTable Creates an HTML table of per sample summary statistics from base-
calling results
Description
Creates an HTML table of per sample summary statistics from basecalling results
18 summaryDemuxTable
Usage
## S4 method for signature 'baseCallQC'
summaryConvStatsTable(object)
## S4 method for signature 'basecallQC'
summaryConvStatsTable(object = "basecallQC",
output = "static")
## S4 method for signature 'list'
summaryConvStatsTable(object = "basecallQC",
output = "static")
Arguments
object A basecallQC object or list from call to baseCallMetrics()
output Whether the report contains frozen or sortable tables. Options are "static" and
"html"
Value
An HTML table for reporting basecalling results.
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
summaryDemuxTable(bclQC,output="static")
summaryDemuxTable Generate an HTML table of per sample summary demultiplexing
statistics
Description
Generate an HTML table of per sample summary demultiplexing statistics
validateBCLSheet 19
Usage
## S4 method for signature 'baseCallQC'
summaryDemuxTable(object)
## S4 method for signature 'basecallQC'
summaryDemuxTable(object = "basecallQC",
output = "static")
## S4 method for signature 'list'
summaryDemuxTable(object = "basecallQC", output = "static")
Arguments
object A basecallQC object or list from call to demultiplexMetrics()
output Whether the report contains frozen or sortable tables. Options are "static" and
"html"
Value
An HTML table for reporting demultiplexing results.
Author(s)
Thomas Carroll
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
outDir <- file.path(fileLocations,"Runs/161105_D00467_0205_AC9L0AANXX/C9L0AANXX/")
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),outDir,verbose=FALSE)
bclQC <- basecallQC(bcl2fastqparams,RunMetaData=NULL,sampleSheet)
summaryDemuxTable(bclQC,output="static")
validateBCLSheet Illumina sample sheet cleaning and updating for bcl2Fastq versions
>= 2.1.7
Description
Parses an Illumina bcl2Fastq sample sheet to create a standardised/updated sample sheet for bcl2Fastq
>= Version 2.1.7
Usage
validateBCLSheet(sampleSheet, param = bcl2fastqparams)
20 validateBCLSheet
Arguments
sampleSheet File location of a sample sheet for Illumina basecalling using bcl2Fastq (See
vignette for more details).
param A BCL2FastQparams object
Value
cleanedSampleSheet A data.frame containing the cleaned sample sheet for Illumina basecalling
using bcl2Fastq versions >= 2.1.7.
Author(s)
Thomas Carroll and Marian Dore
Examples
fileLocations <- system.file("extdata",package="basecallQC")
runXML <- dir(fileLocations,pattern="runParameters.xml",full.names=TRUE)
config <- dir(fileLocations,pattern="config.ini",full.names=TRUE)
sampleSheet <- dir(fileLocations,pattern="*\\.csv",full.names=TRUE)
bcl2fastqparams <- BCL2FastQparams(runXML,config,runDir=getwd(),verbose=FALSE)
cleanedSampleSheet <- validateBCLSheet(sampleSheet,param=bcl2fastqparams)
Index
baseCallMetrics, 2
basecallQC (basecallQC-class), 3
basecallQC-basecallQC
(basecallQC-class), 3
basecallQC-class, 3
BCL2FastQparams
(BCL2FastQparams-class), 4
BCL2FastQparams-BCL2FastQparams
(BCL2FastQparams-class), 4
BCL2FastQparams-class, 4
createBasemasks, 5
createBCLcommand, 6
demultiplexMetrics, 7
demuxBarplot, 8
demuxBarplot,baseCallQC-method
(demuxBarplot), 8
demuxBarplot,basecallQC-method
(demuxBarplot), 8
demuxBarplot,list-method
(demuxBarplot), 8
demuxBarplot.basecallQC (demuxBarplot),
8
demuxBoxplot, 9
demuxBoxplot,baseCallQC-method
(demuxBoxplot), 9
demuxBoxplot,basecallQC-method
(demuxBoxplot), 9
demuxBoxplot,list-method
(demuxBoxplot), 9
demuxBoxplot.basecallQC (demuxBoxplot),
9
indexlengths, 10
indexlengths,BCL2FastQparams-method
(indexlengths), 10
indexlengths,BCL2FastQparams-method,
(indexlengths), 10
indexlengths.bcl2fastqparams
(indexlengths), 10
interOpsReport, 10
makeFQTable, 12
passFilterBar, 13
passFilterBar,baseCallQC-method
(passFilterBar), 13
passFilterBar,basecallQC-method
(passFilterBar), 13
passFilterBar,list-method
(passFilterBar), 13
passFilterBar.basecallQC
(passFilterBar), 13
passFilterBoxplot, 14
passFilterBoxplot,baseCallQC-method
(passFilterBoxplot), 14
passFilterBoxplot,basecallQC-method
(passFilterBoxplot), 14
passFilterBoxplot,list-method
(passFilterBoxplot), 14
passFilterBoxplot.basecallQC
(passFilterBoxplot), 14
passFilterTilePlot, 15
passFilterTilePlot,baseCallQC-method
(passFilterTilePlot), 15
passFilterTilePlot,basecallQC-method
(passFilterTilePlot), 15
passFilterTilePlot,list-method
(passFilterTilePlot), 15
passFilterTilePlot.basecallQC
(passFilterTilePlot), 15
readlengths, 16
readlengths,BCL2FastQparams-method
(readlengths), 16
readlengths,BCL2FastQparams-method,
(readlengths), 16
21
22 INDEX
readlengths.bcl2fastqparams
(readlengths), 16
reportBCL, 16
reportBCL,baseCallQC-method
(reportBCL), 16
reportBCL,basecallQC-method
(reportBCL), 16
reportBCL.basecallQC (reportBCL), 16
summaryConvStatsTable, 17
summaryConvStatsTable,baseCallQC-method
(summaryConvStatsTable), 17
summaryConvStatsTable,basecallQC-method
(summaryConvStatsTable), 17
summaryConvStatsTable,list-method
(summaryConvStatsTable), 17
summaryConvStatsTable.basecallQC
(summaryConvStatsTable), 17
summaryDemuxTable, 18
summaryDemuxTable,baseCallQC-method
(summaryDemuxTable), 18
summaryDemuxTable,basecallQC-method
(summaryDemuxTable), 18
summaryDemuxTable,list-method
(summaryDemuxTable), 18
summaryDemuxTable.basecallQC
(summaryDemuxTable), 18
validateBCLSheet, 19