Structural Bioinformatics Library
Template C++ / Python API for developping structural bioinformatics applications.
T_Alignment_engine_sequences_seqan< SequenceType, SeqanSequenceConverter, FreeEndsAlignment, ScoreType, SeqanUnitType, SeqanCustomMatrix, SeqanAlgorithm > Class Template Reference

Engine for making alignments between sequences using Seqan. More...

#include <Alignment_engine_sequences_seqan.hpp>

Public Types

enum  Sequence_length_type
 Enum for selecting the length of the sequence to consider when computing the similarity / identity percentages. More...
 
typedef SequenceType Sequence_or_structure
 Type for representing the sequences or the structures to align.
More...
 
typedef SequenceOrStructure::Alignment_unit Alignment_unit
 Type for a unit (e.g residue or nucleotid) More...
 
typedef SequenceOrStructure::Alignment_unit_name Alignment_unit_name
 Representation of the name of a unit. More...
 
typedef SequenceOrStructure::Alignment_unit_rep Alignment_unit_rep
 Representation of a unit for index purposes. More...
 
typedef AlignerAlgorithm::Score_type Score_type
 Representation of the score of the algorithm. More...
 
typedef std::pair< Alignment_unit_rep, Alignment_unit_repAligned_pair
 Representation of two aligned units. More...
 
typedef std::vector< Aligned_pairAlignment_type
 Representation of an alignment as a sequence of aligned units. More...
 
typedef std::pair< Alignment_unit_name, Alignment_unit_nameName_pair
 Representation of two unit names. More...
 
typedef std::map< Name_pair, double, Is_lower_name_pair > Substitution_matrix
 Representation of the storage of values of the subtitution matrix. More...
 

Public Member Functions

void align (void)
 Overwritten method for setting the seqan substitution matrix before the alignment. More...
 
void load_seqan_substitution_matrix (const std::string &matrix_filename)
 Overwritten method for checking that all pairs in the matrix are correctly loaded. More...
 

Accessors

const T_Aligner_sequence_seqan_wrapper< T_Default_amino_acid_seqan_sequence_converter<>, true, int, seqan::AminoAcid, seqan::CustomAminoAcidMatrix, seqan::Gotoh > & get_aligner (void) const
 
T_Aligner_sequence_seqan_wrapper< T_Default_amino_acid_seqan_sequence_converter<>, true, int, seqan::AminoAcid, seqan::CustomAminoAcidMatrix, seqan::Gotoh > & get_aligner (void)
 
const SequenceType & get_first_sos (void) const
 Returns the first sequence or structure of the alignment. More...
 
const SequenceType & get_second_sos (void) const
 Returns the second sequence or structure of the alignment. More...
 

Algorithm

void align (unsigned verbose, std::ostream &out)
 

Outputs

const Score_typeget_score (void) const
 
Score_typeget_score (void)
 
Alignment_typeget_alignment (void)
 
const Alignment_typeget_alignment (void) const
 

Analysis

void set_substitution_matrix (const Substitution_matrix &matrix)
 The input is a map from pairs of Alignment_unit_name to the associated score. Note that by construction, it is enough to specify non null values on a half of the matrix, diagonal excluded. More...
 
void load_substitution_matrix (const std::string &matrix_filename, unsigned nb_units_names)
 Loads directly the substitution matrix from a file of the format : Alignment_unit_name Alignment_unit_name Value. If duplicated values, the last one is taken. More...
 
void set_blosum_30 (void)
 Set the current matrix as the BLOSUM30 matrix as defined in Seqan 2.0. More...
 
void set_blosum_45 (void)
 Set the current matrix as the BLOSUM45 matrix as defined in Seqan 2.0. More...
 
void set_blosum_62 (void)
 Set the current matrix as the BLOSUM62 matrix as defined in Seqan 2.0. More...
 
void set_blosum_80 (void)
 Set the current matrix as the BLOSUM80 matrix as defined in Seqan 2.0. More...
 
