| trimLRPatterns {Biostrings} | R Documentation |
The trimLRPatterns function trims left and/or right flanking patterns
from sequences.
trimLRPatterns(Lpattern = "", Rpattern = "", subject,
max.Lmismatch = 0, max.Rmismatch = 0,
with.Lindels = FALSE, with.Rindels = FALSE,
Lfixed = TRUE, Rfixed = TRUE, ranges = FALSE)
Lpattern |
The left pattern. |
Rpattern |
The right pattern. |
subject |
An XString object, XStringSet object, or character vector containing the target sequence(s). |
max.Lmismatch |
Either an integer vector of length nLp = nchar(Lpattern)
representing an absolute number of mismatches (or edit distance if
with.Lindels is TRUE) or a single numeric value in the
interval [0, 1) representing a mismatch rate when aligning
terminal substrings (suffixes) of Lpattern with the beginning
(prefix) of subject following the conventions set by
neditStartingAt, isMatchingStartingAt,
etc.
When
Otherwise,
Once the integer vector is constructed using the rules given above, when
For a given element |
max.Rmismatch |
Same as max.Lmismatch but with Rpattern, along with
with.Rindels (below), and its initial segments (prefixes)
substring(Rpattern, 1, i).
For a given element |
with.Lindels |
If TRUE, indels are allowed in the alignments of the suffixes
of Lpattern with the subject, at its beginning.
See the with.indels arguments of the matchPattern
and neditStartingAt functions for detailed information.
|
with.Rindels |
Same as with.Lindels but for alignments of the prefixes of
Rpattern with the subject, at its end.
See the with.indels arguments of the matchPattern
and neditEndingAt functions for detailed information.
|
Lfixed, Rfixed |
Whether IUPAC extended letters in the left or right pattern should
be interpreted as ambiguities (see ?`lowlevel-matching`
for the details).
|
ranges |
If TRUE, then return the ranges to use to trim subject.
If FALSE, then returned the trimmed subject.
|
A new XString object, XStringSet object, or character vector with the "longest" flanking matches removed, as described above.
P. Aboyoun and H. Jaffee
matchPattern,
matchLRPatterns,
lowlevel-matching,
XString-class,
XStringSet-class
Lpattern <- "TTCTGCTTG"
Rpattern <- "GATCGGAAG"
subject <- DNAString("TTCTGCTTGACGTGATCGGA")
subjectSet <- DNAStringSet(c("TGCTTGACGGCAGATCGG", "TTCTGCTTGGATCGGAAG"))
## Only allow for perfect matches on the flanks
trimLRPatterns(Lpattern = Lpattern, subject = subject)
trimLRPatterns(Rpattern = Rpattern, subject = subject)
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subjectSet)
## Allow for perfect matches on the flanking overlaps
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subjectSet,
max.Lmismatch = 0, max.Rmismatch = 0)
## Allow for mismatches on the flanks
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subject,
max.Lmismatch = 0.2, max.Rmismatch = 0.2)
maxMismatches <- as.integer(0.2 * 1:9)
maxMismatches
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subjectSet,
max.Lmismatch = maxMismatches, max.Rmismatch = maxMismatches)
## Produce ranges that can be an input into other functions
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subjectSet,
max.Lmismatch = 0, max.Rmismatch = 0, ranges = TRUE)
trimLRPatterns(Lpattern = Lpattern, Rpattern = Rpattern, subject = subject,
max.Lmismatch = 0.2, max.Rmismatch = 0.2, ranges = TRUE)