NCBI C++ ToolKit
Classes | Macros | Typedefs | Functions | Variables
blast_engine.c File Reference

Function calls to actually perform a BLAST search (high level). More...

#include <algo/blast/core/blast_engine.h>
#include <algo/blast/core/blast_util.h>
#include <algo/blast/core/blast_aalookup.h>
#include <algo/blast/core/blast_aascan.h>
#include <algo/blast/core/blast_nalookup.h>
#include <algo/blast/core/blast_nascan.h>
#include <algo/blast/core/blast_sw.h>
#include <algo/blast/core/aa_ungapped.h>
#include <algo/blast/core/na_ungapped.h>
#include <algo/blast/core/phi_lookup.h>
#include <algo/blast/core/phi_gapalign.h>
#include <algo/blast/core/phi_extend.h>
#include <algo/blast/core/link_hsps.h>
#include <algo/blast/core/blast_gapalign.h>
#include <algo/blast/core/blast_parameters.h>
#include <algo/blast/core/blast_setup.h>
#include <algo/blast/core/blast_traceback.h>
#include <algo/blast/core/mb_indexed_lookup.h>
#include <algo/blast/core/gencode_singleton.h>
#include "blast_gapalign_priv.h"
#include "jumper.h"
+ Include dependency graph for blast_engine.c:

Go to the source code of this file.

Go to the SVN repository for this file.

Classes

struct  BlastCoreAuxStruct
 Structure to be passed to s_BlastSearchEngineCore, containing pointers to various preallocated structures and arrays. More...
 
struct  SubjectSplitStruct
 Structure used for subject sequence split. More...
 

Macros

#define CONV_NUCL2PROT_COORDINATES(length)   (length) / CODON_LENGTH
 Converts nucleotide coordinates to protein. More...
 

Typedefs

typedef struct BlastCoreAuxStruct BlastCoreAuxStruct
 Structure to be passed to s_BlastSearchEngineCore, containing pointers to various preallocated structures and arrays. More...
 
typedef struct SubjectSplitStruct SubjectSplitStruct
 Structure used for subject sequence split. More...
 

Functions

static BlastCoreAuxStructs_BlastCoreAuxStructFree (BlastCoreAuxStruct *aux_struct)
 Deallocates all memory in BlastCoreAuxStruct. More...
 
static void s_BackupSubject (BLAST_SequenceBlk *subject, SubjectSplitStruct *backup)
 
static void s_AllocateSeqRange (BLAST_SequenceBlk *subject, SubjectSplitStruct *backup, Int4 num_seq_ranges)
 
static void s_RestoreSubject (BLAST_SequenceBlk *subject, SubjectSplitStruct *backup)
 
static Int2 s_GetNextSubjectChunk (BLAST_SequenceBlk *subject, SubjectSplitStruct *backup, Boolean is_nucleotide, int chunk_overlap)
 
static void s_TranslateHSPsToDNAPCoord (EBlastProgramType program, BlastInitHitList *init_hitlist, const BlastQueryInfo *query_info, Int2 subject_frame, Int4 subject_length, Int4 offset)
 Adjust HSP coordinates for out-of-frame gapped extension. More...
 
static void s_RPSOffsetArrayToContextOffsets (BlastQueryInfo *info, Int4 *new_offsets)
 Set up context offsets for the auxiliary BlastQueryInfo structure that is created for the concatenated database in RPS BLAST search. More...
 