void set_pam_40 (void)
 Set the current matrix as the PAM40 matrix as defined in Seqan 2.0. More...
 
void set_pam_120 (void)
 Set the current matrix as the PAM120 matrix as defined in Seqan 2.0. More...
 
void set_pam_200 (void)
 Set the current matrix as the PAM200 matrix as defined in Seqan 2.0. More...
 
void set_pam_250 (void)
 Set the current matrix as the PAM250 matrix as defined in Seqan 2.0. More...
 
void set_vtml_200 (void)
 Set the current matrix as the VTML200 matrix as defined in Seqan 2.0. More...
 
void set_standard_substitution_matrix (SBL::CSB::Alignment_substitution_matrix_type matrix_type)
 Set the current matrix as one of the standard substitution matrices. More...
 
Substitution_matrixget_substitution_matrix (void)
 
double get_substitution_score (const Alignment_unit_name &p, const Alignment_unit_name &q) const
 
double get_identity_percentage (Sequence_length_type type=ALIGNMENT_SEQUENCE_LENGTH) const
 
double get_similarity_percentage (Sequence_length_type type=ALIGNMENT_SEQUENCE_LENGTH) const
 
void statistics (std::ostream &out, Sequence_length_type type) const
 
void print_alignment_txt (std::ostream &out) const
 Print the alignment in txt format. More...
 
void print_alignment_dot (std::ostream &out) const
 Print the alignment as a dot graph. More...
 

Detailed Description

template<class SequenceType, class SeqanSequenceConverter = T_Default_amino_acid_seqan_sequence_converter<>, bool FreeEndsAlignment = true, class ScoreType = int, class SeqanUnitType = seqan::AminoAcid, class SeqanCustomMatrix = seqan::CustomAminoAcidMatrix, class SeqanAlgorithm = seqan::Gotoh>
class SBL::CSB::T_Alignment_engine_sequences_seqan< SequenceType, SeqanSequenceConverter, FreeEndsAlignment, ScoreType, SeqanUnitType, SeqanCustomMatrix, SeqanAlgorithm >

Engine for making alignments between sequences using Seqan.

It provides a generic interface for wrapping existing algorithms aligning pairs of sequences or structures. It is designed such that it provides all common statistics to both type of algorithms. For a data structure that is more specific to structural alignments, see the class T_Alignement_engine_for_structures

Template Parameters
SequenceTypeRepresentation of a sequence : it requires in particular to define the type Alignment_unit for the base representation of a unit to align (e.g a residue or a nucleotid), the type Alignment_unit_rep for indexing the units in the sequence, the method size() to return the length of the sequence, and the operator [Alignment_unit_rep] for accessing to the corresponding unit of the sequence.
SeqanSequenceConverterConverter between unit names in a sequence defined in the SBL and a unit of a Seqan sequence. It defines the two methods to_seqan_sequence and to_unit_name.
FreeEndsAlignmentSeqan parameter (true by default) (see Simon)
ScoreTypeReturned score type (int by default)
SeqanUnitTypeUnit type in Seqan, that could be amino acid or nucleotid (AminoAcid by default)
SeqanCustomMatrixSeqan matrix defined in this header file for adapting the seqan matrix notation to this engine matrix notation (default is custom for amino acids)
SeqanAlgorithmSeqan algorithm used to solve the alignments (default is Gotoh)

Member Typedef Documentation

◆ Aligned_pair

typedef std::pair<Alignment_unit_rep, Alignment_unit_rep> Aligned_pair
inherited

Representation of two aligned units.

◆ Alignment_type

typedef std::vector<Aligned_pair> Alignment_type
inherited

Representation of an alignment as a sequence of aligned units.

◆ Alignment_unit

typedef SequenceOrStructure::Alignment_unit Alignment_unit
inherited

Type for a unit (e.g residue or nucleotid)

◆ Alignment_unit_name

typedef SequenceOrStructure::Alignment_unit_name Alignment_unit_name
inherited

Representation of the name of a unit.

◆ Alignment_unit_rep

typedef SequenceOrStructure::Alignment_unit_rep Alignment_unit_rep
inherited

