NCBI C++ ToolKit
prelim_stage.cpp
Go to the documentation of this file.

Go to the SVN repository for this file.

1 /* ===========================================================================
2  *
3  * PUBLIC DOMAIN NOTICE
4  * National Center for Biotechnology Information
5  *
6  * This software/database is a "United States Government Work" under the
7  * terms of the United States Copyright Act. It was written as part of
8  * the author's official duties as a United States Government employee and
9  * thus cannot be copyrighted. This software/database is freely available
10  * to the public for use. The National Library of Medicine and the U.S.
11  * Government have not placed any restriction on its use or reproduction.
12  *
13  * Although all reasonable efforts have been taken to ensure the accuracy
14  * and reliability of the software and data, the NLM and the U.S.
15  * Government do not and cannot warrant the performance or results that
16  * may be obtained by using this software or data. The NLM and the U.S.
17  * Government disclaim all warranties, express or implied, including
18  * warranties of performance, merchantability or fitness for any particular
19  * purpose.
20  *
21  * Please cite the author in any work or product based on this material.
22  *
23  * ===========================================================================
24  *
25  * Author: Christiam Camacho
26  *
27  */
28 
29 /** @file prelim_stage.cpp
30  * NOTE: This file contains work in progress and the APIs are likely to change,
31  * please do not rely on them until this notice is removed.
32  */
33 
34 #include <ncbi_pch.hpp>
35 
37 #include <algo/blast/api/uniform_search.hpp> // for CSearchDatabase
41 
42 #include "prelim_search_runner.hpp"
43 #include "blast_aux_priv.hpp"
44 #include "psiblast_aux_priv.hpp"
45 #include "split_query_aux_priv.hpp"
46 #include "blast_seqalign.hpp"
47 #include <sstream>
48 
50 
51 /** @addtogroup AlgoBlast
52  *
53  * @{
54  */
55 
58 BEGIN_SCOPE(blast)
59 
60 
61 
63  CRef<CBlastOptions> options,
64  const CSearchDatabase& dbinfo)
65  : m_QueryFactory(query_factory), m_InternalData(new SInternalData),
66  m_Options(options), m_DbAdapter(NULL), m_DbInfo(&dbinfo)
67 {
69  CRef<TBlastSeqSrc> wrapped_src(new TBlastSeqSrc(seqsrc, BlastSeqSrcFree));
70  x_Init(query_factory, options, CRef<CPssmWithParameters>(), seqsrc);
71 
72  m_InternalData->m_SeqSrc = wrapped_src;
73 }
74 
76  CRef<CBlastOptions> options,
78  size_t num_threads)
79  : m_QueryFactory(query_factory), m_InternalData(new SInternalData),
80  m_Options(options), m_DbAdapter(db), m_DbInfo(NULL)
81 {
82  BlastSeqSrc* seqsrc = db->MakeSeqSrc();
83  x_Init(query_factory, options, CRef<CPssmWithParameters>(), seqsrc,
84  num_threads);
85  m_InternalData->m_SeqSrc.Reset(new TBlastSeqSrc(seqsrc, 0));
86  if (num_threads > 1) {
87  SetNumberOfThreads(num_threads);
88  }
89 }
90 
92  CRef<CBlastOptions> options,
93  BlastSeqSrc* seqsrc,
95  : m_QueryFactory(query_factory), m_InternalData(new SInternalData),
96  m_Options(options), m_DbAdapter(NULL), m_DbInfo(NULL)
97 {
98  x_Init(query_factory, options, pssm, seqsrc);
99  m_InternalData->m_SeqSrc.Reset(new TBlastSeqSrc(seqsrc, 0));
100 }
101 
102 void
104 {
105  const bool was_multithreaded = IsMultiThreaded();
106 
108  if (was_multithreaded != IsMultiThreaded()) {
114 
115  CRef<ILocalQueryData> query_data
117  unique_ptr<const CBlastOptionsMemento> opts_memento
119  if (IsMultiThreaded())
122  }
123 }
124 
125 void
127  CRef<CBlastOptions> options,
129  BlastSeqSrc* seqsrc,
130  size_t num_threads)
131 {
132  CRef<SBlastSetupData> setup_data =
133  BlastSetupPreliminarySearchEx(query_factory, options, pssm, seqsrc,
134  num_threads);
135  m_InternalData = setup_data->m_InternalData;
136  copy(setup_data->m_Masks.begin(), setup_data->m_Masks.end(),
137  back_inserter(m_MasksForAllQueries));
138  m_Messages = setup_data->m_Messages;
139 }
140 
141 int
143 {
144  typedef vector< CRef<CPrelimSearchThread> > TBlastThreads;
145  TBlastThreads the_threads(GetNumberOfThreads());
146 
147  unique_ptr<const CBlastOptionsMemento> opts_memento
149  _TRACE("Launching BLAST with " << GetNumberOfThreads() << " threads");
150 
151  // -RMH- This appears to be a problem right now. When used...this
152  // can cause all the work to go to a single thread! (-MN- This is fixed in SB-768)
154  static_cast<int>(GetNumberOfThreads()));
155 
156  // Create the threads ...
157  NON_CONST_ITERATE(TBlastThreads, thread, the_threads) {
158  thread->Reset(new CPrelimSearchThread(internal_data,
159  opts_memento.get()));
160  if (thread->Empty()) {
161  NCBI_THROW(CBlastSystemException, eOutOfMemory,
162  "Failed to create preliminary search thread");
163  }
164  }
165 
166  // Inform indexing library about the number of concurrent
167  // search threads.
168  //
170 
171  // ... launch the threads ...
172  NON_CONST_ITERATE(TBlastThreads, thread, the_threads) {
173  (*thread)->Run();
174  }
175 
176  // ... and wait for the threads to finish
177  Uint8 retv(0);
178  NON_CONST_ITERATE(TBlastThreads, thread, the_threads) {
179  void * result(0);
180  (*thread)->Join(&result);
181  if (result) {
182  // Thread is not really returning a pointer, it's actually
183  // retruning an int
184  retv = reinterpret_cast<Uint8> (result);
185  }
186  }
187 
189 
190  if (retv) {
191  NCBI_THROW(CBlastException, eCoreBlastError,
192  BlastErrorCode2String((Int2)retv));
193  }
194  return 0;
195 }
196 
199 {
200  if (! BlastSeqSrcGetNumSeqs(m_InternalData->m_SeqSrc->GetPointer())) {
201  string msg =
202  "Filtering resulted in an empty database.";
203 
206  msg);
207  }
208 
210 
214  int retval = 0;
215 
216  unique_ptr<const CBlastOptionsMemento> opts_memento
219  LookupTableOptions * lut_options = opts_memento->m_LutOpts;
220  BlastInitialWordOptions * word_options = opts_memento->m_InitWordOpts;
221 
222  // Query splitting data structure (used only if applicable)
224  CRef<CQuerySplitter> query_splitter = setup_data->m_QuerySplitter;
225  if (query_splitter->IsQuerySplit()) {
226 
227  CRef<CSplitQueryBlk> split_query_blk = query_splitter->Split();
228 
229  for (Uint4 i = 0; i < query_splitter->GetNumberOfChunks(); i++) {
230  try {
231  CRef<IQueryFactory> chunk_qf =
232  query_splitter->GetQueryFactoryForChunk(i);
233  _TRACE("Query chunk " << i << "/" <<
234  query_splitter->GetNumberOfChunks());
235  CRef<SInternalData> chunk_data =
239 
240  CRef<ILocalQueryData> query_data(
241  chunk_qf->MakeLocalQueryData( &*m_Options ) );
242  BLAST_SequenceBlk * chunk_queries =
243  query_data->GetSequenceBlk();
246  chunk_queries, lut_options, word_options );
247 
248  if (IsMultiThreaded()) {
249  x_LaunchMultiThreadedSearch(*chunk_data);
250  } else {
251  retval =
252  CPrelimSearchRunner(*chunk_data, opts_memento.get())();
253  if (retval) {
254  NCBI_THROW(CBlastException, eCoreBlastError,
255  BlastErrorCode2String(retval));
256  }
257  }
258 
259 
260  _ASSERT(chunk_data->m_HspStream->GetPointer());
261  BlastHSPStreamMerge(split_query_blk->GetCStruct(), i,
262  chunk_data->m_HspStream->GetPointer(),
263  m_InternalData->m_HspStream->GetPointer());
264  _ASSERT(m_InternalData->m_HspStream->GetPointer());
265  // free this as the query_splitter keeps a reference to the
266  // chunk factories, which in turn keep a reference to the local
267  // query data.
268  query_data->FlushSequenceData();
269  } catch (const CBlastException& e) {
270  // This error message is safe to ignore for a given chunk,
271  // because the chunks might end up producing a region of
272  // the query for which ungapped Karlin-Altschul blocks
273  // cannot be calculated
274  const string err_msg1("search cannot proceed due to errors "
275  "in all contexts/frames of query "
276  "sequences");
277  const string err_msg2(kBlastErrMsg_CantCalculateUngappedKAParams);
278  if (e.GetMsg().find(err_msg1) == NPOS && e.GetMsg().find(err_msg2) == NPOS) {
279  throw;
280  }
281  }
282  }
283 
284  // Restore the full query sequence for the traceback stage!
285  if (m_InternalData->m_Queries == NULL) {
286  CRef<ILocalQueryData> query_data
288  // Query masking info is calculated as a side-effect
289  CBlastScoreBlk sbp
290  (CSetupFactory::CreateScoreBlock(opts_memento.get(), query_data,
291  NULL, m_Messages, NULL, NULL));
292  m_InternalData->m_Queries = query_data->GetSequenceBlk();
293  }
294  } else {
295 
297  GetDbIndexRunSearchFn()( queries, lut_options, word_options );
298 
299  if (IsMultiThreaded()) {
301  } else {
302  retval = CPrelimSearchRunner(*m_InternalData, opts_memento.get())();
303  if (retval) {
304  NCBI_THROW(CBlastException, eCoreBlastError,
305  BlastErrorCode2String(retval));
306  }
307  }
308  }
309 
310  return m_InternalData;
311 }
312 
313 int
315 {
316  int retval = 0;
317  retval = BlastScoreBlkCheck(m_InternalData->m_ScoreBlk->GetPointer());
318  return retval;
319 }
320 
321 
324  Uint4 max_num_hsps,
325  bool* rm_hsps,
326  vector<bool> *rm_hsps_info) const
327 {
328  bool any_query_hsp_limited = false;
329  unique_ptr<const CBlastOptionsMemento> opts_memento
331 
333  Boolean *removed_hsps = new Boolean [ m_InternalData->m_QueryInfo->num_queries ];
334  SBlastHitsParameters* hit_param = NULL;
335  SBlastHitsParametersNew(opts_memento->m_HitSaveOpts,
336  opts_memento->m_ExtnOpts,
337  opts_memento->m_ScoringOpts,
338  &hit_param);
339  BlastHSPResults* retval =
342  hit_param,
343  max_num_hsps,
344  removed_hsps);
345  if( rm_hsps_info){
346  rm_hsps_info->reserve(m_InternalData->m_QueryInfo->num_queries );
347  for( int query_index = 0 ; query_index < m_InternalData->m_QueryInfo->num_queries ; query_index ++ ){
348  (*rm_hsps_info)[ query_index ] = removed_hsps[query_index] == FALSE ? false : true;
349  if( (*rm_hsps_info)[ query_index ] ) any_query_hsp_limited = true;
350  }
351  }
352  delete [] removed_hsps;
353  if( rm_hsps ) *rm_hsps = any_query_hsp_limited ;
354  // applications assume the HSPLists in the HSPResults are
355  // sorted in order of worsening best e-value
357  return retval;
358 }
359 
360 bool CBlastPrelimSearch::Run( vector<list<CRef<CStd_seg> > > & l )
361 {
362  Run();
363  return x_BuildStdSegList(l);
364 }
365 
366 void s_FixNumIdent(BlastHSPList *hsp_list, bool gapped_calculation)
367 {
368  BlastHSP* hsp;
369  int i;
370 
371  for (i=0; i < hsp_list->hspcnt; i++)
372  {
373  hsp = hsp_list->hsp_array[i];
374  if (gapped_calculation)
375  hsp->num_ident = -1;
376  }
377 }
378 
379 void s_GetBitScores(BlastHitList * hit_list, bool gapped_calculation, const BlastScoreBlk * sbp)
380 {
381 
382  for (int i = 0; i < hit_list->hsplist_count; i++)
383  {
384  BlastHSPList* hsp_list = hit_list->hsplist_array[i];
385  if (!hsp_list)
386  continue;
387 
388  Blast_HSPListGetBitScores(hsp_list, gapped_calculation, sbp);
389  s_FixNumIdent(hsp_list, gapped_calculation);
390  }
391 }
392 
393 // Results is trimmed by Blast Hits Save Options if set
395 {
396  if(m_InternalData->m_HspStream.Empty())
397  {
398  _TRACE("HSP Stream is empty");
399  return false;
400  }
401 
402  if(NULL != m_DbInfo)
403  {
405  }
406 
407  if(m_DbAdapter.Empty())
408  {
409  _TRACE("This method does not support CBlastPrelimSearch constructed with BlastSeqSrc");
410  return false;
411  }
412 
413  BlastHSPStream * hsp_stream = m_InternalData->m_HspStream->GetPointer();
414  if(NULL == hsp_stream)
415  {
416  _TRACE("NULL HSP Stream Pointer");
417  return false;
418  }
419 
420  IBlastSeqInfoSrc * s_seqInfoSrc = m_DbAdapter->MakeSeqInfoSrc();
421  EBlastProgramType program = hsp_stream->program;
422 
425 
426  if(NULL == results.GetPointer())
427  return false;
428 
429  int num_queries = results->num_queries;
430 
431  BlastHitList ** q_list_ptr = results->hitlist_array;
433  l.resize(num_queries);
434  const BlastScoreBlk * sbp = m_InternalData->m_ScoreBlk->GetPointer();
435  bool gapped_cal = m_Options->GetGappedMode();
436  for(int i=0; i < num_queries; i++)
437  {
438  CConstRef<CSeq_loc> query_loc = local_query_data->GetSeq_loc(i);
439  TSeqPos query_length = static_cast<TSeqPos>(local_query_data->GetSeqLength(i));
440  BlastHitList * hit_list = q_list_ptr[i];
441  if(NULL != hit_list)
442  {
443 
444  s_GetBitScores(hit_list, gapped_cal, sbp);
445  BLASTPrelminSearchHitListToStdSeg(program, hit_list, *query_loc, query_length, s_seqInfoSrc, l[i]);
446  }
447  }
448 
449  return true;
450 }
451 END_SCOPE(blast)
453 
454 /* @} */
455 
Auxiliary functions for BLAST.
Declarations for indexed blast databases.
BlastDiagnostics * Blast_DiagnosticsFree(BlastDiagnostics *diagnostics)
Free the BlastDiagnostics structure and all substructures.
Structures and API used for saving BLAST hits.
BlastHSPResults * Blast_HSPResultsFromHSPStreamWithLimitEx(struct BlastHSPStream *hsp_stream, Uint4 num_queries, SBlastHitsParameters *hit_param, Uint4 max_num_hsps, Boolean *removed_hsps)
As Blast_HSPResultsFromHSPStreamWithLimit, except accept and return array of Boolen flags specifying ...
Definition: blast_hits.c:3867
BlastHSPResults * Blast_HSPResultsFree(BlastHSPResults *results)
Deallocate memory for BLAST results.
Definition: blast_hits.c:3358
Int2 Blast_HSPResultsSortByEvalue(BlastHSPResults *results)
Sort each hit list in the BLAST results by best e-value.
Definition: blast_hits.c:3375
Int2 SBlastHitsParametersNew(const BlastHitSavingOptions *hit_options, const BlastExtensionOptions *ext_options, const BlastScoringOptions *scoring_options, SBlastHitsParameters **retval)
Sets up small structures used by blast_hit.c for saving HSPs.
Definition: blast_hits.c:75
Int2 Blast_HSPListGetBitScores(BlastHSPList *hsp_list, Boolean gapped_calculation, const BlastScoreBlk *sbp)
Calculate bit scores from raw scores in an HSP list.
Definition: blast_hits.c:1907
int BlastHSPStreamRegisterMTLock(BlastHSPStream *hsp_stream, MT_LOCK lock)
Attach a mutex lock to a stream to protect multiple access during writing.
int BlastHSPStreamMerge(SSplitQueryBlk *squery_blk, Uint4 chunk_num, BlastHSPStream *hsp_stream, BlastHSPStream *combined_hsp_stream)
Moves the HSPlists from an HSPStream into the list contained by a second HSPStream.
const char * kBlastErrMsg_CantCalculateUngappedKAParams
Definition: blast_message.c:38
@ eBlastSevWarning
Definition: blast_message.h:57
const int kBlastMessageNoContext
Declared in blast_message.h as extern const.
Definition: blast_message.c:36
C++ version of the initialization for the mutex locking interface.
EBlastProgramType
Defines the engine's notion of the different applications of the BLAST algorithm.
Definition: blast_program.h:72
Utility function to convert internal BLAST result structures into objects::CSeq_align_set objects.
Int4 BlastSeqSrcGetNumSeqs(const BlastSeqSrc *seq_src)
Get the number of sequences contained in the sequence source.
Definition: blast_seqsrc.c:177
BlastSeqSrc * BlastSeqSrcFree(BlastSeqSrc *seq_src)
Frees the BlastSeqSrc structure by invoking the destructor function set by the user-defined construct...
Definition: blast_seqsrc.c:112
void BlastSeqSrcSetNumberOfThreads(BlastSeqSrc *seq_src, int nthreads)
Set the number of threads for MT mode.
Definition: blast_seqsrc.c:168
void BlastSeqSrcResetChunkIterator(BlastSeqSrc *seq_src)
Reset the internal "bookmark" of the last chunk for iteration provided by this object.
Definition: blast_seqsrc.c:436
Definitions and prototypes used by blast_stat.c to calculate BLAST statistics.
int BlastScoreBlkCheck(BlastScoreBlk *sbp)
Check that score blk is valid, returns zero if it is.
Definition: blast_stat.c:853
Defines BLAST error codes (user errors included)
Encapsulates ALL the BLAST algorithm's options.
Search class to perform the preliminary stage of the BLAST search.
Wrapper class for BlastScoreBlk .
Definition: blast_aux.hpp:333
Defines system exceptions occurred while running BLAST.
CConstRef –.
Definition: ncbiobj.hpp:1266
Memento class to save, replace out, and restore the effective search space options of the CBlastOptio...
Interface to create a BlastSeqSrc suitable for use in CORE BLAST from a a variety of BLAST database/s...
Functor to run the preliminary stage of the BLAST search.
Thread class to run the preliminary stage of the BLAST search.
CRef –.
Definition: ncbiobj.hpp:618
Blast Search Subject.
Abstract base class to encapsulate retrieval of sequence identifiers.
Source of query sequence data for BLAST Provides an interface for search classes to retrieve sequence...
Definition: query_data.hpp:147
#define false
Definition: bool.h:36
virtual CConstRef< objects::CSeq_loc > GetSeq_loc(size_t index)=0
Get the Seq_loc for the sequence indicated by index.
CStructWrapper< BlastDiagnostics > TBlastDiagnostics
CRef< CBlastOptions > m_Options
int CheckInternalData()
Checks that internal data is valid.
size_t GetNumberOfThreads(void) const
Accessor for the number of threads to use.
void SplitQuery_SetEffectiveSearchSpace(CRef< CBlastOptions > options, CRef< IQueryFactory > full_query_fact, CRef< SInternalData > full_data)
this might supercede the function below...
void BLASTPrelminSearchHitListToStdSeg(EBlastProgramType program, BlastHitList *hit_list, const CSeq_loc &query_loc, TSeqPos query_length, const IBlastSeqInfoSrc *subject_seqinfo, list< CRef< CStd_seg > > &seg_list)
virtual void SetNumberOfThreads(size_t nthreads)
Mutator for the number of threads.
CRef< SBlastSetupData > BlastSetupPreliminarySearchEx(CRef< IQueryFactory > qf, CRef< CBlastOptions > options, CConstRef< CPssmWithParameters > pssm, BlastSeqSrc *seqsrc, size_t num_threads)
Extended interface to set up internal data structures used by the BLAST CORE engine.
static BlastDiagnostics * CreateDiagnosticsStructure()
Create and initialize the BlastDiagnostics structure for single-threaded applications.
SSplitQueryBlk * GetCStruct() const
Returns the C structure managed by objects of this class.
void s_FixNumIdent(BlastHSPList *hsp_list, bool gapped_calculation)
string BlastErrorCode2String(Int2 error_code)
Returns a string containing a human-readable interpretation of the error_code passed as this function...
CRef< CLocalDbAdapter > m_DbAdapter
bool IsQuerySplit() const
Determines whether the query sequence(s) are split or not.
Definition: split_query.hpp:88
const CBlastOptionsMemento * CreateSnapshot() const
Create a snapshot of the state of this object for internal use of its data structures (BLAST C++ APIs...
CRef< SInternalData > Run()
Borrow the internal data and results results.
CRef< TBlastSeqSrc > m_SeqSrc
The source of subject sequence data.
void AddMessageAllQueries(EBlastSeverity severity, int error_id, const string &message)
Add a message for all queries.
Definition: blast_aux.cpp:1056
bool IsMultiThreaded(void) const
Returns true if more than 1 thread is specified.
TSeqLocInfoVector m_MasksForAllQueries
Query masking information.
virtual BLAST_SequenceBlk * GetSequenceBlk()=0
Accessor for the BLAST_SequenceBlk structure.
BlastSeqSrc * MakeSeqSrc()
Retrieves or constructs the BlastSeqSrc.
CRef< ILocalQueryData > MakeLocalQueryData(const CBlastOptions *opts)
Creates and caches an ILocalQueryData.
Definition: query_data.cpp:52
int x_LaunchMultiThreadedSearch(SInternalData &internal_data)
Runs the preliminary search in multi-threaded mode.
CRef< SInternalData > SplitQuery_CreateChunkData(CRef< IQueryFactory > qf, CRef< CBlastOptions > options, CRef< SInternalData > full_data, size_t num_threads)
Function used by search class to retrieve a query factory for a given chunk.
TSearchMessages m_Messages
Warnings and error messages.
void s_GetBitScores(BlastHitList *hit_list, bool gapped_calculation, const BlastScoreBlk *sbp)
static BlastDiagnostics * CreateDiagnosticsStructureMT()
Create and initialize the BlastDiagnostics structure for multi-threaded applications.
virtual void SetNumberOfThreads(size_t nthreads)
@inheritDoc
CRef< SInternalData > m_InternalData
static BlastScoreBlk * CreateScoreBlock(const CBlastOptionsMemento *opts_memento, CRef< ILocalQueryData > query_data, BlastSeqLoc **lookup_segments, TSearchMessages &search_messages, TSeqLocInfoVector *masked_query_regions=NULL, const CBlastRPSInfo *rps_info=NULL)
Initializes the BlastScoreBlk.
static BlastSeqSrc * CreateBlastSeqSrc(const CSearchDatabase &db)
Create a BlastSeqSrc from a CSearchDatabase (uses CSeqDB)
CStructWrapper< BlastSeqSrc > TBlastSeqSrc
CRef< TBlastScoreBlk > m_ScoreBlk
BLAST score block structure.
TData * GetPointer()
The a pointer to the wrapped object.
BLAST_SequenceBlk * m_Queries
The query sequence data, these fields are "borrowed" from the query factory (which owns them)
IBlastSeqInfoSrc * MakeSeqInfoSrc()
Retrieves or constructs the IBlastSeqInfoSrc.
BlastQueryInfo * m_QueryInfo
The query information structure.
DbIndexSetUsingThreadsFnType GetDbIndexSetUsingThreadsFn()
Return the appropriate callback to set the concurrency state in the index structure.
CRef< IQueryFactory > GetQueryFactoryForChunk(Uint4 chunk_num)
Returns a IQueryFactory suitable to be executed by a BLAST search class.
void FlushSequenceData()
Frees the cached sequence data structure (as this is usually the larger data structure).
Definition: query_data.cpp:148
DbIndexRunSearchFnType GetDbIndexRunSearchFn()
Return the appropriate callback to run indexed seed search.
DbIndexSetNumThreadsFnType GetDbIndexSetNumThreadsFn()
Return the appropriate callback to set the number of threads in the index structure.
CRef< IQueryFactory > m_QueryFactory
Query factory is retained to ensure the lifetime of the data (queries) produced by it.
bool x_BuildStdSegList(vector< list< CRef< CStd_seg > > > &list)
CRef< CSplitQueryBlk > Split()
Split the query sequence(s)
CBlastPrelimSearch(CRef< IQueryFactory > query_factory, CRef< CBlastOptions > options, const CSearchDatabase &dbinfo)
Constructor which creates and manages a BLAST database handle for the caller.
const CSearchDatabase * m_DbInfo
CRef< TBlastDiagnostics > m_Diagnostics
Diagnostic output from preliminary and traceback stages.
virtual size_t GetSeqLength(size_t index)=0
Get the length of the sequence indicated by index.
void x_Init(CRef< IQueryFactory > query_factory, CRef< CBlastOptions > options, CConstRef< objects::CPssmWithParameters > pssm, BlastSeqSrc *seqsrc, size_t num_threads=1)
Internal initialization function Initializes internal data structures except the BlastSeqSrc.
Uint4 GetNumberOfChunks() const
Returns the number of chunks the query/queries will be split into.
Definition: split_query.hpp:85
MT_LOCK Blast_CMT_LOCKInit(void)
Initialize a mutex locking mechanism for BLAST.
CRef< TBlastHSPStream > m_HspStream
HSP output of the preliminary stage goes here.
BlastHSPResults * ComputeBlastHSPResults(BlastHSPStream *stream, Uint4 max_num_hsps=0, bool *rm_hsps=NULL, vector< bool > *rm_hsps_info=NULL) const
Return HSPs in a structure other than the HSPStream? Provide conversion? How to combine this with CBl...
bool GetGappedMode() const
Returns true if gapped BLAST is set, false otherwise.
unsigned int TSeqPos
Type for sequence locations and lengths.
Definition: ncbimisc.hpp:875
#define NON_CONST_ITERATE(Type, Var, Cont)
Non constant version of ITERATE macro.
Definition: ncbimisc.hpp:822
#define NULL
Definition: ncbistd.hpp:225
#define _TRACE(message)
Definition: ncbidbg.hpp:122
#define NCBI_THROW(exception_class, err_code, message)
Generic macro to throw an exception, given the exception class, error code and message string.
Definition: ncbiexpt.hpp:704
const string & GetMsg(void) const
Get message string.
Definition: ncbiexpt.cpp:461
TObjectType * GetPointer(void) THROWS_NONE
Get pointer,.
Definition: ncbiobj.hpp:998
void Reset(void)
Reset reference object.
Definition: ncbiobj.hpp:773
bool Empty(void) const THROWS_NONE
Check if CRef is empty – not pointing to any object, which means having a null value.
Definition: ncbiobj.hpp:719
int16_t Int2
2-byte (16-bit) signed integer
Definition: ncbitype.h:100
uint32_t Uint4
4-byte (32-bit) unsigned integer
Definition: ncbitype.h:103
uint64_t Uint8
8-byte (64-bit) unsigned integer
Definition: ncbitype.h:105
#define END_NCBI_SCOPE
End previously defined NCBI scope.
Definition: ncbistl.hpp:103
#define USING_SCOPE(ns)
Use the specified namespace.
Definition: ncbistl.hpp:78
#define END_SCOPE(ns)
End the previously defined scope.
Definition: ncbistl.hpp:75
#define BEGIN_NCBI_SCOPE
Define ncbi namespace.
Definition: ncbistl.hpp:100
#define BEGIN_SCOPE(ns)
Define a new scope.
Definition: ncbistl.hpp:72
#define NPOS
Definition: ncbistr.hpp:133
int i
Uint1 Boolean
bool replacment for C
Definition: ncbi_std.h:94
#define FALSE
bool replacment for C indicating false.
Definition: ncbi_std.h:101
void copy(Njn::Matrix< S > *matrix_, const Njn::Matrix< T > &matrix0_)
Definition: njn_matrix.hpp:613
Defines internal auxiliary functor object to run the preliminary stage of the BLAST search.
NOTE: This file contains work in progress and the APIs are likely to change, please do not rely on th...
Declarations of auxiliary functions/classes for PSI-BLAST.
Auxiliary functions and classes to assist in query splitting.
Structure to hold a sequence.
Definition: blast_def.h:242
Return statistics from the BLAST search.
The structure to hold all HSPs for a given sequence after the gapped alignment.
Definition: blast_hits.h:153
Int4 hspcnt
Number of HSPs saved.
Definition: blast_hits.h:158
BlastHSP ** hsp_array
Array of pointers to individual HSPs.
Definition: blast_hits.h:157
The structure to contain all BLAST results, for multiple queries.
Definition: blast_hits.h:183
Default implementation of BlastHSPStream.
EBlastProgramType program
BLAST program type.
Structure holding all information about an HSP.
Definition: blast_hits.h:126
Int4 num_ident
Number of identical base pairs in this HSP.
Definition: blast_hits.h:128
The structure to contain all BLAST results for one query sequence.
Definition: blast_hits.h:169
BlastHSPList ** hsplist_array
Array of HSP lists for individual database hits.
Definition: blast_hits.h:176
Int4 hsplist_count
Filled size of the HSP lists array.
Definition: blast_hits.h:170
Options needed for initial word finding and processing.
int num_queries
Number of query sequences.
Structure used for scoring calculations.
Definition: blast_stat.h:177
Complete type definition of Blast Sequence Source ADT.
Definition: blast_seqsrc.c:43
Options needed to construct a lookup table Also needed: query sequence and query length.
Keeps prelim_hitlist_size and HitSavingOptions together, mostly for use by hspstream.
Definition: blast_hits.h:57
Return type of BlastSetupPreliminarySearch.
Lightweight wrapper to enclose C structures needed for running the preliminary and traceback stages o...
#define _ASSERT
else result
Definition: token2.c:20
Uniform BLAST Search Interface.
#define const
Definition: zconf.h:232
Modified on Tue Apr 30 06:40:55 2024 by modify_doxy.py rev. 669887