alignSeq | R Documentation |
Perform multiple sequence alignment using one of three methods and output results to the console or as a pdf file. One may perform the alignment of all amino acid or nucleotide sequences in a single repertoire_id. Alternatively, one may search for a given sequence within a list of samples using an edit distance threshold.
alignSeq(
study_table,
repertoire_ids = NULL,
sequence_list = NULL,
edit_distance = 15,
type = "junction",
method = "ClustalOmega",
top = 150
)
study_table |
A tibble consisting of antigen receptor sequences
imported by the LymphoSeq2 function |
repertoire_ids |
A character vector indicating the name of the repertoire_id in the productive sequence list. |
sequence_list |
A character vector of one ore more amino acid or nucleotide CDR3 sequences to search. |
edit_distance |
An integer giving the minimum edit distance that the sequence must be less than or equal to. See details below. |
type |
A character vector indicating whether "junction_aa" or "junction"
sequences should be aligned. If "junction_aa" is specified, then run
|
method |
A character vector indicating the multiple sequence alignment method to be used. Refer to the Bioconductor "msa" package for more details. Options include "ClustalW", "ClustalOmega", and "Muscle". |
top |
The number of top sequences to perform alignment. |
Edit distance is a way of quantifying how dissimilar two sequences are to one another by counting the minimum number of operations required to transform one sequence into the other. For example, an edit distance of 0 means the sequences are identical and an edit distance of 1 indicates that the sequences different by a single amino acid or junction.
Performs a multiple sequence alignment object.
If having trouble saving pdf files, refer to Bioconductor package msa for installation instructions http://bioconductor.org/packages/release/bioc/vignettes/msa/inst/doc/msa.pdf
file_path <- system.file("extdata", "IGH_sequencing", package = "LymphoSeq2")
study_table <- LymphoSeq2::readImmunoSeq(path = file_path, threads = 1)
study_table <- LymphoSeq2::topSeqs(study_table, top = 100)
nucleotide_table <- LymphoSeq2::productiveSeq(study_table, aggregate = "junction")
LymphoSeq2::alignSeq(nucleotide_table,
repertoire_ids = "IGH_MVQ92552A_BL", type = "junction",
method = "ClustalW"
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.