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

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

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

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

Public Member Functions

Constructors, destructor and assignment
 pairwise_alignment_algorithm_banded ()=default
 Defaulted.
 
 pairwise_alignment_algorithm_banded (pairwise_alignment_algorithm_banded const &)=default
 Defaulted.
 
 pairwise_alignment_algorithm_banded (pairwise_alignment_algorithm_banded &&)=default
 Defaulted.
 
pairwise_alignment_algorithm_bandedoperator= (pairwise_alignment_algorithm_banded const &)=default
 Defaulted.
 
pairwise_alignment_algorithm_bandedoperator= (pairwise_alignment_algorithm_banded &&)=default
 Defaulted.
 
 ~pairwise_alignment_algorithm_banded ()=default
 Defaulted.
 
 pairwise_alignment_algorithm_banded (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.
 
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.
 

Protected Types

using base_algorithm_t = pairwise_alignment_algorithm< alignment_configuration_t, policies_t... >
 The alignment configuration traits type with auxiliary information extracted from the configuration type.
 
- Protected Types inherited from seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_t... >
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::forward_range alignment_column_t, std::ranges::input_range cell_index_column_t, typename alphabet1_t , std::ranges::input_range sequence2_t>
void compute_band_column (alignment_column_t &&alignment_column, cell_index_column_t &&cell_index_column, alphabet1_t const &alphabet1, sequence2_t &&sequence2)
 Computes a column of the band that does not start in the first row of the alignment matrix.
 
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 banded alignment.
 
- Protected Member Functions inherited from seqan3::detail::pairwise_alignment_algorithm< alignment_configuration_t, policies_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.
 
void compute_matrix (sequence1_t &&sequence1, sequence2_t &&sequence2, alignment_matrix_t &&alignment_matrix, index_matrix_t &&index_matrix)
 Compute the actual alignment.
 
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.
 
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.
 
::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.
 
 pairwise_alignment_algorithm ()=default
 Defaulted.
 
 pairwise_alignment_algorithm (pairwise_alignment_algorithm const &)=default
 Defaulted.
 
 pairwise_alignment_algorithm (pairwise_alignment_algorithm &&)=default
 Defaulted.
 
 pairwise_alignment_algorithm (alignment_configuration_t const &config)
 Constructs and initialises the algorithm using the alignment configuration.
 
pairwise_alignment_algorithmoperator= (pairwise_alignment_algorithm const &)=default
 Defaulted.
 
pairwise_alignment_algorithmoperator= (pairwise_alignment_algorithm &&)=default
 Defaulted.
 
 ~pairwise_alignment_algorithm ()=default
 Defaulted.
 
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.
 

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_banded< alignment_configuration_t, policies_t >

The alignment algorithm type to compute the banded 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_banded()

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

Constructs and initialises the algorithm using the alignment configuration.

Parameters
configThe configuration passed into the algorithm.

Initialises the algorithm given the user settings from the alignment configuration object.

Exceptions
seqan3::invalid_alignment_configuration.

Member Function Documentation

◆ compute_band_column()

template<typename alignment_configuration_t , typename... policies_t>
template<std::ranges::forward_range alignment_column_t, std::ranges::input_range cell_index_column_t, typename alphabet1_t , std::ranges::input_range sequence2_t>
void seqan3::detail::pairwise_alignment_algorithm_banded< alignment_configuration_t, policies_t >::compute_band_column ( alignment_column_t &&  alignment_column,
cell_index_column_t &&  cell_index_column,
alphabet1_t const &  alphabet1,
sequence2_t &&  sequence2 
)
inlineprotected

Computes a column of the band that does not start in the first row of the alignment matrix.

Template Parameters
alignment_column_tThe type of the alignment column; must model std::ranges::forward_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. Note that the length of sequence2 determines the size of the column.

Implementation with a score matrix using linear memory

When the score matrix uses only linear memory, i.e. only one column is stored, the algorithm reuses the cells of the same column to compute the current one. This means, before the current cell is updated its value is cached and used as the previous diagonal value when computing the next cell below. In the banded case, however, the position of the referenced cell is shifted by one and needs a different handling.

0 1 2 3 4 5 6
_____________
0 |0|0|0|\ |
1 |1|1|1|0|\ |
2 |\|2|2|1|0|\ |
3 | \|3|2|1|0|\|
4 | \|3|2|1|0|
5 | \|3|2|1|
–––––––––––––

The picture above depicts the banded matrix with the indices of the column. As long as the band touches the first row (column 0 - 2), the indices of the actual stored column refer to the same position as the previous column. Hence, to compute these columns the regular algorithm seqan3::detail::pairwise_alignment_algorithm::compute_column can be used. Starting at column 3 the first cell of the band moves downwards, causing a shift in the position of the previous cell. In this state, the value of the current cell represents the previous diagonal value before it is updated. To read the previous horizontal value the next cell below has to be dereferenced. Accordingly, two iterators are used to point to the respective cells in the matrix. The first one points to the current cell (the one that is written to) and the second points to the next cell (the one where the horizontal and vertical scores are read from). After computing the last cell of the column the value of the current iterator can be used to track the score of the cell.

◆ 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_banded< 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 banded 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.

In the banded alignment the iteration of the inner columns is split into two phases. The first phase reuses the unbanded column computation and assumes that the banded score matrix always starts at the beginning of the matrix. In the second phase the special interface seqan3::detail::pairwise_alignment_algorithm_banded::compute_band_column is used to compute the banded column. The current implementation (this might change when we need to work with full-matrices like Waterman-Eggert does) assumes that the first cell of the current score matrix column is always the first cell to compute in every column. The respective seqan3::detail::score_matrix_single_column is resized to the band size plus one additional field to cover the end of the band where the end of the column is not reached yet. This cell will never be written to but only read from, i.e. it represents minus infinity. This allows the algorithm to reuse the standard unbanded cell computation. The following figure depicts the referenced cell of the underlying score matrix (assuming seqan3::detail::score_matrix_single_column):

A G G T C A
0 1 2 3 4 5 6
|–|—|—|—|—|—|—|
0|0|0|0|0|0| | |
A 1|1|1|1|1|1|0| |
C 2|x|2|2|2|2|1|0|
G 3| |x|3|3|3|2|1|
T 4| | |x|4|4|3|2|

The coordinate matrix represents the global matrix index and not the local band coordinate. Data structures that require the coordinate might need to map the global matrix coordinate to their local coordinate:

A G G T C A
0 1 2 3 4 5 6
|–––––|–––––|–––––|–––––|–––––|–––––|–––––|
0|(0,0)|(0,1)|(0,2)|(0,3)|(0,4)| | |
A 1|(1,0)|(1,1)|(1,2)|(1,3)|(1,4)|(1,5)| |
C 2| |(2,1)|(2,2)|(2,3)|(2,4)|(2,5)|(2,6)|
G 3| | |(3,2)|(3,3)|(3,4)|(3,5)|(3,6)|
T 4| | | |(4,3)|(4,4)|(4,5)|(4,6)|

◆ 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_banded< 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