NCBI C++ ToolKit
|
Smith-Waterman gapped alignment, for use with infrastructure of BLAST. More...
Go to the source code of this file.
Go to the SVN repository for this file.
Classes | |
struct | BlastGapSW |
Auxiliary structures for Smith-Waterman alignment with traceback. More... | |
Macros | |
#define | SWAP_SEQS(A, B) {const Uint1 *tmp = (A); (A) = (B); (B) = tmp; } |
swap (pointers to) a pair of sequences More... | |
#define | SWAP_INT(A, B) {Int4 tmp = (A); (A) = (B); (B) = tmp; } |
swap two integers More... | |
Typedefs | |
typedef struct BlastGapSW | BlastGapSW |
Auxiliary structures for Smith-Waterman alignment with traceback. More... | |
Enumerations | |
enum | { EDIT_SUB = eGapAlignSub , EDIT_GAP_IN_A = eGapAlignDel , EDIT_GAP_IN_B = eGapAlignIns , EDIT_OP_MASK = 0x07 , EDIT_START_GAP_A = 0x10 , EDIT_START_GAP_B = 0x20 } |
Values for the editing script operations in traceback. More... | |
Functions | |
static Int4 | s_SmithWatermanScoreOnly (const Uint1 *A, Int4 a_size, const Uint1 *B, Int4 b_size, Int4 gap_open, Int4 gap_extend, BlastGapAlignStruct *gap_align) |
Compute the score of the best local alignment between two protein sequences. More... | |
static Int4 | s_NuclSmithWaterman (const Uint1 *B, Int4 b_size, const Uint1 *A, Int4 a_size, Int4 gap_open, Int4 gap_extend, BlastGapAlignStruct *gap_align) |
Compute the score of the best local alignment between two nucleotide sequences. More... | |
static void | s_GetTraceback (EBlastProgramType program_number, Uint1 *trace, const Uint1 *A, const Uint1 *B, Int4 b_size, Int4 gap_open, Int4 gap_extend, BlastGapAlignStruct *gap_align, Int4 a_end, Int4 b_end, Int4 best_score, BlastHSPList *hsp_list, Boolean swapped, BlastHSP *template_hsp, const BlastScoringOptions *score_options, const BlastHitSavingOptions *hit_options, Int4 start_shift) |
Generate the traceback for a single local alignment between two (unpacked) sequences, then create an HSP for the alignment and add to a list of such HSPs. More... | |
void | SmithWatermanScoreWithTraceback (EBlastProgramType program_number, const Uint1 *A, Int4 a_size, const Uint1 *B, Int4 b_size, BlastHSP *template_hsp, BlastHSPList *hsp_list, const BlastScoringParameters *score_params, const BlastHitSavingParameters *hit_params, BlastGapAlignStruct *gap_align, Int4 start_shift, Int4 cutoff) |
Find all local alignments between two (unpacked) sequences, using the Smith-Waterman algorithm, then save the list of alignments found. More... | |
Int2 | BLAST_SmithWatermanGetGappedScore (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, BLAST_SequenceBlk *subject, BlastGapAlignStruct *gap_align, const BlastScoringParameters *score_params, const BlastExtensionParameters *ext_params, const BlastHitSavingParameters *hit_params, const BlastInitialWordParameters *word_params, BlastInitHitList *init_hitlist, BlastHSPList **hsp_list_ptr, BlastGappedStats *gapped_stats, Boolean *fence_hit) |
Performs score-only Smith-Waterman gapped alignment of the subject sequence with all contexts in the query. More... | |
Smith-Waterman gapped alignment, for use with infrastructure of BLAST.
Definition in file blast_sw.c.
swap two integers
Definition at line 42 of file blast_sw.c.
swap (pointers to) a pair of sequences
Definition at line 39 of file blast_sw.c.
typedef struct BlastGapSW BlastGapSW |
Auxiliary structures for Smith-Waterman alignment with traceback.
anonymous enum |
Values for the editing script operations in traceback.
Enumerator | |
---|---|
EDIT_SUB | Substitution. |
EDIT_GAP_IN_A | Deletion. |
EDIT_GAP_IN_B | Insertion. |
EDIT_OP_MASK | Mask for edit script operations. |
EDIT_START_GAP_A | open a gap in A |
EDIT_START_GAP_B | open a gap in B |
Definition at line 243 of file blast_sw.c.
|
static |
Generate the traceback for a single local alignment between two (unpacked) sequences, then create an HSP for the alignment and add to a list of such HSPs.
program_number | Blast program requesting traceback [in] |
trace | 2-D array of edit actions, size (a_size+1) x (b_size+1) [in] |
A | The first sequence [in] |
B | The second sequence [in] |
b_size | Length of the second sequence [in] |
gap_open | Gap open penalty [in] |
gap_extend | Gap extension penalty [in] |
gap_align | Auxiliary data for gapped alignment (used for score matrix info) [in] |
a_end | The alignment end offset on A (plus one) [in] |
b_end | The alignment end offset on B (plus one) [in] |
best_score | Score of the alignment [in] |
hsp_list | Collection of alignments found so far [in][out] |
swapped | TRUE if A and B were swapped before the alignment was found [in] |
template_hsp | Placeholder alignment, used only to determine contexts and frames [in] |
score_options | Structure containing gap penalties [in] |
hit_options | Structure used for percent identity calculation [in] |
start_shift | Bias to be applied to subject offsets [in] |
Definition at line 278 of file blast_sw.c.
References A, B, Blast_HSPAdjustSubjectOffset(), Blast_HSPFree(), Blast_HSPInit(), Blast_HSPListSaveHSP(), Blast_HSPTestIdentityAndLength(), BlastHSP::context, SBlastScoreMatrix::data, EDIT_GAP_IN_A, EDIT_GAP_IN_B, EDIT_OP_MASK, GapPrelimEditBlock::edit_ops, EDIT_START_GAP_A, EDIT_START_GAP_B, EDIT_SUB, eGapAlignDel, eGapAlignIns, BlastSeg::frame, BlastGapAlignStruct::fwd_prelim_tback, GapEditScriptNew(), GapPrelimEditBlockAdd(), GapPrelimEditBlockReset(), i, BlastScoreBlk::matrix, GapEditScript::num, GapPrelimEditScript::num, GapPrelimEditBlock::num_ops, GapEditScript::op_type, GapPrelimEditScript::op_type, BlastGapAlignStruct::positionBased, BlastScoreBlk::psi_matrix, SPsiBlastScoreMatrix::pssm, BlastHSP::query, BlastGapAlignStruct::sbp, BlastHSP::subject, SWAP_INT, SWAP_SEQS, and trace.
Referenced by SmithWatermanScoreWithTraceback().
|
static |
Compute the score of the best local alignment between two nucleotide sequences.
One of the sequences must be in packed format. For nucleotide Smith-Waterman, the vast majority of the runtime is tied up in this routine.
B | The first sequence (must be in ncbi2na format) [in] |
b_size | Length of the first sequence [in] |
A | The second sequence [in] |
a_size | Length of the second sequence [in] |
gap_open | Gap open penalty [in] |
gap_extend | Gap extension penalty [in] |
gap_align | Auxiliary data for gapped alignment (used for score matrix info) [in] |
Definition at line 164 of file blast_sw.c.
References A, B, BlastGapDP::best, BlastGapDP::best_gap, SBlastScoreMatrix::data, BlastGapAlignStruct::dp_mem, BlastGapAlignStruct::dp_mem_alloc, i, malloc(), BlastScoreBlk::matrix, MAX, NCBI2NA_UNPACK_BASE, BlastGapAlignStruct::sbp, and sfree.
Referenced by BLAST_SmithWatermanGetGappedScore().
|
static |
Compute the score of the best local alignment between two protein sequences.
When using Smith-Waterman, the vast majority of the runtime is tied up in this routine.
A | The first sequence [in] |
a_size | Length of the first sequence [in] |
B | The second sequence [in] |
b_size | Length of the second sequence [in] |
gap_open | Gap open penalty [in] |
gap_extend | Gap extension penalty [in] |
gap_align | Auxiliary data for gapped alignment (used for score matrix info) [in] |
Definition at line 57 of file blast_sw.c.
References A, B, BlastGapDP::best, BlastGapDP::best_gap, SBlastScoreMatrix::data, BlastGapAlignStruct::dp_mem, BlastGapAlignStruct::dp_mem_alloc, i, malloc(), BlastScoreBlk::matrix, MAX, BlastGapAlignStruct::positionBased, BlastScoreBlk::psi_matrix, SPsiBlastScoreMatrix::pssm, BlastGapAlignStruct::sbp, sfree, SWAP_INT, and SWAP_SEQS.
Referenced by BLAST_SmithWatermanGetGappedScore().