Representation of a unit for index purposes.

◆ Name_pair

typedef std::pair<Alignment_unit_name, Alignment_unit_name> Name_pair
inherited

Representation of two unit names.

◆ Score_type

typedef AlignerAlgorithm::Score_type Score_type
inherited

Representation of the score of the algorithm.

◆ Sequence_or_structure

typedef SequenceType Sequence_or_structure
inherited

Type for representing the sequences or the structures to align.

◆ Substitution_matrix

typedef std::map<Name_pair, double , Is_lower_name_pair> Substitution_matrix
inherited

Representation of the storage of values of the subtitution matrix.

Member Enumeration Documentation

◆ Sequence_length_type

enum Sequence_length_type
inherited

Enum for selecting the length of the sequence to consider when computing the similarity / identity percentages.

Member Function Documentation

◆ align()

void align ( void  )
inline

Overwritten method for setting the seqan substitution matrix before the alignment.

◆ get_first_sos()

const SequenceType & get_first_sos ( void  ) const
inlineinherited

Returns the first sequence or structure of the alignment.

◆ get_second_sos()

const SequenceType & get_second_sos ( void  ) const
inlineinherited

Returns the second sequence or structure of the alignment.

◆ load_seqan_substitution_matrix()

void load_seqan_substitution_matrix ( const std::string &  matrix_filename)
inline

Overwritten method for checking that all pairs in the matrix are correctly loaded.

◆ load_substitution_matrix()

void load_substitution_matrix ( const std::string &  matrix_filename,
unsigned  nb_units_names 
)
inlineinherited

Loads directly the substitution matrix from a file of the format : Alignment_unit_name Alignment_unit_name Value. If duplicated values, the last one is taken.

◆ print_alignment_dot()

void print_alignment_dot ( std::ostream &  out) const
inlineinherited

Print the alignment as a dot graph.

◆ print_alignment_txt()

void print_alignment_txt ( std::ostream &  out) const
inlineinherited

Print the alignment in txt format.

◆ set_blosum_30()

void set_blosum_30 ( void  )
inlineinherited

Set the current matrix as the BLOSUM30 matrix as defined in Seqan 2.0.

◆ set_blosum_45()

void set_blosum_45 ( void  )
inlineinherited

Set the current matrix as the BLOSUM45 matrix as defined in Seqan 2.0.

◆ set_blosum_62()

void set_blosum_62 ( void  )
inlineinherited

Set the current matrix as the BLOSUM62 matrix as defined in Seqan 2.0.

◆ set_blosum_80()

void set_blosum_80 ( void  )
inlineinherited

Set the current matrix as the BLOSUM80 matrix as defined in Seqan 2.0.

◆ set_pam_120()

void set_pam_120 ( void  )
inlineinherited

Set the current matrix as the PAM120 matrix as defined in Seqan 2.0.

◆ set_pam_200()

void set_pam_200 ( void  )
inlineinherited

Set the current matrix as the PAM200 matrix as defined in Seqan 2.0.

◆ set_pam_250()

void set_pam_250 ( void  )
inlineinherited

Set the current matrix as the PAM250 matrix as defined in Seqan 2.0.

◆ set_pam_40()

void set_pam_40 ( void  )
inlineinherited

Set the current matrix as the PAM40 matrix as defined in Seqan 2.0.

◆ set_standard_substitution_matrix()

void set_standard_substitution_matrix ( SBL::CSB::Alignment_substitution_matrix_type  matrix_type)
inlineinherited

Set the current matrix as one of the standard substitution matrices.

◆ set_substitution_matrix()

void set_substitution_matrix ( const Substitution_matrix matrix)
inlineinherited

The input is a map from pairs of Alignment_unit_name to the associated score. Note that by construction, it is enough to specify non null values on a half of the matrix, diagonal excluded.

◆ set_vtml_200()

void set_vtml_200 ( void  )
inlineinherited

Set the current matrix as the VTML200 matrix as defined in Seqan 2.0.