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

Go to the SVN repository for this file.

1 
2 /* ===========================================================================
3  *
4  * PUBLIC DOMAIN NOTICE
5  * National Center for Biotechnology Information
6  *
7  * This software/database is a "United States Government Work" under the
8  * terms of the United States Copyright Act. It was written as part of
9  * the author's official duties as a United States Government employee and
10  * thus cannot be copyrighted. This software/database is freely available
11  * to the public for use. The National Library of Medicine and the U.S.
12  * Government have not placed any restriction on its use or reproduction.
13  *
14  * Although all reasonable efforts have been taken to ensure the accuracy
15  * and reliability of the software and data, the NLM and the U.S.
16  * Government do not and cannot warrant the performance or results that
17  * may be obtained by using this software or data. The NLM and the U.S.
18  * Government disclaim all warranties, express or implied, including
19  * warranties of performance, merchantability or fitness for any particular
20  * purpose.
21  *
22  * Please cite the author in any work or product based on this material.
23  *
24  * ===========================================================================
25  *
26  * Author: Christiam Camacho
27  *
28  */
29 
30 /// @file blast_aux_priv.cpp
31 /// Implements various auxiliary (private) functions for BLAST
32 
33 #include <ncbi_pch.hpp>
34 #include "blast_aux_priv.hpp"
41 #include "psiblast_aux_priv.hpp"
42 #include "blast_memento_priv.hpp"
43 
47 
48 /** @addtogroup AlgoBlast
49  *
50  * @{
51  */
52 
56 
59 {
60  _ASSERT(!seqids.empty());
61  CRef<CSeq_loc> retval(new CSeq_loc);
62  retval->SetWhole().Assign(**seqids.begin());
63  return retval;
64 }
65 
66 void
68  const BlastQueryInfo* query_info,
69  TSearchMessages& messages)
70 {
71  if ( !blmsg || !query_info ) {
72  return;
73  }
74 
75  if (messages.size() != (size_t) query_info->num_queries) {
76  messages.resize(query_info->num_queries);
77  }
78 
79  const BlastContextInfo* kCtxInfo = query_info->contexts;
80 
81  // First copy the errors...
82  for (; blmsg; blmsg = blmsg->next)
83  {
84  const int kContext = blmsg->context;
85  _ASSERT(blmsg->message);
86  string msg(blmsg->message);
87 
88  if (kContext != kBlastMessageNoContext) {
89  // applies only to a single query
90  const int kQueryIndex = kCtxInfo[kContext].query_index;
92  kQueryIndex, msg));
93  messages[kCtxInfo[kContext].query_index].push_back(sm);
94  } else {
95  // applies to all queries
98  msg));
99  NON_CONST_ITERATE(TSearchMessages, query_messages, messages) {
100  query_messages->push_back(sm);
101  }
102  }
103 
104 
105  }
106 
107  // ... then remove duplicate error messages
108  messages.RemoveDuplicates();
109 }
110 
111 string
113 {
114  Blast_Message* blast_msg = NULL;
115  Blast_PerrorWithLocation(&blast_msg, error_code, kBlastMessageNoContext);
116  string retval(blast_msg != NULL ? blast_msg->message : kEmptyStr);
117  blast_msg = Blast_MessageFree(blast_msg);
118  return retval;
119 }
120 
123  CRef<CBlastOptions> options,
124  size_t num_threads /* = 1 */)
125 {
126  return BlastSetupPreliminarySearchEx(query_factory, options,
128  NULL, num_threads);
129 }
130 
133  CRef<CBlastOptions> options,
135  BlastSeqSrc* seqsrc,
136  size_t num_threads)
137 {
138  CRef<SBlastSetupData> retval(new SBlastSetupData(qf, options));
139  TSearchMessages m;
140  options->Validate();
141  bool is_multi_threaded = num_threads > 1;
142 
143  // 0. Initialize the megablast database index.
144  if (options->GetUseIndex()) {
147  }
148 
149  // 1. Initialize the query data (borrow it from the factory)
150  CRef<ILocalQueryData> query_data(qf->MakeLocalQueryData(&*options));
151  retval->m_InternalData->m_Queries = query_data->GetSequenceBlk();
152  retval->m_InternalData->m_QueryInfo = query_data->GetQueryInfo();
153  // get any warning messages from instantiating the queries
154  query_data->GetMessages(m);
155  retval->m_Messages.resize(query_data->GetNumQueries());
156  retval->m_Messages.Combine(m);
157 
158  // 2. Take care of any rps information
159  if (Blast_ProgramIsRpsBlast(options->GetProgramType())) {
160  const char* name = BlastSeqSrcGetName(seqsrc);
161  const string rps_dbname(name ? name : "");
162  retval->m_InternalData->m_RpsData =
163  CSetupFactory::CreateRpsStructures(rps_dbname, options);
164  }
165 
166  // 3. Create the options memento
167  auto_ptr<const CBlastOptionsMemento> opts_memento
168  (options->CreateSnapshot());
169 
170  // 4. Create the BlastScoreBlk
171  BlastSeqLoc* lookup_segments = NULL;
172  BlastScoreBlk* sbp = NULL;
173  try {
174  sbp =
175  CSetupFactory::CreateScoreBlock(opts_memento.get(), query_data,
176  &lookup_segments, retval->m_Messages,
177  &retval->m_Masks,
178  retval->m_InternalData->m_RpsData);
179  } catch (CBlastException & e) {
180  const string kCatchThisError (kBlastErrMsg_CantCalculateUngappedKAParams);
181  if(e.GetMsg() == kCatchThisError) {
182  return retval;
183  }
184  else {
185  NCBI_RETHROW(e, CBlastException, eCoreBlastError, e.GetMsg());
186  }
187  }
188  CRef< CBlastSeqLocWrap > lookup_segments_wrap(
189  new CBlastSeqLocWrap( lookup_segments ) );
190  retval->m_InternalData->m_ScoreBlk.Reset
191  (new TBlastScoreBlk(sbp, BlastScoreBlkFree));
192  if (pssm.NotEmpty()) {
193  if (query_data->GetNumQueries() > 1) {
194  NCBI_THROW(CBlastException, eNotSupported,
195  "Multiple queries cannot be specified with a PSSM");
196  }
197  PsiBlastSetupScoreBlock(sbp, pssm, retval->m_Messages, options);
198  }
199 
200  // Initialize callbacks needed for indexed search.
201  //
202  if (options->GetUseIndex()) {
204  }
205 
206  // 5. Create the lookup table
207  if ( !retval->m_QuerySplitter->IsQuerySplit() ) {
208  LookupTableWrap* lut =
209  CSetupFactory::CreateLookupTable(query_data, opts_memento.get(),
210  sbp, lookup_segments_wrap,
211  retval->m_InternalData->m_RpsData,
212  seqsrc,
213  num_threads);
214  retval->m_InternalData->m_LookupTable.Reset
216  }
217 
218  // 6. Create diagnostics
219  BlastDiagnostics* diags = is_multi_threaded
222  retval->m_InternalData->m_Diagnostics.Reset
224 
225  // 7. Create the HSP stream
226  BlastHSPStream* hsp_stream =
227  CSetupFactory::CreateHspStream(opts_memento.get(),
228  query_data->GetNumQueries(),
229  CSetupFactory::CreateHspWriter(opts_memento.get(),
230  retval->m_InternalData->m_Queries,
231  query_data->GetQueryInfo()));
232 
233  if (is_multi_threaded)
235 
236  // 8. Register a traceback HSP Pipe(s)
237  BlastHSPStreamRegisterPipe(hsp_stream,
238  CSetupFactory::CreateHspPipe(opts_memento.get(),
239  query_data->GetQueryInfo()),
241 
242  retval->m_InternalData->m_HspStream.Reset
243  (new TBlastHSPStream(hsp_stream, BlastHSPStreamFree));
244 
245  // 8. Get errors/warnings
246  query_data->GetMessages(m);
247  retval->m_Messages.Combine(m);
248 
249  if (retval->m_QuerySplitter->IsQuerySplit()) {
250  // We don't need the full sequence for the preliminary stage, so we
251  // free it and NULL out references to it (this MUST be restored prior
252  // to the traceback stage)
253  query_data->FlushSequenceData();
254  retval->m_InternalData->m_Queries = NULL;
255  }
256 
257  retval->m_InternalData->m_FnInterrupt = NULL;
259  return retval;
260 }
261 
262 
263 void
265  const vector< CConstRef<CSeq_id> >& query_ids,
266  const BlastScoreBlk* sbp,
267  const BlastQueryInfo* qinfo,
268  const TSeqAlignVector& alignments,
269  const EResultType result_type,
271 {
272  retval.clear();
273 
274  if (Blast_ProgramIsPhiBlast(program)) {
275  CRef<CBlastAncillaryData> s(new CBlastAncillaryData(program, 0, sbp,
276  qinfo));
277 
278  for(unsigned i = 0; i < alignments.size(); i++) {
279  retval.push_back(s);
280  }
281  } else {
282  if (result_type == ncbi::blast::eSequenceComparison) {
283  const size_t num_subjects = alignments.size()/query_ids.size();
284  for(size_t i = 0; i < alignments.size(); i += num_subjects) {
286  (new CBlastAncillaryData(program, i/num_subjects, sbp,
287  qinfo));
288  for (size_t j = 0; j < num_subjects; j++) {
289  retval.push_back(s);
290  }
291  }
292  } else {
293  for(size_t i = 0; i < alignments.size(); i++) {
295  sbp,
296  qinfo));
297  retval.push_back(s);
298  }
299  }
300  }
301 
302 }
303 
306  const BlastScoreBlk* sbp,
307  const BlastQueryInfo* qinfo,
308  EBlastProgramType program,
309  const TSeqAlignVector& alignments,
310  TSearchMessages& messages,
311  const vector<TSeqLocInfoVector>& subj_masks,
312  const TSeqLocInfoVector* query_masks,
313  const EResultType result_type)
314 {
315  const bool is_phi = !!Blast_ProgramIsPhiBlast(program);
316 
317  // Collect query Seq-locs
318 
319  vector< CConstRef<CSeq_id> > qlocs;
320 
321  if (is_phi) {
322  qlocs.assign(alignments.size(), query_ids.front());
323  } else {
324  if (result_type == ncbi::blast::eSequenceComparison)
325  {
326  const size_t num_subjects = alignments.size()/query_ids.size();
327  for (size_t i = 0; i < alignments.size(); i += num_subjects) {
328  for (size_t j = 0; j < num_subjects; j++) {
329  qlocs.push_back(query_ids[i/num_subjects]);
330  }
331  }
332  }
333  else
334  copy(query_ids.begin(), query_ids.end(), back_inserter(qlocs));
335  }
336 
337  // Collect ancillary data
338 
339  CSearchResultSet::TAncillaryVector ancillary_data;
340  BuildBlastAncillaryData(program, query_ids, sbp, qinfo, alignments,
341  result_type, ancillary_data);
342 
343  // The preliminary stage also produces errors and warnings; they
344  // should be copied from that code to this class somehow, and
345  // returned here if they have not been returned or reported yet.
346 
347  if (messages.size() < alignments.size()) {
348  messages.resize(alignments.size());
349  }
350 
351  // N.B.: the number of query masks for bl2seq will be adjusted in
352  // CSearchResultSet::SetFilteredQueryRegions
353  const SPHIQueryInfo* phi_query_info = is_phi ? qinfo->pattern_info : NULL;
354  CRef<CSearchResultSet> retval(new CSearchResultSet(qlocs, alignments,
355  messages,
356  ancillary_data,
357  query_masks,
358  result_type,
359  phi_query_info));
360  if (subj_masks.size() == retval->size()) {
361  for (CSearchResultSet::size_type i = 0; i < retval->size(); i++) {
362  (*retval)[i].SetSubjectMasks(subj_masks[i]);
363  }
364  }
365  return retval;
366 }
367 
371  bool assume_both_strands)
372 {
373  if (sloc_in.Empty() ||
374  sloc_in->Which() == CSeq_loc::e_not_set ||
375  sloc_in->IsEmpty() ||
376  sloc_in->IsNull()) {
377  return TMaskedQueryRegions();
378  }
379 
380  CConstRef<CSeq_loc> sloc = sloc_in;
381 
382  if (sloc_in->IsInt()) {
384  iv( const_cast<CSeq_interval *>(& sloc_in->GetInt()) );
385 
386  CRef<CSeq_loc> nsloc(new CSeq_loc);
387  nsloc->SetPacked_int().Set().push_back(iv);
388 
389  sloc.Reset(&*nsloc);
390  }
391 
392  if (! sloc->IsPacked_int()) {
393  NCBI_THROW(CBlastException, eNotSupported,
394  "Unsupported Seq-loc type used for mask");
395  }
396 
397  const objects::CPacked_seqint & psi = sloc->GetPacked_int();
398 
400 
401  ITERATE(list< CRef< objects::CSeq_interval > >, iter, psi.Get()) {
402  objects::CSeq_interval * iv =
403  const_cast<objects::CSeq_interval*>(& (**iter));
404 
405  if (Blast_QueryIsProtein(prog)) {
406  int fr = (int) CSeqLocInfo::eFrameNotSet;
407  mqr.push_back(CRef<CSeqLocInfo>(new CSeqLocInfo(iv, fr)));
408  } else {
409  bool do_pos = false;
410  bool do_neg = false;
411 
412  if (iv->CanGetStrand()) {
413  switch(iv->GetStrand()) {
414  case objects::eNa_strand_plus:
415  do_pos = true;
416  break;
417 
418  case objects::eNa_strand_minus:
419  do_neg = true;
420  break;
421 
422  case objects::eNa_strand_both:
423  do_pos = true;
424  do_neg = true;
425  break;
426 
427  default:
428  NCBI_THROW(CBlastException, eNotSupported,
429  "Unsupported strand type used for query");
430  }
431  } else {
432  // intervals with no strand assignment will use both.
433  do_pos = do_neg = true;
434  }
435 
436  // deliberately override the strand option above, if so requested
437  if (assume_both_strands) {
438  do_pos = do_neg = true;
439  }
440 
441  if (do_pos) {
442  int fr = (int) CSeqLocInfo::eFramePlus1;
443  mqr.push_back(CRef<CSeqLocInfo>(new CSeqLocInfo(iv, fr)));
444  }
445 
446  // No reversal is done here. Tt seems that the code (in core)
447  // that applies the mask reverses it. Whether this is an
448  // accidental or designed is not clear to me, but for now this
449  // will remain as is.
450 
451  if (do_neg) {
452  int fr = (int) CSeqLocInfo::eFrameMinus1;
453  mqr.push_back(CRef<CSeqLocInfo>(new CSeqLocInfo(iv, fr)));
454  }
455  }
456  }
457 
458  return mqr;
459 }
460 
463 {
464  if (sloc.empty()) {
465  return CRef<objects::CSeq_loc>();
466  }
467 
470  if (psi.NotEmpty()) {
471  retval.Reset(new objects::CSeq_loc);
472  retval->SetPacked_int(*psi);
473  }
474  return retval;
475 }
476 
479 
480 /* @} */
TSearchMessages m_Messages
void SetUpDbIndexCallbacks(void)
void SetPacked_int(TPacked_int &v)
Definition: Seq_loc.hpp:968
vector< value_type >::size_type size_type
size_type type definition
bool GetUseIndex() const
static LookupTableWrap * CreateLookupTable(CRef< ILocalQueryData > query_data, const CBlastOptionsMemento *opts_memento, BlastScoreBlk *score_blk, CRef< CBlastSeqLocWrap > lookup_segments, const CBlastRPSInfo *rps_info=NULL, BlastSeqSrc *seqsrc=NULL, size_t num_threads=1)
Initialize the lookup table.
#define NON_CONST_ITERATE(Type, Var, Cont)
Non constant version of ITERATE macro.
Definition: ncbimisc.hpp:783
Class used to return ancillary data from a blast search, i.e.
CConstRef –.
Definition: ncbiobj.hpp:1192
void BuildBlastAncillaryData(EBlastProgramType program, const vector< CConstRef< CSeq_id > > &query_ids, const BlastScoreBlk *sbp, const BlastQueryInfo *qinfo, const TSeqAlignVector &alignments, const EResultType result_type, CSearchResultSet::TAncillaryVector &retval)
Builds an CSearchResultSet::TAncillaryVector.
MT_LOCK Blast_CMT_LOCKInit(void)
Initialize a mutex locking mechanism for BLAST.
Declares the BLAST exception class.
static const unsigned char msg[]
Definition: ccm.c:378
TMaskedQueryRegions PackedSeqLocToMaskedQueryRegions(CConstRef< objects::CSeq_loc > sloc, EBlastProgramType program, bool assume_both_strands=false)
Auxiliary function to convert a Seq-loc describing masked query regions to a TMaskedQueryRegions obje...
static BlastDiagnostics * CreateDiagnosticsStructureMT()
Create and initialize the BlastDiagnostics structure for multi-threaded applications.
bool NotEmpty(void) const THROWS_NONE
Check if CConstRef is not empty – pointing to an object and has a non-null value.
Definition: ncbiobj.hpp:1296
LookupTableWrap * LookupTableWrapFree(LookupTableWrap *lookup)
Deallocate memory for the lookup table.
Definition: lookup_wrap.c:196
#define END_SCOPE(ns)
End the previously defined scope.
Definition: ncbistl.hpp:73
struct Blast_Message * next
next message in this list
Definition: blast_message.h:71
CRef< TBlastHSPStream > m_HspStream
HSP output of the preliminary stage goes here.
CRef< CBlastRPSInfo > m_RpsData
The RPS-BLAST related data.
const TPacked_int & GetPacked_int(void) const
Get the variant data.
Definition: Seq_loc_.cpp:216
structure for seqloc info
Definition: seqlocinfo.hpp:48
int num_queries
Number of query sequences.
CStructWrapper< BlastHSPStream > TBlastHSPStream
Search Results for All Queries.
Complete type definition of Blast Sequence Source ADT.
Definition: blast_seqsrc.c:43
const char * kBlastErrMsg_CantCalculateUngappedKAParams
Definition: blast_message.c:38
CRef< objects::CSeq_loc > MaskedQueryRegionsToPackedSeqLoc(const TMaskedQueryRegions &sloc)
Interface to build a CSeq-loc from a TMaskedQueryRegion; note that conversion conversion in this dire...
void Combine(const TSearchMessages &other_msgs)
Combine another set of search messages with this one.
Definition: blast_aux.cpp:1028
CRef< CQuerySplitter > m_QuerySplitter
static BlastHSPStream * CreateHspStream(const CBlastOptionsMemento *opts_memento, size_t number_of_queries, BlastHSPWriter *writer)
Create and initialize the BlastHSPStream structure.
Defines BLAST error codes (user errors included)
Wrapper for BlastSeqLoc structure.
Auxiliary functions for BLAST.
Structure used for scoring calculations.
Definition: blast_stat.h:177
void Blast_Message2TSearchMessages(const Blast_Message *blmsg, const BlastQueryInfo *query_info, TSearchMessages &messages)
Converts the Blast_Message structure into a TSearchMessages object.
void PsiBlastSetupScoreBlock(BlastScoreBlk *score_blk, CConstRef< objects::CPssmWithParameters > pssm, TSearchMessages &messages, CConstRef< CBlastOptions > options)
Setup CORE BLAST score block structure with data from the scoremat PSSM.
EResultType
Specifies the style of Seq-aligns that should be built from the internal BLAST data structures...
static BlastDiagnostics * CreateDiagnosticsStructure()
Create and initialize the BlastDiagnostics structure for single-threaded applications.
static const unsigned char iv[]
Definition: ccm.c:367
BlastContextInfo * contexts
Information per context.
vector< TMaskedQueryRegions > TSeqLocInfoVector
Collection of masked regions for all queries in a BLAST search.
Definition: seqlocinfo.hpp:139
Declarations of auxiliary functions/classes for PSI-BLAST.
CRef< TLookupTableWrap > m_LookupTable
Lookup table, usually only needed in the preliminary stage of the search, but for PHI-BLAST it's also...
#define NULL
Definition: ncbistd.hpp:225
BlastQueryInfo * m_QueryInfo
The query information structure.
void Reset(void)
Reset reference object.
Definition: ncbiobj.hpp:1343
#define kEmptyStr
Definition: ncbistr.hpp:120
Blast_Message * Blast_MessageFree(Blast_Message *blast_msg)
Deallocates message memory.
Definition: blast_message.c:80
BLAST_SequenceBlk * m_Queries
The query sequence data, these fields are "borrowed" from the query factory (which owns them) ...
vector< CRef< objects::CSeq_align_set > > TSeqAlignVector
Vector of Seq-align-sets.
Boolean Blast_ProgramIsPhiBlast(EBlastProgramType p)
Returns true if program is PHI-BLAST (i.e.
Definition: blast_program.c:70
#define ITERATE(Type, Var, Cont)
ITERATE macro to sequence through container elements.
Definition: ncbimisc.hpp:776
int i
CRef< CSearchResultSet > BlastBuildSearchResultSet(const vector< CConstRef< CSeq_id > > &query_ids, const BlastScoreBlk *sbp, const BlastQueryInfo *qinfo, EBlastProgramType program, const TSeqAlignVector &alignments, TSearchMessages &messages, const vector< TSeqLocInfoVector > &subj_masks, const TSeqLocInfoVector *query_masks, const EResultType result_type)
Build a CSearchResultSet from internal BLAST data structures.
CRef< TBlastScoreBlk > m_ScoreBlk
BLAST score block structure.
Declaration of ADT to save and retrieve lists of HSPs in the BLAST engine.
Wrapper class for SBlastProgress .
Definition: blast_aux.hpp:331
void copy(Njn::Matrix< S > *matrix_, const Njn::Matrix< T > &matrix0_)
Definition: njn_matrix.hpp:612
Int4 query_index
Index of query (same for all frames)
#define END_NCBI_SCOPE
End previously defined NCBI scope.
Definition: ncbistl.hpp:101
static void InitializeMegablastDbIndex(CRef< CBlastOptions > options)
Initialize a megablast BLAST database index.
CRef< objects::CPacked_seqint > ConvertToCPacked_seqint() const
Converts this object to a CPacked_seqint (this is the convention used to encode masking locations in ...
Definition: seqlocinfo.cpp:127
BlastScoreBlk * BlastScoreBlkFree(BlastScoreBlk *sbp)
Deallocates BlastScoreBlk as well as all associated structures.
Definition: blast_stat.c:965
Error or Warning Message from search.
EBlastSeverity severity
severity code
Definition: blast_message.h:72
int BlastHSPStreamRegisterPipe(BlastHSPStream *hsp_stream, BlastHSPPipe *pipe, EBlastStage stage)
Insert the user-specified pipe to the *end* of the pipeline.
CRef< SInternalData > m_InternalData
vector< CRef< CBlastAncillaryData > > TAncillaryVector
typedef for a vector of CRef
Used to hold a set of positions, mostly used for filtering.
Definition: blast_def.h:204
CStructWrapper< LookupTableWrap > TLookupTableWrap
Return statistics from the BLAST search.
Declaration of ADT to retrieve sequences for the BLAST engine.
static CRef< CBlastRPSInfo > CreateRpsStructures(const string &rps_dbname, CRef< CBlastOptions > options)
Initializes RPS-BLAST data structures.
#define USING_SCOPE(ns)
Use the specified namespace.
Definition: ncbistl.hpp:76
static BlastHSPPipe * CreateHspPipe(const CBlastOptionsMemento *opts_memento, BlastQueryInfo *query_info)
Create a pipe to be registered for use by stream.
Boolean Blast_QueryIsProtein(EBlastProgramType p)
Returns true if the query is protein.
Definition: blast_program.c:40
BlastDiagnostics * Blast_DiagnosticsFree(BlastDiagnostics *diagnostics)
Free the BlastDiagnostics structure and all substructures.
bool IsNull(void) const THROWS_NONE
Check if pointer is null – same effect as Empty().
Definition: ncbiobj.hpp:1305
Defines a concrete strategy for the IBlastSeqInfoSrc interface for sequence identifiers retrieval fro...
EBlastProgramType
Defines the engine's notion of the different applications of the BLAST algorithm. ...
Definition: blast_program.h:72
#define Blast_PerrorWithLocation(msg, error_code, context)
Convenient define to call the function Blast_PerrorEx.
int context
Context, allows us to print message for query number.
Definition: blast_message.h:75
bool Empty(void) const THROWS_NONE
Check if CConstRef is empty – not pointing to any object which means having a null value...
Definition: ncbiobj.hpp:1289
Collection of masked regions for a single query sequence.
Definition: seqlocinfo.hpp:112
Seq-aligns in the BLAST 2 Sequence style (one alignment per query-subject pair)
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.
bool Validate() const
Validate the options.
Structure to hold the a message from the core of the BLAST engine.
Definition: blast_message.h:70
Classes that capture the state of the BLAST options (or subsets of options) and restore them later (u...
Default implementation of BlastHSPStream.
const char * BlastSeqSrcGetName(const BlastSeqSrc *seq_src)
Get the Blast Sequence source name (e.g.
Definition: blast_seqsrc.c:235
In PHI BLAST, structure containing information about all pattern occurrences in query.
Definition: blast_def.h:300
static BlastHSPWriter * CreateHspWriter(const CBlastOptionsMemento *opts_memento, BLAST_SequenceBlk *query, BlastQueryInfo *query_info)
Create a writer to be registered for use by stream.
char * message
User message to be saved.
Definition: blast_message.h:73
EBlastProgramType GetProgramType() const
Returns the CORE BLAST notion of program type.
Wrapper structure for different types of BLAST lookup tables.
Definition: lookup_wrap.h:50
CRef –.
Definition: ncbiobj.hpp:616
static char * prog
Definition: mdb_load.c:33
void RemoveDuplicates()
Find and remove redundant messages.
Definition: blast_aux.cpp:1043
The context related information.
unsigned int
A callback function used to compare two keys in a database.
Definition: types.hpp:1153
User-defined methods of the data storage class.
typedef for the messages for an entire BLAST search, which could be comprised of multiple query seque...
CConstRef< objects::CSeq_loc > CreateWholeSeqLocFromIds(const list< CRef< objects::CSeq_id > > seqids)
Create a single CSeq_loc of type whole from the first id in the list.
CRef< CSBlastProgress > m_ProgressMonitor
The user data structure to aid in progress monitoring.
const CBlastOptionsMemento * CreateSnapshot() const
Create a snapshot of the state of this object for internal use of its data structures (BLAST C++ APIs...
Return type of BlastSetupPreliminarySearch.
Traceback stage.
Definition: blast_def.h:330
#define _ASSERT
int BlastHSPStreamRegisterMTLock(BlastHSPStream *hsp_stream, MT_LOCK lock)
Attach a mutex lock to a stream to protect multiple access during writing.
CRef< TBlastDiagnostics > m_Diagnostics
Diagnostic output from preliminary and traceback stages.
size_type size() const
Identical to GetNumResults, provided to facilitate STL-style iteration.
Definitions and functions associated with the BlastQueryInfo structure.
#define BEGIN_SCOPE(ns)
Define a new scope.
Definition: ncbistl.hpp:70
#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:546
Boolean Blast_ProgramIsRpsBlast(EBlastProgramType p)
Returns true if program is RPS-BLAST (i.e.
Definition: blast_program.c:73
CStructWrapper< BlastScoreBlk > TBlastScoreBlk
CRef< SBlastSetupData > BlastSetupPreliminarySearch(CRef< IQueryFactory > query_factory, CRef< CBlastOptions > options, size_t num_threads)
Set up internal data structures used by the BLAST CORE engine.
string BlastErrorCode2String(Int2 error_code)
Returns a string containing a human-readable interpretation of the error_code passed as this function...
#define const
Definition: zconf.h:217
CRef< ILocalQueryData > MakeLocalQueryData(const CBlastOptions *opts)
Creates and caches an ILocalQueryData.
Definition: query_data.cpp:52
TInterruptFnPtr m_FnInterrupt
The interrupt callback.
struct SPHIQueryInfo * pattern_info
Counts of PHI BLAST pattern occurrences, used in PHI BLAST only.
bool NotEmpty(void) const THROWS_NONE
Check if CRef is not empty – pointing to an object and has a non-null value.
Definition: ncbiobj.hpp:709
The query related information.
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.
BlastHSPStream * BlastHSPStreamFree(BlastHSPStream *hsp_stream)
Frees the BlastHSPStream structure by invoking the destructor function set by the user-defined constr...
const int kBlastMessageNoContext
Declared in blast_message.h as extern const.
Definition: blast_message.c:36
No variant selected.
Definition: Seq_loc_.hpp:97
C++ version of the initialization for the mutex locking interface.
const string & GetMsg(void) const
Get message string.
Definition: ncbiexpt.cpp:443
bool IsPacked_int(void) const
Check if variant Packed_int is selected.
Definition: Seq_loc_.hpp:534
#define BEGIN_NCBI_SCOPE
Define ncbi namespace.
Definition: ncbistl.hpp:98
void ClearDbIndexCallbacks(void)
CStructWrapper< BlastDiagnostics > TBlastDiagnostics
#define NCBI_RETHROW(prev_exception, exception_class, err_code, message)
Generic macro to re-throw an exception.
Definition: ncbiexpt.hpp:579
const unsigned int kQueryIndex
Index into multiple sequence alignment structure for the query sequence.
signed short Int2
Alias for signed short.
Definition: ncbitype.h:118
void Reset(void)
Reset reference object.
Definition: ncbiobj.hpp:756
bool IsQuerySplit() const
Determines whether the query sequence(s) are split or not.
Definition: split_query.hpp:88
TSeqLocInfoVector m_Masks
Modified on Tue Jul 25 19:56:51 2017 by modify_doxy.py rev. 533848