static Int2 s_BlastSearchEngineOneContext (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, BLAST_SequenceBlk *subject, Int4 orig_length, LookupTableWrap *lookup, BlastGapAlignStruct *gap_align, const BlastScoringParameters *score_params, const BlastInitialWordParameters *word_params, const BlastExtensionParameters *ext_params, const BlastHitSavingParameters *hit_params, BlastDiagnostics *diagnostics, BlastCoreAuxStruct *aux_struct, BlastHSPList **hsp_list_out_ptr, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 Searches only one context of a database sequence, but does all chunks if it is split. More...
 
static void s_BlastSearchEngineCoreCleanUp (EBlastProgramType program_number, BlastQueryInfo *query_info, const BlastQueryInfo *query_info_in, Uint1 *translation_buffer, Uint4 *frame_offsets_a)
 Clean up function for s_BlastSearchEngineCore. More...
 
static Int2 s_Blast_HSPListReapByPrelimEvalue (BlastHSPList *hsp_list, const BlastHitSavingParameters *hit_params)
 Discard the HSPs above the prelim e-value threshold from the HSP list. More...
 
static Int2 s_BlastSearchEngineCore (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info_in, BLAST_SequenceBlk *subject, LookupTableWrap *lookup, BlastGapAlignStruct *gap_align, const BlastScoringParameters *score_params, const BlastInitialWordParameters *word_params, const BlastExtensionParameters *ext_params, const BlastHitSavingParameters *hit_params, const BlastDatabaseOptions *db_options, BlastDiagnostics *diagnostics, BlastCoreAuxStruct *aux_struct, BlastHSPList **hsp_list_out_ptr, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 The core of the BLAST search: comparison between the (concatenated) query against one subject sequence. More...
 
static Int2 s_FillReturnCutoffsInfo (BlastRawCutoffs *return_cutoffs, const BlastScoringParameters *score_params, const BlastInitialWordParameters *word_params, const BlastExtensionParameters *ext_params, const BlastHitSavingParameters *hit_params)
 Fills the output information about the cutoffs uses in a BLAST search. More...
 
static Int2 s_BlastSetUpAuxStructures (const BlastSeqSrc *seq_src, LookupTableWrap *lookup_wrap, const BlastInitialWordParameters *word_params, const BlastExtensionOptions *ext_options, const BlastHitSavingOptions *hit_options, BLAST_SequenceBlk *query, const BlastQueryInfo *query_info, BlastCoreAuxStruct **aux_struct_ptr)
 Setup of the auxiliary BLAST structures; also calculates internally used parameters from options. More...
 
static Int2 s_RPSPreliminarySearchEngine (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, const BlastScoringParameters *score_params, LookupTableWrap *lookup_wrap, BlastCoreAuxStruct *aux_struct, const BlastInitialWordParameters *word_params, const BlastExtensionParameters *ext_params, BlastGapAlignStruct *gap_align, const BlastHitSavingParameters *hit_params, BlastHSPStream *hsp_stream, BlastDiagnostics *diagnostics, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 Performs the preliminary stage of an RPS BLAST search, after all set up has already been done. More...
 
static void s_AdjustSubjectForTranslatedSraSearch (BlastHSPList *hsp_list, Uint1 offset, Int4 length)
 
static void s_AdjustSubjectForSraSearch (BlastHSPList *hsp_list, Uint1 offset)
 
static Int4 s_GetMinimumSubjSeqLen (LookupTableWrap *lookup_wrap)
 
Int4 BLAST_PreliminarySearchEngine (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, BlastGapAlignStruct *gap_align, BlastScoringParameters *score_params, LookupTableWrap *lookup_wrap, const BlastInitialWordOptions *word_options, BlastExtensionParameters *ext_params, BlastHitSavingParameters *hit_params, BlastEffectiveLengthsParameters *eff_len_params, const PSIBlastOptions *psi_options, const BlastDatabaseOptions *db_options, BlastHSPStream *hsp_stream, BlastDiagnostics *diagnostics, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 Perform the preliminary stage of the BLAST search. More...
 
Int2 Blast_RunPreliminarySearch (EBlastProgramType program, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, const BlastScoringOptions *score_options, BlastScoreBlk *sbp, LookupTableWrap *lookup_wrap, const BlastInitialWordOptions *word_options, const BlastExtensionOptions *ext_options, const BlastHitSavingOptions *hit_options, const BlastEffectiveLengthsOptions *eff_len_options, const PSIBlastOptions *psi_options, const BlastDatabaseOptions *db_options, BlastHSPStream *hsp_stream, BlastDiagnostics *diagnostics)
 The high level function performing the BLAST search against a BLAST database after all the setup has been done. More...
 
Int2 Blast_RunPreliminarySearchWithInterrupt (EBlastProgramType program, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, const BlastScoringOptions *score_options, BlastScoreBlk *sbp, LookupTableWrap *lookup_wrap, const BlastInitialWordOptions *word_options, const BlastExtensionOptions *ext_options, const BlastHitSavingOptions *hit_options, const BlastEffectiveLengthsOptions *eff_len_options, const PSIBlastOptions *psi_options, const BlastDatabaseOptions *db_options, BlastHSPStream *hsp_stream, BlastDiagnostics *diagnostics, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 Same as above, with support for user interrupt function. More...
 
static void s_BlastRunFullSearchCleanUp (BlastGapAlignStruct *gap_align, BlastScoringParameters *score_params, BlastExtensionParameters *ext_params, BlastHitSavingParameters *hit_params, BlastEffectiveLengthsParameters *eff_len_params)
 Function to deallocate data structures allocated in Blast_RunFullSearch. More...
 
Int4 Blast_RunFullSearch (EBlastProgramType program_number, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, BlastScoreBlk *sbp, const BlastScoringOptions *score_options, LookupTableWrap *lookup_wrap, const BlastInitialWordOptions *word_options, const BlastExtensionOptions *ext_options, const BlastHitSavingOptions *hit_options, const BlastEffectiveLengthsOptions *eff_len_options, const PSIBlastOptions *psi_options, const BlastDatabaseOptions *db_options, BlastHSPStream *hsp_stream, const BlastRPSInfo *rps_info, BlastDiagnostics *diagnostics, BlastHSPResults **results, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
 The high level function performing the BLAST search against a BLAST database after all the setup has been done. More...
 

Variables

const int kBlastMajorVersion = 2
 Major version. More...
 
const int kBlastMinorVersion = 16
 Minor version. More...
 
const int kBlastPatchVersion = 1
 Patch version. More...
 
const Int2 SUBJECT_SPLIT_DONE = 0
 return value indicating hitting the end More...
 
const Int2 SUBJECT_SPLIT_OK = 1
 return value indicating OK More...
 
const Int2 SUBJECT_SPLIT_NO_RANGE = 2
 return value indicating all masked out More...
 

Detailed Description

Function calls to actually perform a BLAST search (high level).

The hierarchy of function calls, starting from the top level in the BLAST core, is described below.

Preliminary stage of the BLAST search:

   Blast_RunPreliminarySearch
       BLAST_GapAlignSetUp
       BLAST_PreliminarySearchEngine
           if (RPS BLAST) {
               s_RPSPreliminarySearchEngine
                   s_BlastSearchEngineCore
           } else {
               for (all sequences in the database)
                   s_BlastSearchEngineCore
           }

Full BLAST search, including preliminary search and traceback:

   Blast_RunFullSearch
       BLAST_GapAlignSetUp
       BLAST_PreliminarySearchEngine
       BLAST_ComputeTraceback

Definition in file blast_engine.c.

Macro Definition Documentation

◆ CONV_NUCL2PROT_COORDINATES

#define CONV_NUCL2PROT_COORDINATES (   length)    (length) / CODON_LENGTH

Converts nucleotide coordinates to protein.

Definition at line 79 of file blast_engine.c.

Typedef Documentation

◆ BlastCoreAuxStruct

Structure to be passed to s_BlastSearchEngineCore, containing pointers to various preallocated structures and arrays.

◆ SubjectSplitStruct

Structure used for subject sequence split.

Function Documentation

◆ BLAST_PreliminarySearchEngine()

Int4 BLAST_PreliminarySearchEngine ( EBlastProgramType  program_number,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
const BlastSeqSrc seq_src,
BlastGapAlignStruct gap_align,
BlastScoringParameters score_params,
LookupTableWrap lookup_wrap,
const BlastInitialWordOptions word_options,
BlastExtensionParameters ext_params,
BlastHitSavingParameters hit_params,
BlastEffectiveLengthsParameters eff_len_params,
const PSIBlastOptions psi_options,
const BlastDatabaseOptions db_options,
BlastHSPStream hsp_stream,
BlastDiagnostics diagnostics,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)

Perform the preliminary stage of the BLAST search.

Parameters
program_numberType of BLAST program [in]
queryThe query sequence [in]
query_infoAdditional query information [in]
seq_srcStructure containing BLAST database [in]
gap_alignStructure containing scoring block and memory allocated for gapped alignment. [in]
score_paramsHit scoring parameters [in]
lookup_wrapThe lookup table, constructed earlier [in]
word_optionsOptions for processing initial word hits [in]
ext_paramsParameters for the gapped extension [in]
hit_paramsParameters for saving the HSPs [in]
eff_len_paramsParameters for setting effective lengths [in]
psi_optionsOptions specific to PSI-BLAST [in]
db_optionsOptions for handling BLAST database [in]
hsp_streamPlaceholder for saving HSP lists [in]
diagnosticsReturn statistics containing numbers of hits on different stages of the search. Statistics saved only for the allocated parts of the structure. [in] [out]
interrupt_searchfunction callback to allow interruption of BLAST search [in, optional]
progress_infocontains information about the progress of the current BLAST search [in|out]

Definition at line 1335 of file blast_engine.c.

References ASSERT, BLAST_SequenceBlk::bases_offset, Blast_HSPListFree(), Blast_HSPListGetBitScores(), Blast_HSPListGetEvalues(), Blast_HSPListReapByQueryCoverage(), Blast_HSPListReapByRawScore(), Blast_HSPListReevaluateUngapped(), BLAST_LinkHsps(), BLAST_OneSubjectUpdateParameters(), Blast_ProgramIsMapping(), Blast_ProgramIsNucleotide(), Blast_ProgramIsRpsBlast(), BLAST_SEQSRC_EOF, BLAST_SEQSRC_ERROR, Blast_SubjectIsTranslated(), BLASTERR_INTERRUPTED, BLASTERR_SEQSRC, BlastHSPStreamWrite(), BlastInitialWordParametersFree(), BlastInitialWordParametersNew(), BlastLinkHSPParametersUpdate(), BlastSeqSrcGetAvgSeqLen(), BlastSeqSrcGetNumSeqs(), BlastSeqSrcGetSequence(), BlastSeqSrcGetTotLen(), BlastSeqSrcIteratorFree(), BlastSeqSrcIteratorNewEx(), BlastSeqSrcIteratorNext(), BlastSeqSrcReleaseSequence(), BlastSequenceBlkFree(), CalculateLinkHSPCutoffs(), LookupTableWrap::check_index_oid, CODON_LENGTH, BlastDiagnostics::cutoffs, DoAnchoredSearch(), eBlastEncodingProtein, BlastSeqSrcGetSeqArg::encoding, LookupTableWrap::end_search_indication, eNoResults, ePrelimSearch, FALSE, BlastScoringOptions::gapped_calculation, BLAST_SequenceBlk::gen_code_string, GenCodeSingletonFind(), BlastDatabaseOptions::genetic_code, BlastHitList::heapified, BlastHSPResults::hitlist_array, BlastHSPList::hspcnt, if(), LAST_VOL_IDX_INIT, BLAST_SequenceBlk::length, BlastHitSavingParameters::link_hsp_params, BlastHitList::low_score, BlastHitSavingParameters::low_score, BlastHitSavingOptions::low_score_perc, BlastScoreBlk::matrix_only_scoring, MAX, NULL, BlastHSPResults::num_queries, BlastSeqSrcGetSeqArg::oid, BlastExtensionParameters::options, BlastHitSavingParameters::options, BlastScoringParameters::options, query, BlastHSPStream::results, s_AdjustSubjectForSraSearch(), s_AdjustSubjectForTranslatedSraSearch(), s_Blast_HSPListReapByPrelimEvalue(), s_BlastCoreAuxStructFree(), s_BlastSearchEngineCore(), s_BlastSetUpAuxStructures(), s_FillReturnCutoffsInfo(), s_GetMinimumSubjSeqLen(), s_RPSPreliminarySearchEngine(), BlastGapAlignStruct::sbp, BlastSeqSrcGetSeqArg::seq, SBlastProgress::stage, and TRUE.

Referenced by Blast_RunFullSearch(), and Blast_RunPreliminarySearchWithInterrupt().

◆ Blast_RunFullSearch()

Int4 Blast_RunFullSearch ( EBlastProgramType  program_number,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
const BlastSeqSrc seq_src,
BlastScoreBlk sbp,
const BlastScoringOptions score_options,
LookupTableWrap lookup_wrap,
const BlastInitialWordOptions word_options,
const BlastExtensionOptions ext_options,
const BlastHitSavingOptions hit_options,
const BlastEffectiveLengthsOptions eff_len_options,
const PSIBlastOptions psi_options,
const BlastDatabaseOptions db_options,
BlastHSPStream hsp_stream,
const BlastRPSInfo rps_info,
BlastDiagnostics diagnostics,
BlastHSPResults **  results,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)

The high level function performing the BLAST search against a BLAST database after all the setup has been done.

Parameters
program_numberType of BLAST program [in]
queryThe query sequence [in]
query_infoAdditional query information [in]
seq_srcStructure containing BLAST database [in]
sbpScoring and statistical parameters [in]
score_optionsHit scoring options [in]
lookup_wrapThe lookup table, constructed earlier [in]
word_optionsOptions for processing initial word hits [in]
ext_optionsOptions and parameters for the gapped extension [in]
hit_optionsOptions for saving the HSPs [in]
eff_len_optionsOptions for setting effective lengths [in]
psi_optionsOptions specific to PSI-BLAST [in]
db_optionsOptions for handling BLAST database [in]
hsp_streamStructure for streaming results [in] [out]
rps_infoRPS BLAST auxiliary data structure [in]
diagnosticsReturn statistics containing numbers of hits on different stages of the search [out]
resultsResults of the BLAST search [out]
interrupt_searchfunction callback to allow interruption of BLAST search [in, optional]
progress_infocontains information about the progress of the current BLAST search [in|out]

Definition at line 1727 of file blast_engine.c.

References BLAST_ComputeTraceback(), BLAST_GapAlignSetUp(), BLAST_PreliminarySearchEngine(), Blast_ProgramIsPhiBlast(), BlastHSPStreamClose(), BlastUngappedStats::lookup_hits, LookupTableWrap::lut, NULL, SPHIPatternSearchBlk::num_patterns_db, query, results, s_BlastRunFullSearchCleanUp(), and BlastDiagnostics::ungapped_stat.

◆ Blast_RunPreliminarySearch()

Int2 Blast_RunPreliminarySearch ( EBlastProgramType  program,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
const BlastSeqSrc seq_src,
const BlastScoringOptions score_options,
BlastScoreBlk sbp,
LookupTableWrap lookup_wrap,
const BlastInitialWordOptions word_options,
const BlastExtensionOptions ext_options,
const BlastHitSavingOptions hit_options,
const BlastEffectiveLengthsOptions eff_len_options,
const PSIBlastOptions psi_options,
const BlastDatabaseOptions db_options,
BlastHSPStream hsp_stream,
BlastDiagnostics diagnostics 
)

The high level function performing the BLAST search against a BLAST database after all the setup has been done.

Parameters
program_numberType of BLAST program [in]
queryThe query sequence [in]
query_infoAdditional query information [in]
seq_srcStructure containing BLAST database [in]
score_optionsHit scoring options [in]
sbpScoring and statistical parameters [in]
lookup_wrapThe lookup table, constructed earlier [in]
word_optionsOptions for processing initial word hits [in]
ext_optionsOptions and parameters for the gapped extension [in]
hit_optionsOptions for saving the HSPs [in]
eff_len_optionsOptions for setting effective lengths [in]
psi_optionsOptions specific to PSI-BLAST [in]
db_optionsOptions for handling BLAST database [in]
hsp_streamStructure for streaming results [in] [out]
diagnosticsReturn statistics containing numbers of hits on different stages of the search [out]

Definition at line 1610 of file blast_engine.c.

References Blast_RunPreliminarySearchWithInterrupt(), NULL, and query.

◆ Blast_RunPreliminarySearchWithInterrupt()

Int2 Blast_RunPreliminarySearchWithInterrupt ( EBlastProgramType  program,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
const BlastSeqSrc seq_src,
const BlastScoringOptions score_options,
BlastScoreBlk sbp,
LookupTableWrap lookup_wrap,
const BlastInitialWordOptions word_options,
const BlastExtensionOptions ext_options,
const BlastHitSavingOptions hit_options,
const BlastEffectiveLengthsOptions eff_len_options,
const PSIBlastOptions psi_options,
const BlastDatabaseOptions db_options,
BlastHSPStream hsp_stream,
BlastDiagnostics diagnostics,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)

Same as above, with support for user interrupt function.

Parameters
program_numberType of BLAST program [in]
queryThe query sequence [in]
query_infoAdditional query information [in]
seq_srcStructure containing BLAST database [in]
score_optionsHit scoring options [in]
sbpScoring and statistical parameters [in]
lookup_wrapThe lookup table, constructed earlier [in]
word_optionsOptions for processing initial word hits [in]
ext_optionsOptions and parameters for the gapped extension [in]
hit_optionsOptions for saving the HSPs [in]
eff_len_optionsOptions for setting effective lengths [in]
psi_optionsOptions specific to PSI-BLAST [in]
db_optionsOptions for handling BLAST database [in]
hsp_streamStructure for streaming results [in] [out]
diagnosticsReturn statistics containing numbers of hits on different stages of the search [out]
interrupt_searchUser defined function to interrupt search [in]
progress_infoUser supplied data structure to aid interrupt [in]

Definition at line 1633 of file blast_engine.c.

References Blast_DiagnosticsFree(), Blast_DiagnosticsInit(), Blast_DiagnosticsUpdate(), BLAST_GapAlignSetUp(), BLAST_GapAlignStructFree(), BLAST_PreliminarySearchEngine(), BlastEffectiveLengthsParametersFree(), BlastExtensionParametersFree(), BlastHitSavingParametersFree(), BlastScoringParametersFree(), NULL, query, and BlastGapAlignStruct::sbp.

Referenced by Blast_RunPreliminarySearch(), and CPrelimSearchRunner::operator()().

◆ s_AdjustSubjectForSraSearch()

static void s_AdjustSubjectForSraSearch ( BlastHSPList hsp_list,
Uint1  offset 
)
static

◆ s_AdjustSubjectForTranslatedSraSearch()

static void s_AdjustSubjectForTranslatedSraSearch ( BlastHSPList hsp_list,
Uint1  offset,
Int4  length 
)
static

◆ s_AllocateSeqRange()

static void s_AllocateSeqRange ( BLAST_SequenceBlk subject,
SubjectSplitStruct backup,
Int4  num_seq_ranges 
)
static

◆ s_BackupSubject()

static void s_BackupSubject ( BLAST_SequenceBlk subject,
SubjectSplitStruct backup 
)
static

◆ s_Blast_HSPListReapByPrelimEvalue()

static Int2 s_Blast_HSPListReapByPrelimEvalue ( BlastHSPList hsp_list,
const BlastHitSavingParameters hit_params 
)
static

Discard the HSPs above the prelim e-value threshold from the HSP list.

Parameters
hsp_listList of HSPs for one subject sequence [in] [out]
hit_paramsParameters block containing the prelim e-value [in]

Definition at line 643 of file blast_engine.c.

References ASSERT, Blast_HSPFree(), BlastHSP::evalue, BlastHSPList::hsp_array, BlastHSPList::hspcnt, NULL, and BlastHitSavingParameters::prelim_evalue.

Referenced by BLAST_PreliminarySearchEngine(), and s_BlastSearchEngineCore().

◆ s_BlastCoreAuxStructFree()

static BlastCoreAuxStruct* s_BlastCoreAuxStructFree ( BlastCoreAuxStruct aux_struct)
static

◆ s_BlastRunFullSearchCleanUp()

static void s_BlastRunFullSearchCleanUp ( BlastGapAlignStruct gap_align,
BlastScoringParameters score_params,
BlastExtensionParameters ext_params,
BlastHitSavingParameters hit_params,
BlastEffectiveLengthsParameters eff_len_params 
)
static

Function to deallocate data structures allocated in Blast_RunFullSearch.

Definition at line 1710 of file blast_engine.c.

References BLAST_GapAlignStructFree(), BlastEffectiveLengthsParametersFree(), BlastExtensionParametersFree(), BlastHitSavingParametersFree(), BlastScoringParametersFree(), NULL, and BlastGapAlignStruct::sbp.

Referenced by Blast_RunFullSearch().

◆ s_BlastSearchEngineCore()

static Int2 s_BlastSearchEngineCore ( EBlastProgramType  program_number,
BLAST_SequenceBlk query,
BlastQueryInfo query_info_in,
BLAST_SequenceBlk subject,
LookupTableWrap lookup,
BlastGapAlignStruct gap_align,
const BlastScoringParameters score_params,
const BlastInitialWordParameters word_params,
const BlastExtensionParameters ext_params,
const BlastHitSavingParameters hit_params,
const BlastDatabaseOptions db_options,
BlastDiagnostics diagnostics,
BlastCoreAuxStruct aux_struct,
BlastHSPList **  hsp_list_out_ptr,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)
static

The core of the BLAST search: comparison between the (concatenated) query against one subject sequence.

Translation of the subject sequence into 6 frames is done inside, if necessary. If subject sequence is too long, it can be split into several chunks.

Parameters
program_numberBLAST program type [in]
queryQuery sequence structure [in]
query_info_inQuery information [in]
subjectSubject sequence structure [in]
lookupLookup table [in]
gap_alignStructure for gapped alignment information [in]
score_paramsScoring parameters [in]
word_paramsInitial word finding and ungapped extension parameters [in]
ext_paramsGapped extension parameters [in]
hit_paramsHit saving parameters [in]
db_optionsDatabase options [in]
diagnosticsHit counts and other diagnostics [in] [out]
aux_structStructure containing different auxiliary data and memory for the preliminary stage of the BLAST search [in]
hsp_list_out_ptrList of HSPs found for a given subject sequence [in]
interrupt_searchfunction callback to allow interruption of BLAST search [in, optional]
progress_infocontains information about the progress of the current BLAST search [in|out]

Definition at line 702 of file blast_engine.c.

References BLAST_ContextToFrame(), BLAST_GetAllTranslations(), Blast_HSPListAppend(), Blast_HSPListFree(), Blast_HSPListGetEvalues(), Blast_HSPListReapByRawScore(), BLAST_LinkHsps(), Blast_ProgramIsNucleotide(), Blast_ProgramIsPhiBlast(), Blast_ProgramIsRpsBlast(), Blast_SubjectIsTranslated(), BLASTERR_INTERRUPTED, BlastHspNumMax(), BlastQueryInfoNew(), context, ContextOffsetsToOffsetArray(), CONV_NUCL2PROT_COORDINATES, eBlastEncodingNcbi2na, eBlastTypeBlastx, eBlastTypeMapping, eBlastTypeRpsBlast, eBlastTypeRpsTblastn, eNoSubjMasking, FALSE, SubjectSplitStruct::full_range, BlastScoringOptions::gapped_calculation, BlastDiagnostics::gapped_stat, BlastScoreBlk::gbp, BlastGappedStats::good_extensions, BlastHSPList::hspcnt, i, BlastScoringOptions::is_ooframe, SSeqRange::left, BlastHitSavingParameters::link_hsp_params, lookup(), BlastScoreBlk::matrix_only_scoring, NULL, BlastRPSLookupTable::num_profiles, SubjectSplitStruct::num_seq_ranges, BlastGappedStats::num_seqs_passed, BlastHitSavingParameters::options, BlastScoringParameters::options, query, SSeqRange::right, BlastRPSLookupTable::rps_seq_offsets, s_AllocateSeqRange(), s_BackupSubject(), s_Blast_HSPListReapByPrelimEvalue(), s_BlastSearchEngineCoreCleanUp(), s_BlastSearchEngineOneContext(), s_RestoreSubject(), s_RPSOffsetArrayToContextOffsets(), BlastGapAlignStruct::sbp, BlastScoringParameters::scale_factor, SubjectSplitStruct::seq_ranges, SubjectSplitStruct::sequence, subject, and TRUE.

Referenced by BLAST_PreliminarySearchEngine(), and s_RPSPreliminarySearchEngine().

◆ s_BlastSearchEngineCoreCleanUp()

static void s_BlastSearchEngineCoreCleanUp ( EBlastProgramType  program_number,
BlastQueryInfo query_info,
const BlastQueryInfo query_info_in,
Uint1 translation_buffer,
Uint4 frame_offsets_a 
)
static

Clean up function for s_BlastSearchEngineCore.

Parameters
program_numberBLAST program type [in]
query_infoQuery information structure local to s_BlastSearchEngineCore, which may or may not be deallocated [in]
query_info_inQuery information [in]
translation_bufferbuffer containing translated sequence data [in]
frame_offsets_aFIXME

Definition at line 613 of file blast_engine.c.

References BlastQueryInfoFree(), eBlastTypeRpsTblastn, and sfree.

Referenced by s_BlastSearchEngineCore().

◆ s_BlastSearchEngineOneContext()

static Int2 s_BlastSearchEngineOneContext ( EBlastProgramType  program_number,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
BLAST_SequenceBlk subject,
Int4  orig_length,
LookupTableWrap lookup,
BlastGapAlignStruct gap_align,
const BlastScoringParameters score_params,
const BlastInitialWordParameters word_params,
const BlastExtensionParameters ext_params,
const BlastHitSavingParameters hit_params,
BlastDiagnostics diagnostics,
BlastCoreAuxStruct aux_struct,
BlastHSPList **  hsp_list_out_ptr,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)
static

Searches only one context of a database sequence, but does all chunks if it is split.

Parameters
program_numberBLAST program type [in]
queryQuery sequence structure [in]
query_infoQuery information [in]
subjectSubject sequence structure [in]
orig_lengthoriginal length of query before translation [in]
lookupLookup table [in]
gap_alignStructure for gapped alignment information [in]
score_paramsScoring parameters [in]
word_paramsInitial word finding and ungapped extension parameters [in]
ext_paramsGapped extension parameters [in]
hit_paramsHit saving parameters [in]
diagnosticsHit counts and other diagnostics [in] [out]
aux_structStructure containing different auxiliary data and memory for the preliminary stage of the BLAST search [in]
hsp_list_out_ptrList of HSPs found for a given subject sequence [out]
interrupt_searchfunction callback to allow interruption of BLAST search [in, optional]
progress_infocontains information about the progress of the current BLAST search [in|out]

Definition at line 423 of file blast_engine.c.

References ASSERT, BLAST_GetUngappedHSPList(), Blast_HSPListAdjustOddBlastnScores(), Blast_HSPListAdjustOffsets(), Blast_HSPListFree(), Blast_HSPListPurgeHSPsWithCommonEndpoints(), Blast_HSPListsMerge(), Blast_HSPListSortByScore(), Blast_HSPListSubjectBestHit(), Blast_ProgramIsMapping(), Blast_ProgramIsNucleotide(), Blast_SubjectIsTranslated(), BLASTERR_INTERRUPTED, BlastHspNumMax(), BlastInitHitListReset(), SBlastScoreMatrix::data, DBSEQ_CHUNK_OVERLAP, eBlastTypeRpsTblastn, BlastCoreAuxStruct::ewp, BlastScoringOptions::gapped_calculation, BlastDiagnostics::gapped_stat, BlastCoreAuxStruct::GetGappedScore, GetOffsetArraySize(), SubjectSplitStruct::hard_ranges, SubjectSplitStruct::hm_index, BlastHitSavingOptions::hsp_filt_opt, BlastHSPList::hspcnt, BlastCoreAuxStruct::init_hitlist, INT4_MIN, BlastScoringOptions::is_ooframe, BlastCoreAuxStruct::JumperGapped, SSeqRange::left, lookup(), BlastCoreAuxStruct::mapper_wordhits, BlastScoreBlk::matrix, BlastQueryInfo::max_length, NULL, SubjectSplitStruct::offset, BlastCoreAuxStruct::offset_pairs, BlastHSPList::oid, BlastHitSavingParameters::options, BlastScoringParameters::options, BlastGapAlignStruct::positionBased, BlastScoreBlk::psi_matrix, SPsiBlastScoreMatrix::pssm, query, s_BackupSubject(), s_GetNextSubjectChunk(), s_RestoreSubject(), s_TranslateHSPsToDNAPCoord(), BlastGapAlignStruct::sbp, SubjectSplitStruct::sequence, subject, BlastHSPFilteringOptions::subject_besthit_opts, SUBJECT_SPLIT_DONE, SUBJECT_SPLIT_NO_RANGE, SUBJECT_SPLIT_OK, BlastInitHitList::total, TRUE, BlastDiagnostics::ungapped_stat, and BlastCoreAuxStruct::WordFinder.

Referenced by s_BlastSearchEngineCore().

◆ s_BlastSetUpAuxStructures()

static Int2 s_BlastSetUpAuxStructures ( const BlastSeqSrc seq_src,
LookupTableWrap lookup_wrap,
const BlastInitialWordParameters word_params,
const BlastExtensionOptions ext_options,
const BlastHitSavingOptions hit_options,
BLAST_SequenceBlk query,
const BlastQueryInfo query_info,
BlastCoreAuxStruct **  aux_struct_ptr 
)
static

Setup of the auxiliary BLAST structures; also calculates internally used parameters from options.

Parameters
seq_srcSequence source information, with callbacks to get sequences, their lengths, etc. [in]
lookup_wrapLookup table, already constructed. [in]
word_paramsParameters for initial word finding and ungapped extension. [in]
ext_optionsoptions for gapped extension. [in]
hit_optionsoptions for saving hits. [in]
queryThe query sequence block [in]
aux_struct_ptrPlaceholder joining various auxiliary memory structures [out]

Definition at line 971 of file blast_engine.c.

References ASSERT, BLAST_GetGappedScore(), BLAST_InitHitListNew(), BLAST_SmithWatermanGetGappedScore(), BlastAaWordFinder(), BlastChooseNaExtend(), BlastChooseNucleotideScanSubject(), BlastChooseProteinScanSubject(), BlastExtendWordNew(), BlastNaWordFinder(), BlastRPSWordFinder(), calloc(), eAaLookupTable, eCompressedAaLookupTable, eIndexedMBLookupTable, eJumperWithTraceback, ePhiLookupTable, ePhiNaLookupTable, BlastExtensionOptions::ePrelimGapExt, eRPSLookupTable, eSmithWatermanScoreOnly, BlastCoreAuxStruct::ewp, BlastCoreAuxStruct::GetGappedScore, GetOffsetArraySize(), BlastCoreAuxStruct::init_hitlist, BlastCoreAuxStruct::JumperGapped, JumperNaWordFinder(), LookupTableWrap::lut_type, malloc(), BlastCoreAuxStruct::mapper_wordhits, MapperWordHitsNew(), MB_IndexedWordFinder(), NULL, BlastQueryInfo::num_queries, BlastCoreAuxStruct::offset_pairs, PHI_MAX_HIT, PHIBlastWordFinder(), PHIGetGappedScore(), query, LookupTableWrap::read_indexed_db, ShortRead_IndexedWordFinder(), and BlastCoreAuxStruct::WordFinder.

Referenced by BLAST_PreliminarySearchEngine().

◆ s_FillReturnCutoffsInfo()

static Int2 s_FillReturnCutoffsInfo ( BlastRawCutoffs return_cutoffs,
const BlastScoringParameters score_params,
const BlastInitialWordParameters word_params,
const BlastExtensionParameters ext_params,
const BlastHitSavingParameters hit_params 
)
static

Fills the output information about the cutoffs uses in a BLAST search.

Parameters
return_cutoffsStructure for saving cutoffs information [in] [out]
score_paramsScoring parameters, containing the scaling factor [in]
word_paramsInitial word parameters [in]
ext_paramsGapped extension parameters [in]
hit_paramsHit saving parameters [in]

Definition at line 928 of file blast_engine.c.

References BlastRawCutoffs::cutoff_score, BlastInitialWordParameters::cutoff_score_min, BlastHitSavingParameters::cutoff_score_min, BlastExtensionParameters::gap_x_dropoff, BlastExtensionParameters::gap_x_dropoff_final, if(), BlastScoringParameters::scale_factor, BlastRawCutoffs::ungapped_cutoff, BlastRawCutoffs::x_drop_gap, BlastRawCutoffs::x_drop_gap_final, BlastRawCutoffs::x_drop_ungapped, and BlastInitialWordParameters::x_dropoff_max.

Referenced by BLAST_PreliminarySearchEngine(), and s_RPSPreliminarySearchEngine().

◆ s_GetMinimumSubjSeqLen()

static Int4 s_GetMinimumSubjSeqLen ( LookupTableWrap lookup_wrap)
static

◆ s_GetNextSubjectChunk()

static Int2 s_GetNextSubjectChunk ( BLAST_SequenceBlk subject,
SubjectSplitStruct backup,
Boolean  is_nucleotide,
int  chunk_overlap 
)
static

◆ s_RestoreSubject()

static void s_RestoreSubject ( BLAST_SequenceBlk subject,
SubjectSplitStruct backup 
)
static

◆ s_RPSOffsetArrayToContextOffsets()

static void s_RPSOffsetArrayToContextOffsets ( BlastQueryInfo info,
Int4 new_offsets 
)
static

Set up context offsets for the auxiliary BlastQueryInfo structure that is created for the concatenated database in RPS BLAST search.

Calls the public function OffsetArrayToContextOffsets with a blastp program, because subjects are protein sequences. This guarantees that all frames are set to 0.

Parameters
infoThe preallocated structure [in] [out]
new_offsetsThe array context offsets to fill [in]

Definition at line 392 of file blast_engine.c.

References eBlastTypeBlastp, info, and OffsetArrayToContextOffsets().

Referenced by s_BlastSearchEngineCore().

◆ s_RPSPreliminarySearchEngine()

static Int2 s_RPSPreliminarySearchEngine ( EBlastProgramType  program_number,
BLAST_SequenceBlk query,
BlastQueryInfo query_info,
const BlastSeqSrc seq_src,
const BlastScoringParameters score_params,
LookupTableWrap lookup_wrap,
BlastCoreAuxStruct aux_struct,
const BlastInitialWordParameters word_params,
const BlastExtensionParameters ext_params,
BlastGapAlignStruct gap_align,
const BlastHitSavingParameters hit_params,
BlastHSPStream hsp_stream,
BlastDiagnostics diagnostics,
TInterruptFnPtr  interrupt_search,
SBlastProgress progress_info 
)
static

Performs the preliminary stage of an RPS BLAST search, after all set up has already been done.

Parameters
program_numberType of BLAST program [in]
queryThe query sequence [in]
query_infoAdditional query information [in]
seq_srcStructure containing BLAST database [in]
score_paramsHit scoring parameters [in]
lookup_wrapThe lookup table, constructed earlier [in]
aux_structWrapper for auxiliary structures used in preliminary search [in]
word_paramsParameters for processing initial word hits [in]
ext_paramsParameters for the gapped extension [in]
gap_alignStructure containing scoring block and memory allocated for gapped alignment. [in]
hit_paramsParameters for saving the HSPs [in]
hsp_streamPlaceholder for saving HSP lists [in]
diagnosticsReturn statistics containing numbers of hits on different stages of the search. Statistics saved only for the allocated parts of the structure. [in] [out]
interrupt_searchfunction callback to allow interruption of BLAST search [in, optional]
progress_infocontains information about the progress of the current BLAST search [in|out]

Definition at line 1096 of file blast_engine.c.

References Blast_GetOneQueryStructs(), Blast_HSPListFree(), Blast_ProgramIsRpsBlast(), BLASTERR_INTERRUPTED, BlastExtendWordFree(), BlastExtendWordNew(), BlastHSPStreamWrite(), BlastQueryInfoFree(), BlastSeqSrcGetNumSeqs(), BlastSeqSrcGetTotLen(), BlastSequenceBlkFree(), BlastDiagnostics::cutoffs, BlastCoreAuxStruct::ewp, FALSE, BlastHSPList::hspcnt, if(), INT4_MAX, lookup(), LookupTableWrap::lut, NULL, BlastQueryInfo::num_queries, BlastGapAlignStruct::positionBased, query, BlastHSPList::query_index, RPSPsiMatrixAttach(), RPSPsiMatrixDetach(), s_BlastSearchEngineCore(), s_FillReturnCutoffsInfo(), BlastGapAlignStruct::sbp, and TRUE.

Referenced by BLAST_PreliminarySearchEngine().

◆ s_TranslateHSPsToDNAPCoord()

static void s_TranslateHSPsToDNAPCoord ( EBlastProgramType  program,
BlastInitHitList init_hitlist,
const BlastQueryInfo query_info,
Int2  subject_frame,
Int4  subject_length,
Int4  offset 
)
static

Adjust HSP coordinates for out-of-frame gapped extension.

Parameters
programOne of blastx or tblastn [in]
init_hitlistList of hits after ungapped extension [in]
query_infoQuery information containing context offsets; needed for blastx only [in]
subject_frameFrame of the subject sequence; tblastn only [in]
subject_lengthLength of the original nucleotide subject sequence; tblastn only [in]
offsetShift in the subject sequence protein coordinates [in]

Definition at line 324 of file blast_engine.c.

References Blast_InitHitListSortByScore(), BSearchContextInfo(), CODON_LENGTH, BlastQueryInfo::contexts, eBlastTypeBlastx, BlastInitHitList::init_hsp_array, offset, BlastInitHSP::offsets, BlastOffsetPair::q_off, BlastUngappedData::q_start, BlastOffsetPair::qs_offsets, BlastContextInfo::query_offset, BlastOffsetPair::s_off, BlastUngappedData::s_start, BlastInitHitList::total, and BlastInitHSP::ungapped_data.

Referenced by s_BlastSearchEngineOneContext().

Variable Documentation

◆ kBlastMajorVersion

const int kBlastMajorVersion = 2

Major version.

Definition at line 81 of file blast_engine.c.

◆ kBlastMinorVersion

const int kBlastMinorVersion = 16

Minor version.

Definition at line 82 of file blast_engine.c.

◆ kBlastPatchVersion

const int kBlastPatchVersion = 1

Patch version.

Definition at line 83 of file blast_engine.c.

◆ SUBJECT_SPLIT_DONE

const Int2 SUBJECT_SPLIT_DONE = 0

return value indicating hitting the end

Definition at line 217 of file blast_engine.c.

Referenced by s_BlastSearchEngineOneContext(), and s_GetNextSubjectChunk().

◆ SUBJECT_SPLIT_NO_RANGE

const Int2 SUBJECT_SPLIT_NO_RANGE = 2

return value indicating all masked out

Definition at line 219 of file blast_engine.c.

Referenced by s_BlastSearchEngineOneContext(), and s_GetNextSubjectChunk().

◆ SUBJECT_SPLIT_OK

const Int2 SUBJECT_SPLIT_OK = 1

return value indicating OK

Definition at line 218 of file blast_engine.c.

Referenced by s_BlastSearchEngineOneContext(), and s_GetNextSubjectChunk().

Modified on Fri Sep 20 14:57:48 2024 by modify_doxy.py rev. 669887