Function
globalAlignment
Computes the best global pairwise alignment.
globalAlignment(align,          scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag])
globalAlignment(gapsH, gapsV,   scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag])
globalAlignment(frags, strings, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag])
globalAlignment(alignmentGraph, scoringScheme, [alignConfig,] [lowerDiag, upperDiag,] [algorithmTag])
Include Headers
seqan/align.h
Parameters
align
An Align object that stores the alignment. The number of rows must be 2 and the sequences must have already been set. row(align, 0) is the horizontal one in the alignment matrix alignment, row(align, 1) is the vertical one.
Types: Align
gapsH
Horizontal gapped sequence in alignment matrix.
Types: Gaps
gapsV
Vertical gapped sequence in alignment matrix.
Types: Gaps
frags
String of Fragment objects. The sequence with id 0 is the horizontal one, the sequence with id 1 is the vertical one.
alignmentGraph
Alignment Graph object to store the alignment in.
Remarks: The underlying StringSet must be an Owner StringSet.
strings
A StringSet containing two sequences.
Types: StringSet
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.
Types: Score
alignConfig
The AlignConfig to use for the alignment.
lowerDiag
Optional lower diagonal.
Types: int
upperDiag
Optional upper diagonal.
Types: int
algorithmTag
The Tag for picking the alignment algorithm.
Remarks
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 Alignment Graph, 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.
Return Values
An integer with the alignment score, as given by the Value metafunction of the Score type.
Examples
Global alignment of two sequences using an Align object and the Needleman-Wunsch algorithm.
Dna5String seqH = "CGATT";
Dna5String seqV = "CGAAATT";
 
Align<Dna5String> align;
resize(rows(align), 2);
assignSource(row(align, 0), seqH);
assignSource(row(align, 0), seqV);
Score<int, Simple> scoringScheme(2, -1, -2);
AlignConfig<> alignConfig;
 
int result = globalAlignment(align, scoringScheme, alignConfig,
                             NeedlemanWunsch());
Global banded alignment of two sequences using two Gaps objects and the Gotoh algorithm.
Dna5String seqH = "CGATT";
Gaps<Dna5String, ArrayGaps> gapsH(seqH);
DnaString seqV = "CGAAATT";
Gaps<Dna5String, AnchorGaps<> > gapsV(seqV);
 
Score<int, Simple> scoringScheme(5, -3, -1, -5);
AlignConfig<> alignConfig;
 
int result = globalAlignment(gapsH, gapsV, scoringScheme, alignConfig, -2, 2);
Example Programs
SeqAn - Sequence Analysis Library - www.seqan.de
 

Page built @2013/07/11 09:12:17