filterAndTrim: Filter and trim fastq file(s).

Description Usage Arguments Details Value See Also Examples

View source: R/filter.R

Description

Filters and trims an input fastq file(s) (can be compressed) based on several user-definable criteria, and outputs fastq file(s) (compressed by default) containing those trimmed reads which passed the filters. Corresponding forward and reverse fastq file(s) can be provided as input, in which case filtering is performed on the forward and reverse reads independently, and both reads must pass for the read pair to be output.

Usage

1
2
3
4
5
filterAndTrim(fwd, filt, rev = NULL, filt.rev = NULL, compress = TRUE,
  truncQ = 2, truncLen = 0, trimLeft = 0, maxLen = Inf, minLen = 20,
  maxN = 0, minQ = 0, maxEE = Inf, rm.phix = TRUE, primer.fwd = NULL,
  matchIDs = FALSE, id.sep = "\\s", id.field = NULL,
  multithread = FALSE, n = 1e+05, OMP = TRUE, verbose = FALSE)

Arguments

fwd

(Required). character. The path(s) to the input fastq file(s).

filt

(Required). character. The path(s) to the output filtered file(s) corresponding to the fwd input files. If containing directory does not exist, it will be created.

rev

(Optional). Default NULL. The path(s) to the input reverse fastq file(s) from paired-end sequence data corresponding to those provided to the fwd argument. If NULL, the fwd files are processed as single-reads.

filt.rev

(Optional). Default NULL, but required if rev is provided. The path(s) to the output fastq file(s) corresponding to the rev input. Can also provide a directory, which if not existing will be created (how to differentiate between dir/file in len1 case?).

compress

(Optional). Default TRUE. If TRUE, the output fastq file(s) are gzipped.

FILTERING AND TRIMMING PARAMETERS ———

Note: When filtering paired reads... If a length 1 vector is provided, the same parameter value is used for the forward and reverse reads. If a length 2 vector is provided, the first value is used for the forward reads, and the second for the reverse reads.

truncQ

(Optional). Default 2. Truncate reads at the first instance of a quality score less than or equal to truncQ.

truncLen

(Optional). Default 0 (no truncation). Truncate reads after truncLen bases. Reads shorter than this are discarded.

trimLeft

(Optional). Default 0. The number of nucleotides to remove from the start of each read. If both truncLen and trimLeft are provided, filtered reads will have length truncLen-trimLeft.

maxLen

(Optional). Default Inf (no maximum). Remove reads with length greater than maxLen. maxLen is enforced before trimming and truncation.

minLen

(Optional). Default 20. Remove reads with length less than minLen. minLen is enforced after trimming and truncation.

maxN

(Optional). Default 0. After truncation, sequences with more than maxN Ns will be discarded. Note that dada does not allow Ns.

minQ

(Optional). Default 0. After truncation, reads contain a quality score less than minQ will be discarded.

maxEE

(Optional). Default Inf (no EE filtering). After truncation, reads with higher than maxEE "expected errors" will be discarded. Expected errors are calculated from the nominal definition of the quality score: EE = sum(10^(-Q/10))

rm.phix

(Optional). Default TRUE. If TRUE, discard reads that match against the phiX genome, as determined by isPhiX.

primer.fwd

(Optional). Default NULL. Paired-read filtering only. A character string defining the forward primer. Only allows unambiguous nucleotides. The primer will be compared to the first len(primer.fwd) nucleotides at the start of the read. If there is not an exact match, the read is filtered out. For paired reads, the reverse read is also interrogated, and if the primer is detected on the reverse read, the forward/reverse reads are swapped.

matchIDs

(Optional). Default FALSE. Paired-read filtering only. Whether to enforce matching between the id-line sequence identifiers of the forward and reverse fastq files. If TRUE, only paired reads that share id fields (see below) are output. If FALSE, no read ID checking is done. Note: matchIDs=FALSE essentially assumes matching order between forward and reverse reads. If that matched order is not present future processing steps may break (in particular mergePairs).

id.sep

(Optional). Default "\s" (white-space). Paired-read filtering only. The separator between fields in the id-line of the input fastq files. Passed to the strsplit.

id.field

(Optional). Default NULL (automatic detection). Paired-read filtering only. The field of the id-line containing the sequence identifier. If NULL (the default) and matchIDs is TRUE, the function attempts to automatically detect the sequence identifier field under the assumption of Illumina formatted output.

multithread

(Optional). Default is FALSE. If TRUE, input files are filtered in parallel via mclapply. If an integer is provided, it is passed to the mc.cores argument of mclapply. Note that the parallelization here is by forking, and each process is loading another fastq file into memory. If memory is an issue, execute in a clean environment and reduce the chunk size n and/or the number of threads.

n

(Optional). Default 1e5. The number of records (reads) to read in and filter at any one time. This controls the peak memory requirement so that very large fastq files are supported. See FastqStreamer for details.

OMP

(Optional). Default TRUE. Whether or not to use OMP multithreading when calling FastqStreamer. Should be set to FALSE if calling this function within a parallelized chunk of code. If multithread=TRUE, this argument will be coerced to FALSE.

verbose

(Optional). Default FALSE. Whether to output status messages.

Details

filterAndTrim is a multithreaded convenience interface for the fastqFilter and fastqPairedFilter filtering functions. Note that error messages and tracking are not handled gracefully when using the multithreading functionality. If errors arise, it is recommended to re-run without multithreading to troubleshoot the issue.

Value

Integer matrix. Returned invisibly (i.e. only if assigned to something). Rows correspond to the input files, columns record the reads.in and reads.out after filtering.

See Also

fastqFilter fastqPairedFilter FastqStreamer

Examples

1
2
3
4
5
testFastqs = c(system.file("extdata", "sam1F.fastq.gz", package="dada2"), 
               system.file("extdata", "sam2F.fastq.gz", package="dada2"))
filtFastqs <- c(tempfile(fileext=".fastq.gz"), tempfile(fileext=".fastq.gz"))
filterAndTrim(testFastqs, filtFastqs, maxN=0, maxEE=2, verbose=TRUE)
filterAndTrim(testFastqs, filtFastqs, truncQ=2, truncLen=200, rm.phix=TRUE)

Bioconductor-mirror/dada2 documentation built on June 1, 2017, 7:31 a.m.