SeqAn3 3.4.0-rc.1
The Modern C++ library for sequence analysis.
Loading...
Searching...
No Matches
seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t > Class Template Reference

The alignment algorithm type to compute standard pairwise alignment using dynamic programming. More...

#include <seqan3/alignment/pairwise/detail/pairwise_alignment_algorithm.hpp>

+ Inheritance diagram for seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >:

Public Member Functions

template<indexed_sequence_pair_range indexed_sequence_pairs_t, typename callback_t >
requires traits_type
::is_vectorised &&std::invocable< callback_t, alignment_result_type > auto operator() (indexed_sequence_pairs_t &&indexed_sequence_pairs, callback_t &&callback)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
Constructors, destructor and assignment
 pairwise_alignment_algorithm ()=default
 Defaulted.
 
 pairwise_alignment_algorithm (pairwise_alignment_algorithm const &)=default
 Defaulted.
 
 pairwise_alignment_algorithm (pairwise_alignment_algorithm &&)=default
 Defaulted.
 
pairwise_alignment_algorithmoperator= (pairwise_alignment_algorithm const &)=default
 Defaulted.
 
pairwise_alignment_algorithmoperator= (pairwise_alignment_algorithm &&)=default
 Defaulted.
 
 ~pairwise_alignment_algorithm ()=default
 Defaulted.
 
 pairwise_alignment_algorithm (alignment_configuration_t const &config)
 Constructs and initialises the algorithm using the alignment configuration.
 
Invocation
template<indexed_sequence_pair_range indexed_sequence_pairs_t, typename callback_t >
requires std::invocable<callback_t, alignment_result_type>
void operator() (indexed_sequence_pairs_t &&indexed_sequence_pairs, callback_t &&callback)
 Computes the pairwise sequence alignment for the given range over indexed sequence pairs.
 

Protected Types

using alignment_result_type = typename traits_type::alignment_result_type
 The configured alignment result type.
 
using score_type = typename traits_type::score_type
 The configured score type.
 
using traits_type = alignment_configuration_traits< alignment_configuration_t >
 The alignment configuration traits type with auxiliary information extracted from the configuration type.
 

Protected Member Functions

template<std::ranges::input_range alignment_column_t, std::ranges::input_range cell_index_column_t, typename alphabet1_t , std::ranges::input_range sequence2_t>
requires semialphabet<alphabet1_t> || simd_concept<alphabet1_t>
void compute_column (alignment_column_t &&alignment_column, cell_index_column_t &&cell_index_column, alphabet1_t const &alphabet1, sequence2_t &&sequence2)
 Initialise any column of the alignment matrix except the first one.
 
template<std::ranges::forward_range sequence1_t, std::ranges::forward_range sequence2_t, std::ranges::input_range alignment_matrix_t, std::ranges::input_range index_matrix_t>
requires std::ranges::forward_range<std::ranges::range_reference_t<alignment_matrix_t>> && std::ranges::forward_range<std::ranges::range_reference_t<index_matrix_t>>
void compute_matrix (sequence1_t &&sequence1, sequence2_t &&sequence2, alignment_matrix_t &&alignment_matrix, index_matrix_t &&index_matrix)
 Compute the actual alignment.
 
template<typename simd_sequence_t , std::ranges::forward_range sequence_collection_t, arithmetic padding_symbol_t>
requires std::ranges::output_range<simd_sequence_t, score_type>
void convert_batch_of_sequences_to_simd_vector (simd_sequence_t &simd_sequence, sequence_collection_t &sequences, padding_symbol_t const &padding_symbol)
 Converts a batch of sequences to a sequence of simd vectors.
 
template<std::ranges::input_range alignment_column_t, std::ranges::input_range cell_index_column_t, std::ranges::input_range sequence2_t>
void initialise_column (alignment_column_t &&alignment_column, cell_index_column_t &&cell_index_column, sequence2_t &&sequence2)
 Initialise the first column of the alignment matrix.
 

Detailed Description

template<typename alignment_configuration_t, typename... policies_t>
requires is_type_specialisation_of_v<alignment_configuration_t, configuration>
class seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >

The alignment algorithm type to compute standard pairwise alignment using dynamic programming.

