fn() globalAlignmentComputes the best global pairwise alignment.
Defined in | <seqan/align.h> |
---|---|
Signature |
TScoreVal globalAlignment(align, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag]);
TScoreVal globalAlignment(gapsH, gapsV, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag]);
TScoreVal globalAlignment(frags, strings, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag]);
TScoreVal globalAlignment(alignGraph, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag]);
|
Parameters
align
|
The Align object to use for storing the pairwise alignment. |
---|---|
gapsH
|
The Gaps object for the first row (horizontal in the DP matrix). |
gapsV
|
The Gaps object for the second row (vertical in the DP matrix). |
frags
|
String of Fragment objects to store alignment in. |
strings
|
StringSet of length two with the strings to align. |
alignGraph
|
Alignment Graph for the resulting alignment. Must be initialized with two strings. |
scoringScheme
|
The scoring scheme to use for the alignment. Note that the user is responsible for ensuring that the scoring scheme is compatible with algorithmTag. |
alignConfig
|
AlignConfig instance to use for the alignment configuration. |
lowerDiag
|
Optional lower diagonal (int). |
upperDiag
|
Optional upper diagonal (int). |
algorithmTag
|
Tag to select the alignment algorithm (see AlignmentAlgorithmTags). |
Return Values
TScoreVal |
Score value of the resulting alignment (Metafunction: Value of the type of scoringScheme). |
---|
Detailed Description
There exist multiple overloads for this function with four configuration dimensions.
First, you can select whether begin and end gaps are free in either sequence using alignConfig.
Second, you can select the type of the target storing the alignment. This can be either an Align object, two Gaps objects, a AlignmentGraph, or a string of Fragment objects. Align objects provide an interface to tabular alignments with the restriction of all rows having the same type. Using two Gaps objects has the advantage that you an align sequences with different types, for example DnaString and Dna5String. Alignment Graphs provide a graph-based representation of segment-based colinear alignments. Using Fragment strings is useful for collecting many pairwise alignments, for example in the construction of Alignment Graphs for multiple-sequence alignments (MSA).
Third, you can optionally give a band for the alignment using lowerDiag and upperDiag. The center diagonal has index 0, the ith diagonal below has index -i, the ith above has index i.
Fourth, you can select the algorithm to use with algorithmTag. This can be one of NeedlemanWunsch and Gotoh. The Needleman-Wunsch algorithm supports scoring schemes with linear gap costs only while Gotoh's algorithm also allows affine gap costs.
The available alignment algorithms all have some restrictions. Gotoh's algorithm can handle arbitrary substitution and affine gap scores. Needleman-Wunsch is limited to linear gap scores. The implementation of Hirschberg's algorithm is further limited that it does not support alignConfig objects or banding. The implementation of the Myers-Hirschberg algorithm further limits this to only support edit distance (as scores, matches are scored with 0, mismatches are scored with -1).
The examples below show some common use cases.
Examples
Global alignment of two sequences using an Align object and the Needleman-Wunsch algorithm.
#include <iostream> #include <seqan/basic.h> #include <seqan/align.h> #include <seqan/stream.h> // for printint strings using namespace seqan; int main() { Dna5String seqH = "CGATT"; Dna5String seqV = "CGAAATT"; Align<Dna5String> align; resize(rows(align), 2); assignSource(row(align, 0), seqH); assignSource(row(align, 1), seqV); Score<int, Simple> scoringScheme(2, -1, -2, -1); AlignConfig<> alignConfig; int result = globalAlignment(align, scoringScheme, alignConfig); std::cout << "Score: " << result << "\n"; std::cout << "The resulting alignment is\n" << align << "\n"; return 0; }
The output is as follows:
Score: 8 The resulting alignment is 0 . CG--ATT || ||| CGAAATT
Global banded alignment of two sequences using two Gaps objects and the Gotoh algorithm.
#include <iostream> #include <seqan/basic.h> #include <seqan/align.h> #include <seqan/sequence.h> #include <seqan/stream.h> // for printint strings using namespace seqan; int main() { Dna5String seqH = "CGATT"; Dna5String seqV = "CGAAATT"; Align<Dna5String> align; resize(rows(align), 2); assignSource(row(align, 0), seqH); assignSource(row(align, 1), seqV); Score<int, Simple> scoringScheme(2, -1, -2); AlignConfig<> alignConfig; int lDiag = -2; int uDiag = 2; int result = globalAlignment(align, scoringScheme, alignConfig, lDiag, uDiag); std::cout << "Score: " << result << "\n"; std::cout << "The resulting alignment is\n" << align << "\n"; return 0; }
The output is as follows:
Score: 6 The resulting alignment is 0 . CG--ATT || ||| CGAAATT
http://seqan.readthedocs.org/en/develop/Tutorial/PairwiseSequenceAlignment.html
References
- Needleman SB, Wunsch CD: A general method applicable to the search for similarities in the amino acid sequence of two proteins. J Mol Biol 1970, 48(3): 443-53.
- Gotoh O: An improved algorithm for matching biological sequences. J Mol Biol 1982, 162(3):705-8