checkMatch | R Documentation |
checkMatch()
investigates how well oligos or assays match with their
intended target sequences
within a multiple DNA sequence alignment.
checkMatch(x, target)
## S4 method for signature 'RprimerOligo'
checkMatch(x, target)
## S4 method for signature 'RprimerAssay'
checkMatch(x, target)
x |
An |
target |
A |
The output provides information on the proportion and names of target
sequences
that match perfectly as well as with one, two, three, or four or more
mismatches to the oligo within the intended oligo binding region
in the input alignment (on-target match).
It also gives the proportion and names of target sequences that
match with a maximum of two mismatches to at least one sequence variant of
the oligo outside the
intended oligo binding region (off-target match).
The function is a wrapper to Biostrings::vcountPDict()
(Pages et al., 2020).
An RprimerMatchOligo
or RprimerMatchAssay
object, depending
on whether an RprimerOligo
or RprimerAssay
object was used
as input.
RprimerMatchOligo
objects contain the following information:
The oligo sequence in IUPAC format.
Proportion of target sequences that matches perfectly to the oligo within the intended binding region.
Names of all sequences that matches perfectly.
Proportion of target sequences with one mismatch to the oligo within the intended binding region.
Names of all sequences that matches with one mismatch.
Proportion of target sequences with two mismatches to the oligo within the intended binding region.
Names of all sequences that matches with two mismatches.
Proportion of target sequences with three mismatches to the oligo within the intended binding region.
Names of all sequences that matches with three mismatches.
Proportion of target sequences with four or more mismatches to the oligo within the intended binding region.
Names of all sequences that matches with four or more mismatches.
Proportion of target sequences with maximum two mismatches to at least one site outside the intended oligo binding region in the input alignment.
Names of all off-target matching sequences.
RprimerMatchAssay
objects contain the following information:
The forward primer sequence in IUPAC format.
Proportion of target sequences that matches perfectly with the forward primer withing the intended binding region.
Names of all sequences that matches perfectly.
Proportion of target sequences with one mismatch to the forward primer within the intended binding region.
Names of all sequences that matches with one mismatch.
Proportion of target sequences with two mismatches to the forward primer within the intended binding region.
Names of all sequences that matches with two mismatches.
Proportion of target sequences with three mismatches to the forward primer within the intended binding region.
Names of all sequences that matches with three mismatches.
Proportion of target sequences with four or more mismatches to the forward primer within the intended binding region.
Names of all sequences that matches with four or more mismatches.
Proportion of target sequences with maximum two mismatches to at least one site outside the intended forward primer binding region in the input alignment.
Names of all off-target matching sequences.
The reverse primer sequence in IUPAC format.
Proportion of target sequences that matches perfectly with the reverse primer withing the intended binding region.
Names of all sequences that matches perfectly.
Proportion of target sequences with one mismatch to the reverse primer within the intended binding region.
Names of all sequences that matches with one mismatch.
Proportion of target sequences with two mismatches to the reverse primer within the intended binding region.
Names of all sequences that matches with two mismatches.
Proportion of target sequences with three mismatches to the reverse primer within the intended binding region.
Names of all sequences that matches with three mismatches.
Proportion of target sequences with four or more mismatches to the reverse primer within the intended binding region.
Names of all sequences that matches with four or more mismatches.
Proportion of target sequences with maximum two mismatches to at least one site outside the intended reverse primer binding region in the input alignment.
Names of all off-target matching sequences.
If the input assay contains probes, the following information is also added:
The probe sequence in IUPAC format.
Proportion of target sequences that matches perfectly with the probe withing the intended binding region.
Names of all sequences that matches perfectly.
Proportion of target sequences with one mismatch to the probe within the intended binding region.
Names of all sequences that matches with one mismatch.
Proportion of target sequences with two mismatches to the probe within the intended binding region.
Names of all sequences that matches with two mismatches.
Proportion of target sequences with three mismatches to the probe within the intended binding region.
Names of all sequences that matches with three mismatches.
Proportion of target sequences with four or more mismatches to the probe within the intended binding region.
Names of all sequences that matches with four or more mismatches.
Proportion of target sequences with maximum two mismatches to at least one site outside the intended probe binding region in the input alignment.
Names of all off-target matching sequences.
checkMatch(RprimerOligo)
:
checkMatch(RprimerAssay)
:
There are a few limitations with this function, which is important to be aware of:
False negatives or positives may occur due to poorly aligned sequences
The output does not tell which strand (minus or plus) the oligo matches to. This is important to consider when assessing off-target matches to single-stranded targets
Ambiguous bases and gaps in the target sequences are identified as mismatches
The function checks strictly on- and off-target, and may therefore miss off-target matches that partially overlap the intended target
Pages, H., Aboyoun, P., Gentleman R., and DebRoy S. (2020). Biostrings: Efficient manipulation of biological strings. R package version 2.57.2.
#### RprimerOligo objects
data("exampleRprimerOligo")
data("exampleRprimerAlignment")
x <- exampleRprimerOligo[1:2, ]
target <- exampleRprimerAlignment
checkMatch(x, target)
#### RprimerAssay objects
data("exampleRprimerAssay")
data("exampleRprimerAlignment")
x <- exampleRprimerAssay[1:2, ]
target <- exampleRprimerAlignment
checkMatch(x, target)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.