Template Parameters
alignment_configuration_tThe configuration type; must be of type seqan3::configuration.
policies_tVariadic template argument for the different policies of this alignment algorithm.

Configuration

The first template argument is the type of the alignment configuration. The alignment configuration was used to configure the alignment algorithm type within the seqan3::detail::alignment_configurator. The algorithm computes a column based dynamic programming matrix given two sequences. After the computation a user defined callback function is invoked with the computed seqan3::alignment_result.

Constructor & Destructor Documentation

◆ pairwise_alignment_algorithm()

template<typename alignment_configuration_t , typename... policies_t>
seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::pairwise_alignment_algorithm ( alignment_configuration_t const &  config)
inline

Constructs and initialises the algorithm using the alignment configuration.

Parameters
configThe configuration passed into the algorithm.

Initialises the base policies of the alignment algorithm.

Member Function Documentation

◆ compute_column()

template<typename alignment_configuration_t , typename... policies_t>
template<std::ranges::input_range alignment_column_t, std::ranges::input_range cell_index_column_t, typename alphabet1_t , std::ranges::input_range sequence2_t>
requires semialphabet<alphabet1_t> || simd_concept<alphabet1_t>
void seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::compute_column ( alignment_column_t &&  alignment_column,
cell_index_column_t &&  cell_index_column,
alphabet1_t const &  alphabet1,
sequence2_t &&  sequence2 
)
inlineprotected

Initialise any column of the alignment matrix except the first one.

Template Parameters
alignment_column_tThe type of the alignment column; must model std::ranges::input_range.
cell_index_column_tThe type of the indexed column; must model std::ranges::input_range.
alphabet1_tThe type of the current symbol of sequence1.
sequence2_tThe type of the second sequence; must model std::ranges::input_range.
Parameters
[in]alignment_columnThe current alignment matrix column to compute.
[in]cell_index_columnThe current index matrix column to get the respective cell indices.
[in]alphabet1The current symbol of sequence1.
[in]sequence2The second sequence to align against alphabet1.

Computes the alignment for the given alignment matrix column. The function splits the computation of the column into three phases: the initialisation phase, the iteration phase, and the final phase. In the initialisation phase the first cell of the column is computed and in the iteration phase all remaining cells are computed. In the final phase the last cell is possibly evaluated for a new alignment optimum.

◆ compute_matrix()

template<typename alignment_configuration_t , typename... policies_t>
template<std::ranges::forward_range sequence1_t, std::ranges::forward_range sequence2_t, std::ranges::input_range alignment_matrix_t, std::ranges::input_range index_matrix_t>
requires std::ranges::forward_range<std::ranges::range_reference_t<alignment_matrix_t>> && std::ranges::forward_range<std::ranges::range_reference_t<index_matrix_t>>
void seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::compute_matrix ( sequence1_t &&  sequence1,
sequence2_t &&  sequence2,
alignment_matrix_t &&  alignment_matrix,
index_matrix_t &&  index_matrix 
)
inlineprotected

Compute the actual alignment.

Template Parameters
sequence1_tThe type of the first sequence; must model std::ranges::forward_range.
sequence2_tThe type of the second sequence; must model std::ranges::forward_range.
alignment_matrix_tThe type of the alignment matrix; must model std::ranges::input_range and its std::ranges::range_reference_t type must model std::ranges::forward_range.
index_matrix_tThe type of the index matrix; must model std::ranges::input_range and its std::ranges::range_reference_t type must model std::ranges::forward_range.
Parameters
[in]sequence1The first sequence to compute the alignment for.
[in]sequence2The second sequence to compute the alignment for.
[in]alignment_matrixThe alignment matrix to compute.
[in]index_matrixThe index matrix corresponding to the alignment matrix.

◆ convert_batch_of_sequences_to_simd_vector()

template<typename alignment_configuration_t , typename... policies_t>
template<typename simd_sequence_t , std::ranges::forward_range sequence_collection_t, arithmetic padding_symbol_t>
requires std::ranges::output_range<simd_sequence_t, score_type>
void seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::convert_batch_of_sequences_to_simd_vector ( simd_sequence_t &  simd_sequence,
sequence_collection_t &  sequences,
padding_symbol_t const &  padding_symbol 
)
inlineprotected

Converts a batch of sequences to a sequence of simd vectors.

Template Parameters
simd_sequence_tThe type of the simd sequence; must model std::ranges::output_range for the score_type.
sequence_collection_tThe type of the collection containing the sequences; must model std::ranges::forward_range.
padding_symbol_tThe type of the padding symbol.
Parameters
[out]simd_sequenceThe transformed simd sequence.
[in]sequencesThe batch of sequences to transform.
[in]padding_symbolThe symbol that should be appended during the transformation/

Expects that the size of the collection is less or equal than the number of alignments that can be computed within one simd vector (simd_traits<score_type>::length). Applies an Array-of-Structures (AoS) to Structure-of-Arrays (SoA) transformation by storing one column of the collection as a simd vector. The resulting simd sequence has the size of the longest sequence in the collection. For all sequences with a smaller size the padding symbol will be appended during the simd transformation to fill up the remaining size difference.

◆ initialise_column()

template<typename alignment_configuration_t , typename... policies_t>
template<std::ranges::input_range alignment_column_t, std::ranges::input_range cell_index_column_t, std::ranges::input_range sequence2_t>
void seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::initialise_column ( alignment_column_t &&  alignment_column,
cell_index_column_t &&  cell_index_column,
sequence2_t &&  sequence2 
)
inlineprotected

Initialise the first column of the alignment matrix.

Template Parameters
alignment_column_tThe type of the alignment column; must model std::ranges::input_range.
cell_index_column_tThe type of the indexed column; must model std::ranges::input_range.
sequence2_tThe type of the second sequence; must model std::ranges::input_range.
Parameters
[in]alignment_columnThe current alignment matrix column to compute.
[in]cell_index_columnThe current index matrix column to get the respective cell indices.
[in]sequence2The second sequence used to determine the size of the column.

The first column of the alignment matrix does not require any character comparisons of the sequences that shall be aligned. The second sequence is thus only needed to determine the size of the column. The computation of the column is split into three phases: the initialisation phase, the iteration phase, and the final phase. In the initialisation phase the first cell of the column is computed and in the iteration phase all remaining cells are computed. In the final phase the last cell is possibly evaluated for a new alignment optimum.

◆ operator()()

template<typename alignment_configuration_t , typename... policies_t>
template<indexed_sequence_pair_range indexed_sequence_pairs_t, typename callback_t >
requires std::invocable<callback_t, alignment_result_type>
void seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t >::operator() ( indexed_sequence_pairs_t &&  indexed_sequence_pairs,
callback_t &&  callback 
)
inline

Computes the pairwise sequence alignment for the given range over indexed sequence pairs.

Template Parameters
indexed_sequence_pairs_tThe type of indexed_sequence_pairs; must model seqan3::detail::indexed_sequence_pair_range.
callback_tThe type of the callback function that is called with the alignment result; must model std::invocable with seqan3::alignment_result as argument.
Parameters
[in]indexed_sequence_pairsA range over indexed sequence pairs to be aligned.
[in]callbackThe callback function to be invoked with each computed alignment result.
Exceptions
std::bad_allocduring allocation of the alignment matrices or seqan3::invalid_alignment_configuration if an invalid configuration for the given sequences is detected.

Uses the standard dynamic programming algorithm to compute the pairwise sequence alignment for each sequence pair. The space and runtime complexities depend on the selected configurations (see below). For every computed alignment the given callback is invoked with the respective alignment result.

Exception

Strong exception guarantee. Might throw std::bad_alloc or seqan3::invalid_alignment_configuration.

Thread-safety

Calls to this functions in a concurrent environment are not thread safe. Instead use a copy of the alignment algorithm type.

Complexity

The following table lists the runtime and space complexities for the banded and unbanded algorithm dependent on the given seqan3::align_cfg::output_* per sequence pair. Let n be the length of the first sequence, m be the length of the second sequence and k be the size of the band.

unbanded banded
runtime \( O(n*m) \) \( O(n*k) \)
space (score only) \( O(m) \) \( O(k) \)
space (end positions) \( O(m) \) \( O(k) \)
space (begin positions) \( O(n*m) \) \( O(n*k) \)
space (alignment) \( O(n*m) \) \( O(n*k) \)

The documentation for this class was generated from the following file:
Hide me