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

Go to the SVN repository for this file.

1 /* $Id: traceback_unit_test.cpp 97411 2022-07-14 19:13:48Z merezhuk $
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: Tom Madden
27 *
28 * File Description:
29 * Unit test module for traceback calculation
30 *
31 * ===========================================================================
32 */
33 #include <ncbi_pch.hpp>
34 #include <corelib/test_boost.hpp>
35 
36 #include <corelib/ncbitime.hpp>
38 #include <objmgr/scope.hpp>
40 #include <objmgr/util/sequence.hpp>
41 
46 #include <serial/serial.hpp>
47 #include <serial/iterator.hpp>
48 #include <serial/objostr.hpp>
49 
52 #include <blast_objmgr_priv.hpp>
53 
57 
70 
75 
78 
79 #include <common/test_data_path.h>
80 #include <corelib/ncbifile.hpp>
81 #include <algorithm>
82 
83 #include "test_objmgr.hpp"
84 #include "blast_test_util.hpp"
85 
86 #ifdef _OPENMP
87 #include <omp.h>
88 #endif
89 
90 using namespace std;
91 using namespace ncbi;
92 using namespace ncbi::objects;
93 using namespace ncbi::blast;
94 
96 
97 public:
98  BlastScoreBlk* m_ScoreBlk;
99  BlastSeqLoc* m_LookupSegments;
100  Blast_Message* m_BlastMessage;
101  BlastScoringParameters* m_ScoreParams;
102  BlastExtensionParameters* m_ExtParams;
103  BlastHitSavingParameters* m_HitParams;
104  BlastEffectiveLengthsParameters* m_EffLenParams;
105  BlastGapAlignStruct* m_GapAlign;
106 
107 
109  m_ScoreBlk=NULL;
110  m_LookupSegments=NULL;
111  m_BlastMessage=NULL;
112  m_ScoreParams = NULL;
113  m_ExtParams=NULL;
114  m_HitParams=NULL;
115  m_EffLenParams=NULL;
116  m_GapAlign=NULL;
117  }
118 
120  m_LookupSegments = BlastSeqLocFree(m_LookupSegments);
121  m_BlastMessage = Blast_MessageFree(m_BlastMessage);
122  m_ScoreBlk = BlastScoreBlkFree(m_ScoreBlk);
123  m_ScoreParams = BlastScoringParametersFree(m_ScoreParams);
124  m_HitParams = BlastHitSavingParametersFree(m_HitParams);
125  m_ExtParams = BlastExtensionParametersFree(m_ExtParams);
126  m_EffLenParams = BlastEffectiveLengthsParametersFree(m_EffLenParams);
127  m_GapAlign = BLAST_GapAlignStructFree(m_GapAlign);
128  }
129 
135 
136  BlastHSPWriter* writer = BlastHSPWriterNew(&writer_info, NULL, NULL);
137  BOOST_REQUIRE(writer_info == NULL);
138  return BlastHSPStreamNew(opt.GetProgramType(), opt.GetExtnOpts(),
139  FALSE, 1, writer);
140  }
141 
142  void x_SetupMain(const CBlastOptions &opt,
143  const CBLAST_SequenceBlk &query_blk,
144  const CBlastQueryInfo &query_info) {
146  opt.GetQueryOpts(),
147  opt.GetScoringOpts(),
148  query_blk, query_info, 1.0, &m_LookupSegments, NULL,
149  &m_ScoreBlk, &m_BlastMessage, &BlastFindMatrixPath);
150  }
151 
153  const BlastSeqSrc* seq_src,
154  const CBlastQueryInfo &query_info) {
155  BLAST_GapAlignSetUp(opt.GetProgramType(), seq_src,
156  opt.GetScoringOpts(),
157  opt.GetEffLenOpts(),
158  opt.GetExtnOpts(),
159  opt.GetHitSaveOpts(),
160  query_info, m_ScoreBlk, &m_ScoreParams,
161  &m_ExtParams, &m_HitParams, &m_EffLenParams,
162  &m_GapAlign);
163  }
164 
166  BlastHSPStream *hsp_stream,
167  const CBLAST_SequenceBlk &query_blk,
168  const CBlastQueryInfo &query_info,
169  const BlastSeqSrc* seq_src,
170  BlastHSPResults** results) {
172  hsp_stream, query_blk, query_info, seq_src,
173  m_GapAlign, m_ScoreParams, m_ExtParams, m_HitParams, m_EffLenParams,
174  opt.GetDbOpts(), NULL, NULL, NULL, results, 0, 0);
175  }
176 
177 };
178 
179 BOOST_FIXTURE_TEST_SUITE(traceback, CTracebackTestFixture)
180 
181 /* Checks that HSP data is updated correctly with traceback information. */
182 BOOST_AUTO_TEST_CASE(testHSPUpdateWithTraceback) {
183  const int kOffset=10;
184  BlastHSP* new_hsp = Blast_HSPNew();
185  BlastGapAlignStruct* gap_align =
187  gap_align->query_start = kOffset;
188  gap_align->query_stop = 2*kOffset;
189  gap_align->subject_start = 3*kOffset;
190  gap_align->subject_stop = 4*kOffset;
191  gap_align->score = 10*kOffset;
192  gap_align->edit_script = (GapEditScript*) calloc(1, sizeof(GapEditScript));
193 
194  Blast_HSPUpdateWithTraceback(gap_align, new_hsp);
195 
196  BOOST_REQUIRE_EQUAL((GapEditScript*) 0, gap_align->edit_script); // this was NULL'ed out
197  BOOST_REQUIRE(new_hsp->gap_info); // this got the pointer to edit_script
198  BOOST_REQUIRE_EQUAL(kOffset, new_hsp->query.offset);
199  BOOST_REQUIRE_EQUAL(2*kOffset, new_hsp->query.end);
200  BOOST_REQUIRE_EQUAL(3*kOffset, new_hsp->subject.offset);
201  BOOST_REQUIRE_EQUAL(4*kOffset, new_hsp->subject.end);
202  BOOST_REQUIRE_EQUAL(10*kOffset, new_hsp->score);
203 
204  new_hsp = Blast_HSPFree(new_hsp);
205  BOOST_REQUIRE(new_hsp == NULL);
206  gap_align = BLAST_GapAlignStructFree(gap_align);
207  BOOST_REQUIRE(gap_align == NULL);
208 }
209 
210 BOOST_AUTO_TEST_CASE(testBLASTNTraceBack) {
211  const int k_num_hsps_start = 9;
212  const int k_num_hsps_end = 7;
213 
214  CSeq_id qid("gi|1945388");
215  unique_ptr<SSeqLoc> qsl(
216  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_both));
217  CSeq_id sid("gi|1732684");
218  unique_ptr<SSeqLoc> ssl(
219  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_both));
220 
221  CBl2Seq blaster(*qsl, *ssl, eBlastn);
222 
223  CBlastQueryInfo query_info;
224  CBLAST_SequenceBlk query_blk;
225  TSearchMessages blast_msg;
226 
227  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
229  ENa_strand strand_opt = kOpts.GetStrandOption();
230 
231  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
232  prog, strand_opt, &query_info);
233  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
234  query_info, &query_blk, prog, strand_opt, blast_msg);
235  ITERATE(TSearchMessages, m, blast_msg) {
236  BOOST_REQUIRE(m->empty());
237  }
238 
239  BlastSeqSrc* seq_src =
241  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
244 
245  BlastHSPList* hsp_list =
246  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
247  hsp_list->oid = 0;
248  hsp_list->hspcnt = k_num_hsps_start;
249  hsp_list->allocated = k_num_hsps_start;
250  hsp_list->hsp_max = k_num_hsps_start;
251  hsp_list->do_not_reallocate = FALSE;
252  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
253 
254  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
255 
256  const int query_offset[k_num_hsps_start] = { 6020, 6022, 6622, 6622, 5295, 5199, 7191, 3818, 7408};
257  const int query_end[k_num_hsps_start] = { 6032, 6161, 6730, 6753, 5386, 5219, 7227, 3830, 7419};
258  const int subject_offset[k_num_hsps_start] = { 98, 104, 241, 241, 16, 0, 378, 71, 63};
259  const int subject_end[k_num_hsps_start] = { 110, 241, 350, 376, 107, 20, 415, 83, 74};
260  const int score[k_num_hsps_start] = { 17, 115, 93, 91, 91, 20, 17, 12, 11};
261  const int context[k_num_hsps_start] = { 0, 0, 0, 0, 0, 0, 0, 1, 1};
262  const int subject_frame[k_num_hsps_start] = { 1, 1, 1, 1, 1, 1, 1, 1, 1};
263  const int query_gapped_start[k_num_hsps_start] = { 20, 6035, 6625, 6745, 5295, 5199, 7193, 3819, 7409};
264  const int subject_gapped_start[k_num_hsps_start] = { 115, 116, 244, 368, 16, 0, 380, 72, 64};
265 
266  for (int index=0; index<k_num_hsps_start; index++)
267  {
268  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
269  hsp_list->hsp_array[index]->query.offset =query_offset[index];
270  hsp_list->hsp_array[index]->query.end =query_end[index];
271  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
272  hsp_list->hsp_array[index]->subject.end =subject_end[index];
273  hsp_list->hsp_array[index]->score =score[index];
274  hsp_list->hsp_array[index]->context =context[index];
275  hsp_list->hsp_array[index]->subject.frame =subject_frame[index];
276  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
277  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
278  }
279 
280  // needed after tie-breaking algorithm for scores was changed in
281  // ScoreCompareHSP (blast_hits.c, revision 1.139)
282  Blast_HSPListSortByScore(hsp_list);
283  BlastHSPStreamWrite(hsp_stream, &hsp_list);
284 
285  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
286 
287  // Set "database" length option to the length of subject sequence,
288  // to avoid having to calculate cutoffs and effective lengths twice.
289  Int4 oid = 0;
290  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
291  blaster.SetOptionsHandle().SetDbLength(subj_length);
292 
293  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
294 
295  BlastHSPResults* results = NULL;
296 
297  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
298  hsp_stream, query_blk, query_info, seq_src, &results);
299 
300  BlastHSPStreamFree(hsp_stream);
301 
302  const int query_offset_final[k_num_hsps_end] = { 6022, 6622, 5295, 7191, 5199, 7396, 3818};
303  const int query_end_final[k_num_hsps_end] = { 6161, 6759, 5386, 7231, 5219, 7425, 3830};
304  const int subject_offset_final[k_num_hsps_end] = { 104, 241, 16, 378, 0, 51, 71};
305  const int subject_end_final[k_num_hsps_end] = { 241, 383, 107, 419, 20, 80, 83};
306  const int score_final[k_num_hsps_end] = { 252, 226, 182, 54, 40, 27, 24};
307  const int context_final[k_num_hsps_end] = { 0, 0, 0, 0, 0, 1, 1};
308  const int subject_frame_final[k_num_hsps_end] = { 1, 1, 1, 1, 1, 1, 1};
309  const int query_gapped_start_final[k_num_hsps_end] = { 6035, 6625, 5295, 7193, 5199, 7409, 3819};
310  const int subject_gapped_start_final[k_num_hsps_end] = { 116, 244, 16, 380, 0, 64, 72};
311  const int num_ident_final[k_num_hsps_end] = { 135, 134, 91, 36, 20, 25, 12};
312 
313  // One hsp is dropped when the function runs.
314  BlastHitList* hit_list = results->hitlist_array[0];
315  hsp_list = hit_list->hsplist_array[0];
316 
317  BOOST_REQUIRE(hsp_list != NULL);
318  BOOST_REQUIRE_EQUAL(k_num_hsps_end, hsp_list->hspcnt);
319  for (int index=0; index<k_num_hsps_end; index++)
320  {
321  BlastHSP* tmp_hsp = hsp_list->hsp_array[index];
322  BOOST_REQUIRE_EQUAL(query_offset_final[index], tmp_hsp->query.offset);
323  BOOST_REQUIRE_EQUAL(query_end_final[index], tmp_hsp->query.end);
324  BOOST_REQUIRE_EQUAL(subject_offset_final[index], tmp_hsp->subject.offset);
325  BOOST_REQUIRE_EQUAL(subject_end_final[index], tmp_hsp->subject.end);
326  BOOST_REQUIRE_EQUAL(score_final[index], tmp_hsp->score);
327  BOOST_REQUIRE_EQUAL(context_final[index], (int) tmp_hsp->context);
328  BOOST_REQUIRE_EQUAL(subject_frame_final[index], (int) tmp_hsp->subject.frame);
329  BOOST_REQUIRE_EQUAL(query_gapped_start_final[index], tmp_hsp->query.gapped_start);
330  BOOST_REQUIRE_EQUAL(subject_gapped_start_final[index], tmp_hsp->subject.gapped_start);
331  BOOST_REQUIRE_EQUAL(num_ident_final[index], tmp_hsp->num_ident);
332  }
333 
334  results = Blast_HSPResultsFree(results);
335  BOOST_REQUIRE(results == NULL);
336  seq_src = BlastSeqSrcFree(seq_src);
337  BOOST_REQUIRE(seq_src == NULL);
338 }
339 
340 BOOST_AUTO_TEST_CASE(testBLASTNTraceBackLargeXDrop) {
341  const int k_num_hsps_start = 3;
342  const int k_num_hsps_end = 1;
343 
344  CSeq_id qid("gi|42254502");
345  unique_ptr<SSeqLoc> qsl(
346  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_both));
347  CSeq_id sid("gi|34787366");
348  unique_ptr<SSeqLoc> ssl(
349  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_both));
350 
351  CBlastNucleotideOptionsHandle opts_handle;
352  opts_handle.SetTraditionalBlastnDefaults();
353  opts_handle.SetMatchReward(2);
354  opts_handle.SetGapXDropoffFinal(200);
355 
356  CBl2Seq blaster(*qsl, *ssl, opts_handle);
357 
358  CBlastQueryInfo query_info;
359  CBLAST_SequenceBlk query_blk;
360  TSearchMessages blast_msg;
361 
362  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
364  ENa_strand strand_opt = kOpts.GetStrandOption();
365 
366  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
367  prog, strand_opt, &query_info);
368  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
369  query_info, &query_blk, prog, strand_opt, blast_msg);
370 
371  ITERATE(TSearchMessages, m, blast_msg) {
372  BOOST_REQUIRE(m->empty());
373  }
374 
375  BlastSeqSrc* seq_src =
377  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
378  opts_handle.GetOptions().GetProgramType());
380 
381  BlastHSPList* hsp_list =
382  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
383  hsp_list->oid = 0;
384  hsp_list->hspcnt = k_num_hsps_start;
385  hsp_list->allocated = k_num_hsps_start;
386  hsp_list->hsp_max = k_num_hsps_start;
387  hsp_list->do_not_reallocate = FALSE;
388  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
389 
390  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
391 
392  const int query_offset[k_num_hsps_start] = { 25194, 13986, 22457};
393  const int query_end[k_num_hsps_start] = { 31512, 17712, 25019};
394  const int subject_offset[k_num_hsps_start] = {11211, 0, 8471};
395  const int subject_end[k_num_hsps_start] = { 17529, 3726, 11036};
396  const int score[k_num_hsps_start] = { 12433, 7421, 4870};
397  const int context[k_num_hsps_start] = { 1, 1, 1};
398  const int subject_frame[k_num_hsps_start] = { 1, 1, 1};
399  const int query_gapped_start[k_num_hsps_start] = { 26671, 13986, 23372};
400  const int subject_gapped_start[k_num_hsps_start] = { 12688, 0, 9388};
401 
402  for (int index=0; index<k_num_hsps_start; index++)
403  {
404  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
405  hsp_list->hsp_array[index]->query.offset =query_offset[index];
406  hsp_list->hsp_array[index]->query.end =query_end[index];
407  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
408  hsp_list->hsp_array[index]->subject.end =subject_end[index];
409  hsp_list->hsp_array[index]->score =score[index];
410  hsp_list->hsp_array[index]->context =context[index];
411  hsp_list->hsp_array[index]->subject.frame =subject_frame[index];
412  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
413  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
414  }
415 
416  // needed after tie-breaking algorithm for scores was changed in
417  // ScoreCompareHSP (blast_hits.c, revision 1.139)
418  Blast_HSPListSortByScore(hsp_list);
419  BlastHSPStreamWrite(hsp_stream, &hsp_list);
420 
421  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
422 
423  // Set "database" length option to the length of subject sequence,
424  // to avoid having to calculate cutoffs and effective lengths twice.
425  Int4 oid = 0;
426  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
427  blaster.SetOptionsHandle().SetDbLength(subj_length);
428 
429  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
430 
431  BlastHSPResults* results = NULL;
432 
433  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
434  hsp_stream, query_blk, query_info, seq_src, &results);
435 
436  BlastHSPStreamFree(hsp_stream);
437 
438  const int query_offset_final[k_num_hsps_end] = { 13986};
439  const int query_end_final[k_num_hsps_end] = { 41877};
440  const int subject_offset_final[k_num_hsps_end] = { 0};
441  const int subject_end_final[k_num_hsps_end] = { 27888};
442  const int score_final[k_num_hsps_end] = { 55540};
443  const int context_final[k_num_hsps_end] = { 1};
444  const int subject_frame_final[k_num_hsps_end] = { 1};
445  const int query_gapped_start_final[k_num_hsps_end] = { 26671};
446  const int subject_gapped_start_final[k_num_hsps_end] = { 12688};
447  const int num_ident_final[k_num_hsps_end] = { 27856};
448 
449  // One hsp is dropped when the function runs.
450  BlastHitList* hit_list = results->hitlist_array[0];
451  hsp_list = hit_list->hsplist_array[0];
452 
453  BOOST_REQUIRE(hsp_list != NULL);
454  BOOST_REQUIRE_EQUAL(k_num_hsps_end, hsp_list->hspcnt);
455  for (int index=0; index<k_num_hsps_end; index++)
456  {
457  BlastHSP* tmp_hsp = hsp_list->hsp_array[index];
458  BOOST_REQUIRE_EQUAL(query_offset_final[index], tmp_hsp->query.offset);
459  BOOST_REQUIRE_EQUAL(query_end_final[index], tmp_hsp->query.end);
460  BOOST_REQUIRE_EQUAL(subject_offset_final[index], tmp_hsp->subject.offset);
461  BOOST_REQUIRE_EQUAL(subject_end_final[index], tmp_hsp->subject.end);
462  BOOST_REQUIRE_EQUAL(score_final[index], tmp_hsp->score);
463  BOOST_REQUIRE_EQUAL(context_final[index], (int) tmp_hsp->context);
464  BOOST_REQUIRE_EQUAL(subject_frame_final[index], (int) tmp_hsp->subject.frame);
465  BOOST_REQUIRE_EQUAL(query_gapped_start_final[index], tmp_hsp->query.gapped_start);
466  BOOST_REQUIRE_EQUAL(subject_gapped_start_final[index], tmp_hsp->subject.gapped_start);
467  BOOST_REQUIRE_EQUAL(num_ident_final[index], tmp_hsp->num_ident);
468  }
469 
470  results = Blast_HSPResultsFree(results);
471  BOOST_REQUIRE(results == NULL);
472  seq_src = BlastSeqSrcFree(seq_src);
473  BOOST_REQUIRE(seq_src == NULL);
474 }
475 
476 BOOST_AUTO_TEST_CASE(testBLASTPTraceBack) {
477  const int k_num_hsps_start = 12;
478  const int k_num_hsps_end = 10;
479 
480  CSeq_id qid("gi|42734333");
481  CSeq_id sid("gi|30176631");
482 
483  unique_ptr<SSeqLoc> qsl(
484  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_unknown));
485  unique_ptr<SSeqLoc> ssl(
486  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_unknown));
487 
488  CBlastProteinOptionsHandle opts_handle;
489  CBl2Seq blaster(*qsl, *ssl, opts_handle);
490 
491  CBlastQueryInfo query_info;
492  CBLAST_SequenceBlk query_blk;
493  TSearchMessages blast_msg;
494 
495  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
497  ENa_strand strand_opt = kOpts.GetStrandOption();
498 
499  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
500  prog, strand_opt, &query_info);
501  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
502  query_info, &query_blk, prog, strand_opt, blast_msg);
503  ITERATE(TSearchMessages, m, blast_msg) {
504  BOOST_REQUIRE(m->empty());
505  }
506 
507  BlastSeqSrc* seq_src =
509  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
510  opts_handle.GetOptions().GetProgramType());
512 
513  BlastHSPList* hsp_list =
514  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
515  hsp_list->oid = 0;
516  hsp_list->hspcnt = k_num_hsps_start;
517  hsp_list->allocated = k_num_hsps_start;
518  hsp_list->hsp_max = k_num_hsps_start;
519  hsp_list->do_not_reallocate = FALSE;
520  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
521 
522  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
523 
524  const int query_offset[k_num_hsps_start] = { 0, 3864, 3254, 1828, 2189, 795, 607, 1780, 1363, 2751, 3599, 242};
525  const int query_end[k_num_hsps_start] = { 307, 4287, 3556, 2058, 2269, 914, 741, 1821, 1451, 2810, 3631, 285};
526  const int subject_offset[k_num_hsps_start] = { 1, 2723, 2267, 1028, 1292, 634, 501, 925, 1195, 1795, 477, 1233};
527  const int subject_end[k_num_hsps_start] = { 321, 3171, 2537, 1243, 1371, 749, 618, 966, 1286, 1869, 509, 1276};
528  const int score[k_num_hsps_start] = { 370, 319, 139, 120, 89, 84, 75, 70, 69, 60, 47, 43};
529  const int query_gapped_start[k_num_hsps_start] = { 47, 4181, 3286, 2034, 2228, 871, 632, 1798, 1383, 2759, 3606, 259};
530  const int subject_gapped_start[k_num_hsps_start] = { 48, 3073, 2299, 1219, 1330, 709, 526, 943, 1215, 1803, 484, 1250};
531 
532  for (int index=0; index<k_num_hsps_start; index++)
533  {
534  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
535  hsp_list->hsp_array[index]->query.offset =query_offset[index];
536  hsp_list->hsp_array[index]->query.end =query_end[index];
537  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
538  hsp_list->hsp_array[index]->subject.end =subject_end[index];
539  hsp_list->hsp_array[index]->score =score[index];
540  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
541  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
542  }
543 
544  // needed after tie-breaking algorithm for scores was changed in
545  // ScoreCompareHSP (blast_hits.c, revision 1.139)
546  Blast_HSPListSortByScore(hsp_list);
547  BlastHSPStreamWrite(hsp_stream, &hsp_list);
548 
549  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
550 
551  // Set "database" length option to the length of subject sequence,
552  // to avoid having to calculate cutoffs and effective lengths twice.
553  Int4 oid = 0;
554  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
555  blaster.SetOptionsHandle().SetDbLength(subj_length);
556 
557  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
558 
559  BlastHSPResults* results = NULL;
560 
561  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
562  hsp_stream, query_blk, query_info, seq_src, &results);
563 
564  BlastHSPStreamFree(hsp_stream);
565 
566  const int query_offset_final[k_num_hsps_end] = { 0, 3864, 3254, 1780, 2189, 607, 1363, 2751, 3599, 242};
567  const int query_end_final[k_num_hsps_end] = { 307, 4287, 3556, 2058, 2599, 914, 1451, 2810, 3631, 285};
568  const int subject_offset_final[k_num_hsps_end] = { 1, 2723, 2267, 925, 1292, 501, 1195, 1795, 477, 1233};
569  const int subject_end_final[k_num_hsps_end] = { 321, 3171, 2537, 1243, 1704, 749, 1286, 1869, 509, 1276};
570  const int score_final[k_num_hsps_end] = { 367, 319, 139, 131, 122, 104, 69, 60, 47, 43};
571  const int query_gapped_start_final[k_num_hsps_end] = { 47, 4181, 3286, 2034, 2228, 871, 1383, 2759, 3606, 259};
572  const int subject_gapped_start_final[k_num_hsps_end] = { 48, 3073, 2299, 1219, 1330, 709, 1215, 1803, 484, 1250};
573  const int num_ident_final[k_num_hsps_end] = { 100, 122, 70, 61, 92, 54, 22, 18, 11, 9};
574 
575  BlastHitList* hit_list = results->hitlist_array[0];
576  hsp_list = hit_list->hsplist_array[0];
577 
578  // One hsp is dropped when the function runs.
579  BOOST_REQUIRE(hsp_list != NULL);
580  BOOST_REQUIRE_EQUAL(k_num_hsps_end, hsp_list->hspcnt);
581  for (int index=0; index<k_num_hsps_end; index++)
582  {
583  BlastHSP* tmp_hsp = hsp_list->hsp_array[index];
584  BOOST_REQUIRE_EQUAL(query_offset_final[index], tmp_hsp->query.offset);
585  BOOST_REQUIRE_EQUAL(query_end_final[index], tmp_hsp->query.end);
586  BOOST_REQUIRE_EQUAL(subject_offset_final[index], tmp_hsp->subject.offset);
587  BOOST_REQUIRE_EQUAL(subject_end_final[index], tmp_hsp->subject.end);
588  BOOST_REQUIRE_EQUAL(score_final[index], tmp_hsp->score);
589  BOOST_REQUIRE_EQUAL(query_gapped_start_final[index], tmp_hsp->query.gapped_start);
590  BOOST_REQUIRE_EQUAL(subject_gapped_start_final[index], tmp_hsp->subject.gapped_start);
591  BOOST_REQUIRE_EQUAL(num_ident_final[index], tmp_hsp->num_ident);
592  }
593 
594  results = Blast_HSPResultsFree(results);
595  BOOST_REQUIRE(results == NULL);
596  seq_src = BlastSeqSrcFree(seq_src);
597  BOOST_REQUIRE(seq_src == NULL);
598 }
599 
600 BOOST_AUTO_TEST_CASE(testTBLASTNTraceBack) {
601  const int k_num_hsps_start = 16;
602  const int k_num_hsps_end = 11;
603 
604  CSeq_id qid("gi|42734333");
605  CSeq_id sid("gi|27902043");
606 
607  unique_ptr<SSeqLoc> qsl(
608  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_unknown));
609  unique_ptr<SSeqLoc> ssl(
610  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_both));
611 
612  CTBlastnOptionsHandle opts_handle;
614  opts_handle.SetOptions().SetSegFiltering();
615 
616  CBl2Seq blaster(*qsl, *ssl, opts_handle);
617 
618  CBlastQueryInfo query_info;
619  CBLAST_SequenceBlk query_blk;
620  TSearchMessages blast_msg;
621 
622  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
624  ENa_strand strand_opt = kOpts.GetStrandOption();
625 
626  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
627  prog, strand_opt, &query_info);
628  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
629  query_info, &query_blk, prog, strand_opt, blast_msg);
630  ITERATE(TSearchMessages, m, blast_msg) {
631  BOOST_REQUIRE(m->empty());
632  }
633 
634  BlastSeqSrc* seq_src =
636  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
637  opts_handle.GetOptions().GetProgramType());
639 
640  BlastHSPList* hsp_list =
641  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
642  hsp_list->oid = 0;
643  hsp_list->hspcnt = k_num_hsps_start;
644  hsp_list->allocated = k_num_hsps_start;
645  hsp_list->hsp_max = k_num_hsps_start;
646  hsp_list->do_not_reallocate = FALSE;
647  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
648 
649  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
650 
651  const int query_offset[k_num_hsps_start] = { 149, 606, 656, 313, 221, 701, 57,
652  472, 0, 532, 404, 279, 125, 32, 371, 913};
653  const int query_end[k_num_hsps_start] = { 189, 662, 705, 377, 281, 747, 138,
654  533, 33, 575, 472, 312, 150, 59, 405, 941};
655  const int subject_offset[k_num_hsps_start] = { 58604, 59751, 59831, 58974,
656  58732, 59910, 58411, 59474, 58102, 59566, 59363, 58890, 58552, 58165, 59129, 9193};
657  const int subject_end[k_num_hsps_start] = { 58644, 59807, 59880, 59038, 58792,
658  59956, 58489, 59535, 58135, 59609, 59432, 58923, 58577, 58192, 59163, 9221};
659  const int score[k_num_hsps_start] = { 253, 237, 214, 193, 183, 178, 168, 165, 162,
660  149, 125, 120, 120, 113, 100, 55};
661  const int subject_frame[k_num_hsps_start] = { -3, -2, -2, -2, -3, -2, -1, -3, -2,
662  -1, -2, -2, -3, -3, -3, 3};
663  const int query_gapped_start[k_num_hsps_start] = { 173, 611, 662, 319, 254, 719, 72,
664  491, 11, 554, 438, 286, 131, 39, 379, 929};
665  const int subject_gapped_start[k_num_hsps_start] = { 58629, 59756, 59837, 58980, 58765,
666  59928, 58426, 59493, 58113, 59588, 59399, 58897, 58558, 58172, 59137, 9209};
667 
668  for (int index=0; index<k_num_hsps_start; index++)
669  {
670  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
671  hsp_list->hsp_array[index]->query.offset =query_offset[index];
672  hsp_list->hsp_array[index]->query.end =query_end[index];
673  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
674  hsp_list->hsp_array[index]->subject.end =subject_end[index];
675  hsp_list->hsp_array[index]->score =score[index];
676  hsp_list->hsp_array[index]->subject.frame =subject_frame[index];
677  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
678  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
679  }
680 
681  // needed after tie-breaking algorithm for scores was changed in
682  // ScoreCompareHSP (blast_hits.c, revision 1.139)
683  Blast_HSPListSortByScore(hsp_list);
684  BlastHSPStreamWrite(hsp_stream, &hsp_list);
685 
686  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
687 
688  // Set "database" length option to the length of subject sequence,
689  // to avoid having to calculate cutoffs and effective lengths twice.
690  Int4 oid = 0;
691  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
692  blaster.SetOptionsHandle().SetDbLength(subj_length);
693 
694  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
695 
696  BlastHSPResults* results = NULL;
697 
698  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
699  hsp_stream, query_blk, query_info, seq_src, &results);
700 
701  BlastHSPStreamFree(hsp_stream);
702 
703  const int query_offset_final[k_num_hsps_end] = {606, 125, 279, 57, 472, 0, 532,
704  428, 32, 371, 913};
705  const int query_end_final[k_num_hsps_end] = {747, 281, 377, 138, 533, 33, 575,
706  472, 59, 405, 941};
707  const int subject_offset_final[k_num_hsps_end] = {59751, 58552, 58890, 58411, 59474,
708  58102, 59566, 59389, 58165, 59129, 9193};
709  const int subject_end_final[k_num_hsps_end] = {59956, 58792, 59038, 58489, 59535,
710  58135, 59609, 59432, 58192, 59163, 9221};
711  const int score_final[k_num_hsps_end] = {525, 465, 250, 167, 165, 162, 149,
712  123, 113, 100, 55};
713  const int subject_frame_final[k_num_hsps_end] = {-2, -3, -2, -1, -3, -2,
714  -1, -2, -3, -3, 3};
715  const int query_gapped_start_final[k_num_hsps_end] = {611, 173, 319, 72, 491,
716  11, 554, 438, 39, 379, 929};
717  const int subject_gapped_start_final[k_num_hsps_end] = {59756, 58629,
718  58980, 58426, 59493, 58113, 59588, 59399, 58172, 59137, 9209};
719  const int num_ident_final[k_num_hsps_end] = {116, 105, 54, 44, 27, 31, 29,
720  25, 22, 21, 12};
721  const int nums[k_num_hsps_end] = {1, 1, 1, 1, 2, 2, 2, 1, 2, 1, 1};
722 
723  BlastHitList* hit_list = results->hitlist_array[0];
724  hsp_list = hit_list->hsplist_array[0];
725 
726  BOOST_REQUIRE(hsp_list != NULL);
727  BOOST_REQUIRE_EQUAL(k_num_hsps_end, hsp_list->hspcnt);
728  for (int index=0; index<k_num_hsps_end; index++)
729  {
730  BlastHSP* tmp_hsp = hsp_list->hsp_array[index];
731 
732  BOOST_REQUIRE_EQUAL(query_offset_final[index], tmp_hsp->query.offset);
733  BOOST_REQUIRE_EQUAL(query_end_final[index], tmp_hsp->query.end);
734  BOOST_REQUIRE_EQUAL(subject_offset_final[index], tmp_hsp->subject.offset);
735  BOOST_REQUIRE_EQUAL(subject_end_final[index], tmp_hsp->subject.end);
736  BOOST_REQUIRE_EQUAL(score_final[index], tmp_hsp->score);
737  BOOST_REQUIRE_EQUAL(subject_frame_final[index], (int) tmp_hsp->subject.frame);
738  BOOST_REQUIRE_EQUAL(query_gapped_start_final[index], tmp_hsp->query.gapped_start);
739  BOOST_REQUIRE_EQUAL(subject_gapped_start_final[index], tmp_hsp->subject.gapped_start);
740  BOOST_REQUIRE_EQUAL(num_ident_final[index], tmp_hsp->num_ident);
741  BOOST_REQUIRE_EQUAL(nums[index], tmp_hsp->num);
742  }
743 
744  results = Blast_HSPResultsFree(results);
745  BOOST_REQUIRE(results == NULL);
746  seq_src = BlastSeqSrcFree(seq_src);
747  BOOST_REQUIRE(seq_src == NULL);
748 }
749 
750 BOOST_AUTO_TEST_CASE(testNoHSPEvalueCutoffBeforeLink) {
751  const EProgram kProgram = eTblastn;
752  const EBlastProgramType kProgramType = eBlastTypeTblastn;
753  const int kNumHsps = 3;
754  const int q_offsets[kNumHsps] = { 1, 144, 203 };
755  const int q_ends[kNumHsps] = { 151, 191, 226 };
756  const int q_gapped_starts[kNumHsps] = { 23, 153, 209 };
757  const int s_offsets[kNumHsps] = { 501, 655, 736 };
758  const int s_ends[kNumHsps] = { 648, 702, 756 };
759  const int s_gapped_starts[kNumHsps] = { 523, 664, 742 };
760  const int s_frames[kNumHsps] = { 3, 1, 3 };
761  const int scores[kNumHsps] = { 211, 91, 52 };
762  const Int8 kSearchSp = 20763230804LL;
763  const string kDbName("data/nt.41646578");
764 
765  CRef<CSeq_id> qid(new CSeq_id("gi|129295"));
766 
769 
771  query.AddQuery(Q1);
772 
773  CRef<CSeqDB> seqdb(new CSeqDB(kDbName, CSeqDB::eNucleotide));
774  CBlastSeqSrc seq_src = SeqDbBlastSeqSrcInit(seqdb);
775  CRef<blast::IBlastSeqInfoSrc> seq_info_src(new blast::CSeqDbSeqInfoSrc(seqdb));
776 
777  BlastHSPList* hsp_list = Blast_HSPListNew(kNumHsps);
778  for (int index = 0; index < kNumHsps; ++index) {
779  BlastHSP* hsp = hsp_list->hsp_array[index] = Blast_HSPNew();
780  hsp->score = scores[index];
781  hsp->query.offset = q_offsets[index];
782  hsp->query.end = q_ends[index];
783  hsp->query.gapped_start = q_gapped_starts[index];
784  hsp->subject.offset = s_offsets[index];
785  hsp->subject.end = s_ends[index];
786  hsp->subject.gapped_start = s_gapped_starts[index];
787  hsp->subject.frame = s_frames[index];
788  }
789  hsp_list->hspcnt = kNumHsps;
790 
791  BlastExtensionOptions* ext_options = NULL;
792  BlastExtensionOptionsNew(kProgramType, &ext_options, true);
793 
794  BlastScoringOptions* scoring_options = NULL;
795  BlastScoringOptionsNew(kProgramType, &scoring_options);
796 
797  BlastHitSavingOptions* hit_options = NULL;
798  BlastHitSavingOptionsNew(kProgramType, &hit_options,
799  scoring_options->gapped_calculation);
800 
803  hit_options, ext_options->compositionBasedStats,
804  scoring_options->gapped_calculation));
805 
806  BlastHSPWriter* writer = BlastHSPWriterNew(&writer_info, NULL, NULL);
807  BOOST_REQUIRE(writer_info == NULL);
808  BlastHSPStream* hsp_stream = BlastHSPStreamNew(
809  kProgramType, ext_options, FALSE, 1, writer);
810 
811  hit_options = BlastHitSavingOptionsFree(hit_options);
812  BOOST_REQUIRE(hit_options == NULL);
813  scoring_options = BlastScoringOptionsFree(scoring_options);
814  BOOST_REQUIRE(scoring_options == NULL);
815  ext_options = BlastExtensionOptionsFree(ext_options);
816  BOOST_REQUIRE(ext_options == NULL);
817  // needed after tie-breaking algorithm for scores was changed in
818  // ScoreCompareHSP (blast_hits.c, revision 1.139)
819  Blast_HSPListSortByScore(hsp_list);
820  BlastHSPStreamWrite(hsp_stream, &hsp_list);
821 
822  // Run traceback on this HSP list, without producing a Seq-align.
823 
825 
827  cboh(CBlastOptionsFactory::Create(kProgram));
828 
830  (const_cast<CBlastOptions*>(& cboh->GetOptions()));
831 
832  cboh->SetEffectiveSearchSpace(kSearchSp);
833 
835  (WrapStruct(hsp_stream, BlastHSPStreamFree));
836 
837  CBlastTracebackSearch search(qf, opts, seq_src.Get(), seq_info_src, hsps);
838 
839  CSearchResultSet crs = *search.Run();
840 
841  BOOST_REQUIRE_EQUAL((int) crs.GetNumResults(), 1);
842 
843  BOOST_REQUIRE_EQUAL(kNumHsps, (int)crs[0].GetSeqAlign()->Size());
844 }
845 
846 BOOST_AUTO_TEST_CASE(testFilterBlastResults_QueryCov) {
847  const int k_num_hsps_start = 12;
848  const int k_num_hsps_filtered = 6;
849 
850  CSeq_id qid("gi|42734333");
851  CSeq_id sid("gi|30176631");
852 
853  unique_ptr<SSeqLoc> qsl(
854  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_unknown));
855  unique_ptr<SSeqLoc> ssl(
856  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_unknown));
857 
858  CBlastProteinOptionsHandle opts_handle;
859  CBl2Seq blaster(*qsl, *ssl, opts_handle);
860 
861  CBlastQueryInfo query_info;
862  CBLAST_SequenceBlk query_blk;
863  TSearchMessages blast_msg;
864 
865  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
867  ENa_strand strand_opt = kOpts.GetStrandOption();
868 
869  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
870  prog, strand_opt, &query_info);
871  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
872  query_info, &query_blk, prog, strand_opt, blast_msg);
873  ITERATE(TSearchMessages, m, blast_msg) {
874  BOOST_REQUIRE(m->empty());
875  }
876 
877  BlastSeqSrc* seq_src =
879  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
880  opts_handle.GetOptions().GetProgramType());
882 
883  BlastHSPList* hsp_list =
884  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
885  hsp_list->oid = 0;
886  hsp_list->hspcnt = k_num_hsps_start;
887  hsp_list->allocated = k_num_hsps_start;
888  hsp_list->hsp_max = k_num_hsps_start;
889  hsp_list->do_not_reallocate = FALSE;
890  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
891 
892  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
893 
894  const int query_offset[k_num_hsps_start] = { 0, 3864, 3254, 1828, 2189, 795, 607, 1780, 1363, 2751, 3599, 242};
895  const int query_end[k_num_hsps_start] = { 307, 4287, 3556, 2058, 2269, 914, 741, 1821, 1451, 2810, 3631, 285};
896  const int subject_offset[k_num_hsps_start] = { 1, 2723, 2267, 1028, 1292, 634, 501, 925, 1195, 1795, 477, 1233};
897  const int subject_end[k_num_hsps_start] = { 321, 3171, 2537, 1243, 1371, 749, 618, 966, 1286, 1869, 509, 1276};
898  const int score[k_num_hsps_start] = { 370, 319, 139, 120, 89, 84, 75, 70, 69, 60, 47, 43};
899  const int query_gapped_start[k_num_hsps_start] = { 47, 4181, 3286, 2034, 2228, 871, 632, 1798, 1383, 2759, 3606, 259};
900  const int subject_gapped_start[k_num_hsps_start] = { 48, 3073, 2299, 1219, 1330, 709, 526, 943, 1215, 1803, 484, 1250};
901 
902  for (int index=0; index<k_num_hsps_start; index++)
903  {
904  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
905  hsp_list->hsp_array[index]->query.offset =query_offset[index];
906  hsp_list->hsp_array[index]->query.end =query_end[index];
907  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
908  hsp_list->hsp_array[index]->subject.end =subject_end[index];
909  hsp_list->hsp_array[index]->score =score[index];
910  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
911  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
912  }
913 
914  // needed after tie-breaking algorithm for scores was changed in
915  // ScoreCompareHSP (blast_hits.c, revision 1.139)
916  Blast_HSPListSortByScore(hsp_list);
917  BlastHSPStreamWrite(hsp_stream, &hsp_list);
918 
919  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
920 
921  // Set "database" length option to the length of subject sequence,
922  // to avoid having to calculate cutoffs and effective lengths twice.
923  Int4 oid = 0;
924  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
925  blaster.SetOptionsHandle().SetDbLength(subj_length);
927 
928  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
929 
930  BlastHSPResults* results = NULL;
931 
932  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
933  hsp_stream, query_blk, query_info, seq_src, &results);
934 
935  BlastHSPStreamFree(hsp_stream);
936 
937  BlastHitList* hit_list = results->hitlist_array[0];
938  hsp_list = hit_list->hsplist_array[0];
939 
940  // One hsp is dropped when the function runs.
941  BOOST_REQUIRE(hsp_list != NULL);
942  BOOST_REQUIRE_EQUAL(k_num_hsps_filtered, hsp_list->hspcnt);
943 
944  results = Blast_HSPResultsFree(results);
945  BOOST_REQUIRE(results == NULL);
946  seq_src = BlastSeqSrcFree(seq_src);
947  BOOST_REQUIRE(seq_src == NULL);
948 
949 }
950 
951 BOOST_AUTO_TEST_CASE(testFilterBlastResults_MaxHsps) {
952  const int k_num_hsps_start = 9;
953  const int k_num_hsps_filtered = 4;
954 
955  CSeq_id qid("gi|1945388");
956  unique_ptr<SSeqLoc> qsl(
957  CTestObjMgr::Instance().CreateSSeqLoc(qid, eNa_strand_both));
958  CSeq_id sid("gi|1732684");
959  unique_ptr<SSeqLoc> ssl(
960  CTestObjMgr::Instance().CreateSSeqLoc(sid, eNa_strand_both));
961 
962  CBl2Seq blaster(*qsl, *ssl, eBlastn);
963 
964  CBlastQueryInfo query_info;
965  CBLAST_SequenceBlk query_blk;
966  TSearchMessages blast_msg;
967 
968  const CBlastOptions& kOpts = blaster.GetOptionsHandle().GetOptions();
970  ENa_strand strand_opt = kOpts.GetStrandOption();
971 
972  SetupQueryInfo(const_cast<TSeqLocVector&>(blaster.GetQueries()),
973  prog, strand_opt, &query_info);
974  SetupQueries(const_cast<TSeqLocVector&>(blaster.GetQueries()),
975  query_info, &query_blk, prog, strand_opt, blast_msg);
976  ITERATE(TSearchMessages, m, blast_msg) {
977  BOOST_REQUIRE(m->empty());
978  }
979 
980  BlastSeqSrc* seq_src =
982  const_cast<TSeqLocVector&>(blaster.GetSubjects()),
985 
986  BlastHSPList* hsp_list =
987  (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
988  hsp_list->oid = 0;
989  hsp_list->hspcnt = k_num_hsps_start;
990  hsp_list->allocated = k_num_hsps_start;
991  hsp_list->hsp_max = k_num_hsps_start;
992  hsp_list->do_not_reallocate = FALSE;
993  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
994 
995  BlastHSPStream* hsp_stream = x_MakeStream(blaster.GetOptionsHandle().GetOptions());
996 
997  const int query_offset[k_num_hsps_start] = { 6020, 6022, 6622, 6622, 5295, 5199, 7191, 3818, 7408};
998  const int query_end[k_num_hsps_start] = { 6032, 6161, 6730, 6753, 5386, 5219, 7227, 3830, 7419};
999  const int subject_offset[k_num_hsps_start] = { 98, 104, 241, 241, 16, 0, 378, 71, 63};
1000  const int subject_end[k_num_hsps_start] = { 110, 241, 350, 376, 107, 20, 415, 83, 74};
1001  const int score[k_num_hsps_start] = { 17, 115, 93, 91, 91, 20, 17, 12, 11};
1002  const int context[k_num_hsps_start] = { 0, 0, 0, 0, 0, 0, 0, 1, 1};
1003  const int subject_frame[k_num_hsps_start] = { 1, 1, 1, 1, 1, 1, 1, 1, 1};
1004  const int query_gapped_start[k_num_hsps_start] = { 20, 6035, 6625, 6745, 5295, 5199, 7193, 3819, 7409};
1005  const int subject_gapped_start[k_num_hsps_start] = { 115, 116, 244, 368, 16, 0, 380, 72, 64};
1006 
1007  for (int index=0; index<k_num_hsps_start; index++)
1008  {
1009  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
1010  hsp_list->hsp_array[index]->query.offset =query_offset[index];
1011  hsp_list->hsp_array[index]->query.end =query_end[index];
1012  hsp_list->hsp_array[index]->subject.offset =subject_offset[index];
1013  hsp_list->hsp_array[index]->subject.end =subject_end[index];
1014  hsp_list->hsp_array[index]->score =score[index];
1015  hsp_list->hsp_array[index]->context =context[index];
1016  hsp_list->hsp_array[index]->subject.frame =subject_frame[index];
1017  hsp_list->hsp_array[index]->query.gapped_start =query_gapped_start[index];
1018  hsp_list->hsp_array[index]->subject.gapped_start =subject_gapped_start[index];
1019  }
1020 
1021  // needed after tie-breaking algorithm for scores was changed in
1022  // ScoreCompareHSP (blast_hits.c, revision 1.139)
1023  Blast_HSPListSortByScore(hsp_list);
1024  BlastHSPStreamWrite(hsp_stream, &hsp_list);
1025 
1026  x_SetupMain(blaster.GetOptionsHandle().GetOptions(), query_blk, query_info);
1027 
1028  // Set "database" length option to the length of subject sequence,
1029  // to avoid having to calculate cutoffs and effective lengths twice.
1030  Int4 oid = 0;
1031  Uint4 subj_length = BlastSeqSrcGetSeqLen(seq_src, (void*)&oid);
1032  blaster.SetOptionsHandle().SetDbLength(subj_length);
1033  blaster.SetOptionsHandle().SetMaxHspsPerSubject(k_num_hsps_filtered);
1034 
1035  x_SetupGapAlign(blaster.GetOptionsHandle().GetOptions(), seq_src, query_info);
1036 
1037  BlastHSPResults* results = NULL;
1038 
1039  x_ComputeTracebak(blaster.GetOptionsHandle().GetOptions(),
1040  hsp_stream, query_blk, query_info, seq_src, &results);
1041 
1042  BlastHSPStreamFree(hsp_stream);
1043 
1044  // One hsp is dropped when the function runs.
1045  BlastHitList* hit_list = results->hitlist_array[0];
1046  hsp_list = hit_list->hsplist_array[0];
1047 
1048  BOOST_REQUIRE(hsp_list != NULL);
1049  BOOST_REQUIRE_EQUAL(k_num_hsps_filtered, hsp_list->hspcnt);
1050 
1051  results = Blast_HSPResultsFree(results);
1052  BOOST_REQUIRE(results == NULL);
1053  seq_src = BlastSeqSrcFree(seq_src);
1054  BOOST_REQUIRE(seq_src == NULL);
1055 }
1056 
1057 BlastHSPList * s_GetHSPList(int num_hsps, int oid, int len, int range_len)
1058 {
1059  BlastHSPList* hsp_list = (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
1060  hsp_list->oid = oid;
1061  hsp_list->hspcnt = num_hsps;
1062  hsp_list->allocated = num_hsps;
1063  hsp_list->hsp_max = num_hsps;
1064  hsp_list->do_not_reallocate = FALSE;
1065  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
1066 
1067  vector<int> subject_offset;
1068  for (int i =0; i < num_hsps; i ++) {
1069  subject_offset.push_back(rand()%(len -range_len));
1070  }
1071 
1072  for (int index=0; index < num_hsps; index++) {
1073  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
1074  hsp_list->hsp_array[index]->query.offset = 0;
1075  hsp_list->hsp_array[index]->query.end = 10000;
1076  hsp_list->hsp_array[index]->subject.offset = subject_offset[index];
1077  hsp_list->hsp_array[index]->subject.end = subject_offset[index] + range_len;
1078  hsp_list->hsp_array[index]->score = 500;
1079  hsp_list->hsp_array[index]->context = 0;
1080  hsp_list->hsp_array[index]->subject.frame = 1;
1081  hsp_list->hsp_array[index]->query.gapped_start = 500;
1082  hsp_list->hsp_array[index]->subject.gapped_start =subject_offset[index] + 50;
1083  }
1084 
1085  return hsp_list;
1086 }
1087 
1088 BOOST_AUTO_TEST_CASE(testUnit_TestsDataPath) {
1089  string rel_path = "blast/algo/unit_tests/api/data/long_seqs";
1090  std::replace( rel_path.begin(), rel_path.end(),'/', CDirEntry::GetPathSeparator() ); // SB-3398 fix sub directory path to be platform compatible
1091 
1092  BlastSeqSrc* seqSrc = SeqDbBlastSeqSrcInit( CFile::ConcatPath( NCBI_GetTestDataPath(), rel_path) , false );
1093 
1094  if( seqSrc && BlastSeqSrcGetInitError(seqSrc) ) {
1095  cout << "testUnit_TestsDataPath: BlastSeqSrcGetInitError: "<<BlastSeqSrcGetInitError(seqSrc) << endl;
1096  cout << "testUnit_TestsDataPath: Can't access: "<< ( NCBI_GetTestDataPath() + rel_path ) << endl;
1097  }
1098  BOOST_REQUIRE(seqSrc != NULL);
1099  BOOST_REQUIRE(BlastSeqSrcGetInitError(seqSrc) == NULL);
1100 
1101 }
1102 
1103 BOOST_AUTO_TEST_CASE(testSetupPartialFetching) {
1104  string rel_path = "blast/algo/unit_tests/api/data/long_seqs";
1105  std::replace( rel_path.begin(), rel_path.end(),'/', CDirEntry::GetPathSeparator() ); // SB-3398 fix sub directory path to be platform compatible
1106 
1107  BlastSeqSrc* seqSrc = SeqDbBlastSeqSrcInit( CFile::ConcatPath( NCBI_GetTestDataPath(), rel_path) , false );
1108 
1109  BOOST_REQUIRE(seqSrc != NULL);
1110  BOOST_REQUIRE(BlastSeqSrcGetInitError(seqSrc) == NULL);
1111 
1112  const int k_num_hsps = 6;
1113  BlastHSPList* hsp_list = (BlastHSPList*) calloc(1, sizeof(BlastHSPList));
1114  hsp_list->oid = 4;
1115  hsp_list->hspcnt = k_num_hsps;
1116  hsp_list->allocated = k_num_hsps;
1117  hsp_list->hsp_max = k_num_hsps;
1118  hsp_list->do_not_reallocate = FALSE;
1119  hsp_list->hsp_array = (BlastHSP**) malloc(hsp_list->allocated*sizeof(BlastHSP*));
1120 
1121  const int subject_offset[k_num_hsps] = { 0, 1110000, 98765, 1430868, 2134657, 25285936};
1122 
1123  for (int index=0; index<k_num_hsps; index++) {
1124  hsp_list->hsp_array[index] = (BlastHSP*) calloc(1, sizeof(BlastHSP));
1125  hsp_list->hsp_array[index]->query.offset = 0;
1126  hsp_list->hsp_array[index]->query.end = 10000;
1127  hsp_list->hsp_array[index]->subject.offset = subject_offset[index];
1128  hsp_list->hsp_array[index]->subject.end = subject_offset[index] + 1000;
1129  hsp_list->hsp_array[index]->score = 500;
1130  hsp_list->hsp_array[index]->context = 0;
1131  hsp_list->hsp_array[index]->subject.frame = 1;
1132  hsp_list->hsp_array[index]->query.gapped_start = 500;
1133  hsp_list->hsp_array[index]->subject.gapped_start =subject_offset[index] + 50;
1134  }
1135 
1136 
1138  eBlastTypeBlastn, seqSrc, (const BlastHSPList**) &hsp_list, 1);
1139  Int4 len = BlastSeqSrcGetSeqLen(seqSrc, &(hsp_list->oid));
1140 
1141  BOOST_REQUIRE_EQUAL(ranges->num_ranges, 6);
1142  BOOST_REQUIRE_EQUAL(ranges->ranges[0], 0);
1143  BOOST_REQUIRE_EQUAL(ranges->ranges[11], len);
1144  free(hsp_list);
1145 }
1146 
1147 BOOST_AUTO_TEST_CASE(testPartialFetchingMT) {
1148  string rel_path = "blast/algo/unit_tests/api/data/long_seqs";
1149  std::replace( rel_path.begin(), rel_path.end(),'/', CDirEntry::GetPathSeparator() ); // SB-3398 fix sub directory path to be platform compatible
1150 
1151  BlastSeqSrc* seqSrc = SeqDbBlastSeqSrcInit( CFile::ConcatPath( NCBI_GetTestDataPath(), rel_path) , false );
1152 
1153  BOOST_REQUIRE(seqSrc != NULL);
1154  BOOST_REQUIRE(BlastSeqSrcGetInitError(seqSrc) == NULL);
1155 
1156  const int kMaxNum = 36;
1157  const int kRangeLen = 10000;
1158  int kOid = 4;
1159  Int4 len = BlastSeqSrcGetSeqLen(seqSrc, &kOid);
1160 
1161  BlastSeqSrcGetSeqArg seq_arg_ref;
1162  seq_arg_ref.oid = kOid;
1164  seq_arg_ref.reset_ranges = false;
1165  seq_arg_ref.check_oid_exclusion = true;
1166  seq_arg_ref.seq = NULL;
1167  seq_arg_ref.ranges = NULL;
1168  BlastSeqSrcGetSequence(seqSrc, &seq_arg_ref);
1169 
1170 
1171  BlastSeqSrcGetSeqArg** seq_arg_array = (BlastSeqSrcGetSeqArg**) calloc (kMaxNum, sizeof(BlastSeqSrcGetSeqArg*));
1172  BlastHSPList ** hsp_list_array = (BlastHSPList **) calloc (kMaxNum, sizeof(BlastHSPList *));
1173  int i;
1174 #pragma omp parallel for default(none) num_threads(8) schedule(guided) shared(seqSrc, len, kOid, seq_arg_array, hsp_list_array)
1175  for (i = 0; i < kMaxNum; i++) {
1176  BlastHSPList * hsp_list = s_GetHSPList(i + 2, kOid, len, kRangeLen);
1177  hsp_list_array[i] = hsp_list;
1178  BlastSeqSrc* seq_src = BlastSeqSrcCopy(seqSrc);
1180  seq_arg->oid = hsp_list->oid;
1182  seq_arg->reset_ranges = false;
1183  seq_arg->check_oid_exclusion = true;
1184  seq_arg->seq = NULL;
1185  seq_arg->ranges = BLAST_SetupPartialFetching(eBlastTypeBlastn, seq_src, (const BlastHSPList**) &hsp_list, 1);
1186  seq_arg_array[i] = seq_arg;
1187 
1188  BlastSeqSrcGetSequence(seq_src, seq_arg);
1189  seq_arg->ranges = BlastSeqSrcSetRangesArgFree(seq_arg->ranges);
1190  seq_src = BlastSeqSrcFree(seq_src);
1191  }
1192 
1193 
1194  for(int k=0; k < kMaxNum; k++) {
1195  BlastHSPList * hsp_list = hsp_list_array[k];
1196  BlastSeqSrcGetSeqArg* seq_arg = seq_arg_array[k];
1197  for (int j =0; j < hsp_list->hspcnt; j++){
1198  int b = hsp_list->hsp_array[j]->subject.offset;
1199  int e = hsp_list->hsp_array[j]->subject.end;
1200  int r = e-b;
1201  BOOST_REQUIRE(memcmp((seq_arg->seq->sequence_start + b), (seq_arg_ref.seq->sequence_start + b), r) == 0);
1202  }
1203  seq_arg_array[k]->seq = BlastSequenceBlkFree(seq_arg_array[k]->seq);
1204  free(seq_arg_array[k]);
1205  free (hsp_list_array[k]);
1206  }
1207  free(seq_arg_array);
1208  free(hsp_list_array);
1209 }
1210 
1212 
1213 /*
1214 * ===========================================================================
1215 *
1216 * $Log: traceback-cppunit.cpp,v $
1217 * Revision 1.76 2009/03/13 19:25:56 maning
1218 * Previous commit messed up. Roll back again.
1219 *
1220 * Revision 1.74 2009/03/13 18:45:15 maning
1221 * Roll back to previous version.
1222 *
1223 * Revision 1.72 2008/12/16 21:24:10 madden
1224 * Correct offsets
1225 *
1226 * Revision 1.71 2008/10/31 22:05:27 avagyanv
1227 * Set old default values in options to avoid diffs in unit tests
1228 *
1229 * Revision 1.70 2008/10/27 17:00:13 camacho
1230 * Fix include paths to deprecated headers
1231 *
1232 * Revision 1.69 2008/04/15 13:50:29 madden
1233 * Update tests for svn 124499
1234 *
1235 * Revision 1.68 2008/02/13 21:39:12 camacho
1236 * Re-enable choice to sort by score to meet pre-condition of composition-based
1237 * statistics code.
1238 *
1239 * Revision 1.67 2007/10/22 19:16:10 madden
1240 * BlastExtensionOptionsNew has Boolean gapped arg
1241 *
1242 * Revision 1.66 2007/07/27 18:04:34 papadopo
1243 * change signature of HSPListCollector_Init
1244 *
1245 * Revision 1.65 2007/07/25 12:41:39 madden
1246 * Accomodates changes to blastn type defaults
1247 *
1248 * Revision 1.64 2007/04/16 19:38:26 camacho
1249 * Update code following CSearchResultSet changes
1250 *
1251 * Revision 1.63 2007/04/10 18:24:36 madden
1252 * Remove discontinuous seq-aligns
1253 *
1254 * Revision 1.62 2007/03/22 14:34:44 camacho
1255 * + support for auto-detection of genetic codes
1256 *
1257 * Revision 1.61 2007/03/20 14:54:02 camacho
1258 * changes related to addition of multiple genetic code specification
1259 *
1260 * Revision 1.60 2006/10/16 19:33:22 madden
1261 * Call in BlastHitSavingOptionsFree testNoHSPEvalueCutoffBeforeLink
1262 *
1263 * Revision 1.59 2006/10/10 19:47:01 bealer
1264 * - Remove CDbBlast dependencies in favor of new traceback class.
1265 *
1266 * Revision 1.58 2006/06/29 16:25:24 camacho
1267 * Changed BlastHitSavingOptionsNew signature
1268 *
1269 * Revision 1.57 2006/06/05 13:34:05 madden
1270 * Changes to remove [GS]etMatrixPath and use callback instead
1271 *
1272 * Revision 1.56 2006/03/21 21:01:07 camacho
1273 * + interruptible api support
1274 *
1275 * Revision 1.55 2006/02/15 18:41:11 madden
1276 * Made changes to CBlastTraceBackTest::testBLASTNTraceBack to reflect
1277 * small changes to the starting point due to changes to the
1278 * BLAST_CheckStartForGappedAlignment routine.
1279 * (from Mike Gertz).
1280 *
1281 * Revision 1.54 2005/12/16 20:51:50 camacho
1282 * Diffuse the use of CSearchMessage, TQueryMessages, and TSearchMessages
1283 *
1284 * Revision 1.53 2005/10/14 13:47:32 camacho
1285 * Fixes to pacify icc compiler
1286 *
1287 * Revision 1.52 2005/08/29 15:00:45 camacho
1288 * Update calls to BLAST_MainSetUp to reflect changes in the signature
1289 *
1290 * Revision 1.51 2005/07/19 15:56:39 bealer
1291 * - Undo incorrect Seq-id changes.
1292 *
1293 * Revision 1.50 2005/07/19 13:49:12 madden
1294 * Fixes for dust change
1295 *
1296 * Revision 1.49 2005/07/18 19:17:51 bealer
1297 * - Fix accessions.
1298 *
1299 * Revision 1.48 2005/07/18 17:04:43 bealer
1300 * - Change expired GIs to (hopefully) longer lasting unversioned accessions.
1301 *
1302 * Revision 1.47 2005/06/09 20:37:06 camacho
1303 * Use new private header blast_objmgr_priv.hpp
1304 *
1305 * Revision 1.46 2005/05/25 13:47:34 camacho
1306 * Use CBLAST_SequenceBlk instead of BLAST_SequenceBlk*
1307 *
1308 * Revision 1.45 2005/05/24 20:05:17 camacho
1309 * Changed signature of SetupQueries and SetupQueryInfo
1310 *
1311 * Revision 1.44 2005/05/23 15:53:19 dondosha
1312 * Special case for preliminary hitlist size in RPS BLAST - hence no need for 2 extra parameters in SBlastHitsParametersNew
1313 *
1314 * Revision 1.43 2005/05/17 16:05:23 madden
1315 * Fixes for BLOSUM62 matrix change
1316 *
1317 * Revision 1.42 2005/05/16 12:29:15 madden
1318 * use SBlastHitsParameters in Blast_HSPListCollectorInit and Blast_HSPListCollectorInit[MT]
1319 *
1320 * Revision 1.41 2005/04/27 20:09:38 dondosha
1321 * Extra argument has been added to BLAST_ComputeTraceback for PHI BLAST, pass NULL
1322 *
1323 * Revision 1.40 2005/04/18 14:01:55 camacho
1324 * Updates following BlastSeqSrc reorganization
1325 *
1326 * Revision 1.39 2005/04/06 21:26:37 dondosha
1327 * GapEditBlock structure and redundant fields in BlastHSP have been removed
1328 *
1329 * Revision 1.38 2005/03/31 13:45:58 camacho
1330 * BLAST options API clean-up
1331 *
1332 * Revision 1.37 2005/03/14 19:44:59 madden
1333 * Adjustments for recent dust change
1334 *
1335 * Revision 1.36 2005/03/04 17:20:45 bealer
1336 * - Command line option support.
1337 *
1338 * Revision 1.35 2005/01/27 20:09:54 madden
1339 * Adjustments for dust fixes from Richa
1340 *
1341 * Revision 1.34 2005/01/19 19:07:17 coulouri
1342 * do not cast pointers to smaller integer types
1343 *
1344 * Revision 1.33 2005/01/18 14:55:05 camacho
1345 * Changes for taking new tie-breakers for score comparison into account
1346 *
1347 * Revision 1.32 2005/01/06 15:43:25 camacho
1348 * Make use of modified signature to blast::SetupQueries
1349 *
1350 * Revision 1.31 2004/12/21 17:26:11 dondosha
1351 * BLAST_ComputeTraceback has a new RPSInfo argument; pass it as NULL for non-RPS searches
1352 *
1353 * Revision 1.30 2004/12/08 16:43:55 dondosha
1354 * Removed dead code
1355 *
1356 * Revision 1.29 2004/11/30 17:05:35 dondosha
1357 * Added test to check that HSPs are not removed based on e-value cutoff before linking
1358 *
1359 * Revision 1.28 2004/11/18 21:35:43 dondosha
1360 * Added testHSPUpdateWithTraceback
1361 *
1362 * Revision 1.27 2004/11/17 21:02:01 camacho
1363 * Add error checking to BlastSeqSrc initialization
1364 *
1365 * Revision 1.26 2004/10/19 16:38:14 dondosha
1366 * Changed order of results due to sorting of HSPs by score in BLAST_LinkHsps
1367 *
1368 * Revision 1.25 2004/09/29 17:20:08 papadopo
1369 * shuffle the order of expected HSPs to account for changes in HSP linking
1370 *
1371 * Revision 1.24 2004/09/21 13:54:28 dondosha
1372 * Adjusted results due to fix in the engine that sorts HSPs in HSP lists by e-value
1373 *
1374 * Revision 1.23 2004/07/19 15:04:22 dondosha
1375 * Renamed multiseq_src to seqsrc_multiseq, seqdb_src to seqsrc_seqdb
1376 *
1377 * Revision 1.22 2004/07/06 15:58:45 dondosha
1378 * Use EBlastProgramType enumeration type for program when calling C functions
1379 *
1380 * Revision 1.21 2004/06/28 13:43:51 madden
1381 * Use NULL for unused filter_slp in BLAST_MainSetUp
1382 *
1383 * Revision 1.20 2004/06/24 16:10:26 madden
1384 * Add test testBLASTNTraceBackLargeXDrop that core-dumped before blast_gapalign.c fix from Jason Papadopoulos on 6/24/04
1385 *
1386 * Revision 1.19 2004/06/22 16:46:19 camacho
1387 * Changed the blast_type_* definitions for the EBlastProgramType enumeration.
1388 *
1389 * Revision 1.18 2004/06/08 15:24:34 dondosha
1390 * Use BlastHSPStream interface
1391 *
1392 * Revision 1.17 2004/05/12 12:20:34 madden
1393 * Add (NULL) psi_options to call to BLAST_ComputeTraceback
1394 *
1395 * Revision 1.16 2004/05/07 15:44:42 papadopo
1396 * fill in and use BlastScoringParameters instead of BlastScoringOptions
1397 *
1398 * Revision 1.15 2004/05/05 15:29:39 dondosha
1399 * Renamed functions in blast_hits.h accordance with new convention Blast_[StructName][Task]
1400 *
1401 * Revision 1.14 2004/04/30 16:54:05 dondosha
1402 * Changed a number of function names to have the same conventional Blast_ prefix
1403 *
1404 * Revision 1.13 2004/04/19 13:10:47 madden
1405 * Remove duplicate lines/checks
1406 *
1407 * Revision 1.12 2004/03/23 22:07:01 camacho
1408 * Fix memory leaks
1409 *
1410 * Revision 1.11 2004/03/23 16:10:34 camacho
1411 * Minor changes to CTestObjMgr
1412 *
1413 * Revision 1.10 2004/03/18 15:15:25 dondosha
1414 * Changed duplicate variable names due to compiler warnings
1415 *
1416 * Revision 1.9 2004/03/15 20:02:38 dondosha
1417 * Use BLAST_GapAlignSetUp function to do traceback specific setup; database and two sequences traceback engines merged; corrected order of arguments in BOOST_REQUIRE_EQUAL calls
1418 *
1419 * Revision 1.8 2004/03/11 21:17:16 camacho
1420 * Fix calls to BlastHitSavingParametersNew
1421 *
1422 * Revision 1.7 2004/03/09 18:58:56 dondosha
1423 * Added extension parameters argument to BlastHitSavingParametersNew calls
1424 *
1425 * Revision 1.6 2004/03/01 14:14:50 madden
1426 * Add check for number of identitical letters in final alignment
1427 *
1428 * Revision 1.5 2004/02/27 21:26:56 madden
1429 * Cleanup testBLASTNTraceBack
1430 *
1431 * Revision 1.4 2004/02/27 21:22:25 madden
1432 * Add tblastn test: testTBLASTNTraceBack
1433 *
1434 * Revision 1.3 2004/02/27 20:34:34 madden
1435 * Add setUp and tearDown routines, use tearDown for deallocation
1436 *
1437 * Revision 1.2 2004/02/27 20:14:24 madden
1438 * Add protein-protein test: testBLASTPTraceBack
1439 *
1440 * Revision 1.1 2004/02/27 19:45:35 madden
1441 * Unit test for traceback (mostly BlastHSPListGetTraceback)
1442 *
1443 *
1444 *
1445 * ===========================================================================
1446 */
Declares the CBl2Seq (BLAST 2 Sequences) class.
BlastSeqLoc * BlastSeqLocFree(BlastSeqLoc *loc)
Deallocate all BlastSeqLoc objects in a chain.
Definition: blast_filter.c:737
Structures and functions prototypes used for BLAST gapped extension.
BlastGapAlignStruct * BLAST_GapAlignStructFree(BlastGapAlignStruct *gap_align)
Deallocates memory in the BlastGapAlignStruct structure.
BlastHSP * Blast_HSPNew(void)
Allocate and zeros out memory for an HSP structure.
Definition: blast_hits.c:141
BlastHSPResults * Blast_HSPResultsFree(BlastHSPResults *results)
Deallocate memory for BLAST results.
Definition: blast_hits.c:3364
BlastHSPList * Blast_HSPListNew(Int4 hsp_max)
Creates HSP list structure with a default size HSP array.
Definition: blast_hits.c:1558
BlastHSP * Blast_HSPFree(BlastHSP *hsp)
Deallocate memory for an HSP structure.
Definition: blast_hits.c:130
void Blast_HSPListSortByScore(BlastHSPList *hsp_list)
Sort the HSPs in an HSP list by score.
Definition: blast_hits.c:1374
BlastHSPWriter * BlastHSPWriterNew(BlastHSPWriterInfo **writer_info, BlastQueryInfo *query_info, BLAST_SequenceBlk *query)
A generic function to create writer.
Declaration of ADT to save and retrieve lists of HSPs in the BLAST engine.
int BlastHSPStreamWrite(BlastHSPStream *hsp_stream, BlastHSPList **hsp_list)
Invokes the user-specified write function for this BlastHSPStream implementation.
BlastHSPStream * BlastHSPStreamFree(BlastHSPStream *hsp_stream)
Frees the BlastHSPStream structure by invoking the destructor function set by the user-defined constr...
BlastHSPStream * BlastHSPStreamNew(EBlastProgramType program, const BlastExtensionOptions *extn_opts, Boolean sort_on_read, Int4 num_queries, BlastHSPWriter *writer)
Initialize the HSP stream.
Common definitions for protein and nucleotide lookup tables.
Blast_Message * Blast_MessageFree(Blast_Message *blast_msg)
Deallocates message memory.
Definition: blast_message.c:80
Declares the CBlastNucleotideOptionsHandle class.
Definitions which are dependant on the NCBI C++ Object Manager.
BlastHitSavingOptions * BlastHitSavingOptionsFree(BlastHitSavingOptions *options)
Deallocate memory for BlastHitSavingOptions.
Int2 BlastScoringOptionsNew(EBlastProgramType program, BlastScoringOptions **options)
Allocate memory for BlastScoringOptions and fill with default values.
BlastExtensionOptions * BlastExtensionOptionsFree(BlastExtensionOptions *options)
Deallocate memory for BlastExtensionOptions.
Int2 BlastHitSavingOptionsNew(EBlastProgramType program, BlastHitSavingOptions **options, Boolean gapped_calculation)
Allocate memory for BlastHitSavingOptions.
BlastScoringOptions * BlastScoringOptionsFree(BlastScoringOptions *options)
Deallocate memory for BlastScoringOptions.
Int2 BlastExtensionOptionsNew(EBlastProgramType program, BlastExtensionOptions **options, Boolean gapped)
Allocate memory for BlastExtensionOptions and fill with default values.
Declares the CBlastOptionsHandle and CBlastOptionsFactory classes.
BlastHitSavingParameters * BlastHitSavingParametersFree(BlastHitSavingParameters *parameters)
Deallocate memory for BlastHitSavingOptions*.
BlastEffectiveLengthsParameters * BlastEffectiveLengthsParametersFree(BlastEffectiveLengthsParameters *parameters)
Deallocate memory for BlastEffectiveLengthsParameters*.
BlastExtensionParameters * BlastExtensionParametersFree(BlastExtensionParameters *parameters)
Deallocate memory for BlastExtensionParameters.
BlastScoringParameters * BlastScoringParametersFree(BlastScoringParameters *parameters)
Deallocate memory for BlastScoringParameters.
EBlastProgramType
Defines the engine's notion of the different applications of the BLAST algorithm.
Definition: blast_program.h:72
@ eBlastTypeBlastn
Definition: blast_program.h:74
@ eBlastTypeTblastn
Definition: blast_program.h:77
Declares the CBlastProteinOptionsHandle class.
Defines interface for retrieving sequence identifiers.
Int4 BlastSeqSrcGetSeqLen(const BlastSeqSrc *seq_src, void *oid)
Retrieve sequence length (number of residues/bases)
Definition: blast_seqsrc.c:281
BlastSeqSrc * BlastSeqSrcCopy(const BlastSeqSrc *seq_src)
Copy function: needed to guarantee thread safety.
Definition: blast_seqsrc.c:138
char * BlastSeqSrcGetInitError(const BlastSeqSrc *seq_src)
Function to retrieve NULL terminated string containing the description of an initialization error or ...
Definition: blast_seqsrc.c:159
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
Int2 BlastSeqSrcGetSequence(const BlastSeqSrc *seq_src, BlastSeqSrcGetSeqArg *getseq_arg)
Retrieve an individual sequence.
Definition: blast_seqsrc.c:271
BlastSeqSrcSetRangesArg * BlastSeqSrcSetRangesArgFree(BlastSeqSrcSetRangesArg *arg)
free setrangearg
Definition: blast_seqsrc.c:455
Utilities initialize/setup BLAST.
Int2 BLAST_MainSetUp(EBlastProgramType program_number, const QuerySetUpOptions *qsup_options, const BlastScoringOptions *scoring_options, BLAST_SequenceBlk *query_blk, const BlastQueryInfo *query_info, double scale_factor, BlastSeqLoc **lookup_segments, BlastMaskLoc **mask, BlastScoreBlk **sbpp, Blast_Message **blast_message, GET_MATRIX_PATH get_path)
"Main" setup routine for BLAST.
Definition: blast_setup.c:563
Int2 BLAST_GapAlignSetUp(EBlastProgramType program_number, const BlastSeqSrc *seq_src, const BlastScoringOptions *scoring_options, const BlastEffectiveLengthsOptions *eff_len_options, const BlastExtensionOptions *ext_options, const BlastHitSavingOptions *hit_options, BlastQueryInfo *query_info, BlastScoreBlk *sbp, BlastScoringParameters **score_params, BlastExtensionParameters **ext_params, BlastHitSavingParameters **hit_params, BlastEffectiveLengthsParameters **eff_len_params, BlastGapAlignStruct **gap_align)
Set up the auxiliary structures for gapped alignment / traceback only.
Definition: blast_setup.c:888
BlastScoreBlk * BlastScoreBlkFree(BlastScoreBlk *sbp)
Deallocates BlastScoreBlk as well as all associated structures.
Definition: blast_stat.c:965
Functions to do gapped alignment with traceback.
Int2 Blast_HSPUpdateWithTraceback(BlastGapAlignStruct *gap_align, BlastHSP *hsp)
Modifies the HSP data after the final gapped alignment.
BlastSeqSrcSetRangesArg * BLAST_SetupPartialFetching(EBlastProgramType program_number, BlastSeqSrc *seq_src, const BlastHSPList **hsp_list, Int4 num_hsplists)
Attempts to set up partial fetching, if it fails (e.g.
EBlastEncoding Blast_TracebackGetEncoding(EBlastProgramType program_number)
Get the subject sequence encoding type for the traceback, given a program number.
Int2 BLAST_ComputeTraceback(EBlastProgramType program_number, BlastHSPStream *hsp_stream, BLAST_SequenceBlk *query, BlastQueryInfo *query_info, const BlastSeqSrc *seq_src, BlastGapAlignStruct *gap_align, BlastScoringParameters *score_params, const BlastExtensionParameters *ext_params, BlastHitSavingParameters *hit_params, BlastEffectiveLengthsParameters *eff_len_params, const BlastDatabaseOptions *db_options, const PSIBlastOptions *psi_options, const BlastRPSInfo *rps_info, SPHIPatternSearchBlk *pattern_blk, BlastHSPResults **results, TInterruptFnPtr interrupt_search, SBlastProgress *progress_info)
Given the preliminary alignment results from a database search, redo the gapped alignment with traceb...
Definitions of special type used in BLAST.
EProgram
This enumeration is to evolve into a task/program specific list that specifies sets of default parame...
Definition: blast_types.hpp:56
@ eBlastn
Nucl-Nucl (traditional blastn)
Definition: blast_types.hpp:58
@ eTblastn
Protein-Translated nucl.
Definition: blast_types.hpp:61
BLAST_SequenceBlk * BlastSequenceBlkFree(BLAST_SequenceBlk *seq_blk)
Deallocate memory for a sequence block.
Definition: blast_util.c:245
BOOST_AUTO_TEST_SUITE_END() static int s_GetSegmentFlags(const CBioseq &bioseq)
Declares the CBlastxOptionsHandle class.
Wrapper class for BLAST_SequenceBlk .
Definition: blast_aux.hpp:309
Runs the BLAST algorithm between 2 sequences.
Definition: bl2seq.hpp:58
Handle to the nucleotide-nucleotide options to the BLAST algorithm.
Encapsulates ALL the BLAST algorithm's options.
Handle to the protein-protein options to the BLAST algorithm.
Wrapper class for BlastQueryInfo .
Definition: blast_aux.hpp:311
Query Vector.
Definition: sseqloc.hpp:276
Wrapper class for BlastSeqSrc .
Definition: blast_aux.hpp:350
NCBI C++ Object Manager dependant implementation of IQueryFactory.
Search Results for All Queries.
CSeqDB.
Definition: seqdb.hpp:161
@ eNucleotide
Definition: seqdb.hpp:175
Handle to the protein-translated nucleotide options to the BLAST algorithm.
CRef< blast::CBlastSearchQuery > CreateBlastSearchQuery(objects::CSeq_id &id, objects::ENa_strand s=objects::eNa_strand_unknown)
static CTestObjMgr & Instance()
Definition: test_objmgr.cpp:69
static BlastHSPStream * x_MakeStream(const CBlastOptions &opt)
void x_SetupMain(const CBlastOptions &opt, const CBLAST_SequenceBlk &query_blk, const CBlastQueryInfo &query_info)
void x_SetupGapAlign(const CBlastOptions &opt, const BlastSeqSrc *seq_src, const CBlastQueryInfo &query_info)
void x_ComputeTracebak(const CBlastOptions &opt, BlastHSPStream *hsp_stream, const CBLAST_SequenceBlk &query_blk, const CBlastQueryInfo &query_info, const BlastSeqSrc *seq_src, BlastHSPResults **results)
typedef for the messages for an entire BLAST search, which could be comprised of multiple query seque...
@ eNoCompositionBasedStats
Don't use composition based statistics.
Declares the CDiscNucleotideOptionsHandle class.
static const int kOffset
BlastHitSavingOptions * GetHitSaveOpts() const
Returns BlastHitSavingOptions for eLocal objects, NULL for eRemote.
void SetCompositionBasedStats(ECompoAdjustModes mode)
BlastSeqSrc * Get() const
Definition: blast_aux.hpp:350
const CBlastOptionsHandle & GetOptionsHandle() const
Retrieve the options handle.
Definition: bl2seq.hpp:285
void SetQueryCovHspPerc(double p)
Sets QueryCovHspPerc.
void SetDbLength(Int8 len)
Sets DbLength.
void SetMatchReward(int r)
Sets MatchReward.
BlastExtensionOptions * GetExtnOpts() const
Returns BlastExtensionOptions for eLocal objects, NULL for eRemote.
BlastSeqSrc * SeqDbBlastSeqSrcInit(const string &dbname, bool is_prot, Uint4 first_seq=0, Uint4 last_seq=0, Int4 mask_algo_id=-1, ESubjectMaskingType mask_type=eNoSubjMasking)
Initialize the sequence source structure.
CStructWrapper< TData > * WrapStruct(TData *obj, TData *(*del)(TData *))
Auxiliary function to create a CStructWrapper for a pointer to an object.
void SetupQueries(TSeqLocVector &queries, BlastQueryInfo *qinfo, BLAST_SequenceBlk **seqblk, EBlastProgramType prog, objects::ENa_strand strand_opt, TSearchMessages &messages)
Populates BLAST_SequenceBlk with sequence data for use in CORE BLAST.
objects::ENa_strand GetStrandOption() const
void SetEffectiveSearchSpace(Int8 eff)
Sets EffectiveSearchSpace.
CRef< CSearchResultSet > Run()
Run the traceback search.
static CBlastOptionsHandle * Create(EProgram program, EAPILocality locality=CBlastOptions::eLocal)
Creates an options handle object configured with default options for the requested program,...
CBlastOptions & SetOptions()
Returns a reference to the internal options class which this object is a handle for.
void SetTraditionalBlastnDefaults()
Sets TraditionalBlastnDefaults.
BlastSeqSrc * MultiSeqBlastSeqSrcInit(TSeqLocVector &seq_vector, EBlastProgramType program, bool dbscan_mode=false)
Initialize the sequence source structure.
BlastEffectiveLengthsOptions * GetEffLenOpts() const
Returns BlastEffectiveLengthsOptions for eLocal objects, NULL for eRemote.
void SetMaxHspsPerSubject(int m)
Sets MaxHspPerSubjectQueryPair.
EBlastProgramType GetProgramType() const
Returns the CORE BLAST notion of program type.
const CBlastOptions & GetOptions() const
Return the object which this object is a handle for.
BlastScoringOptions * GetScoringOpts() const
Returns BlastScoringOptions for eLocal objects, NULL for eRemote.
char * BlastFindMatrixPath(const char *matrix_name, Boolean is_prot)
Returns the path to a specified matrix.
void SetGapXDropoffFinal(double x)
Sets GapXDropoffFinal.
QuerySetUpOptions * GetQueryOpts() const
Returns QuerySetUpOptions for eLocal objects, NULL for eRemote.
const TSeqLocVector & GetQueries() const
Retrieve a vector of query sequences.
Definition: bl2seq.hpp:244
CBlastOptionsHandle & SetOptionsHandle()
Set the options handle.
Definition: bl2seq.hpp:278
void SetupQueryInfo(TSeqLocVector &queries, EBlastProgramType prog, objects::ENa_strand strand_opt, BlastQueryInfo **qinfo)
Allocates the query information structure and fills the context offsets, in case of multiple queries,...
BlastDatabaseOptions * GetDbOpts() const
Returns BlastDatabaseOptions for eLocal objects, NULL for eRemote.
void SetSegFiltering(bool val=true)
size_type GetNumResults() const
Return the number of results contained by this object.
const TSeqLocVector & GetSubjects() const
Retrieve a vector of subject sequences.
Definition: bl2seq.hpp:272
#define ITERATE(Type, Var, Cont)
ITERATE macro to sequence through container elements.
Definition: ncbimisc.hpp:815
#define NULL
Definition: ncbistd.hpp:225
static char GetPathSeparator(void)
Get path separator symbol specific for the current platform.
Definition: ncbifile.cpp:433
static string ConcatPath(const string &first, const string &second)
Concatenate two parts of the path for the current OS.
Definition: ncbifile.cpp:776
int32_t Int4
4-byte (32-bit) signed integer
Definition: ncbitype.h:102
uint32_t Uint4
4-byte (32-bit) unsigned integer
Definition: ncbitype.h:103
int64_t Int8
8-byte (64-bit) signed integer
Definition: ncbitype.h:104
ENa_strand
strand of nucleic acid
Definition: Na_strand_.hpp:64
@ eNa_strand_unknown
Definition: Na_strand_.hpp:65
@ eNa_strand_both
in forward orientation
Definition: Na_strand_.hpp:68
Implementation of a number of BlastHSPWriters to save hits from a BLAST search, and subsequently retu...
BlastHSPCollectorParams * BlastHSPCollectorParamsNew(const BlastHitSavingOptions *hit_options, Int4 compositionBasedStats, Boolean gapped_calculation)
Sets up parameter set for use by collector.
BlastHSPWriterInfo * BlastHSPCollectorInfoNew(BlastHSPCollectorParams *params)
WriterInfo to create a default writer: the collecter.
int i
int len
Utility functions for lookup table generation.
static char * prog
Definition: mdb_load.c:33
void CheckForBlastSeqSrcErrors(const BlastSeqSrc *seqsrc)
Magic spell ;-) needed for some weird compilers... very empiric.
#define FALSE
bool replacment for C indicating false.
Definition: ncbi_std.h:101
Defines classes: CDirEntry, CFile, CDir, CSymLink, CMemoryFile, CFileUtil, CFileLock,...
Defines: CTimeFormat - storage class for time format.
double r(size_t dimension_, const Int4 *score_, const double *prob_, double theta_)
The Object manager core.
NOTE: This file contains work in progress and the APIs are likely to change, please do not rely on th...
Utilities to develop and debug unit tests for BLAST.
NOTE: This file contains work in progress and the APIs are likely to change, please do not rely on th...
Defines a concrete strategy for the IBlastSeqInfoSrc interface for sequence identifiers retrieval fro...
Implementation of the BlastSeqSrc interface for a vector of sequence locations.
Implementation of the BlastSeqSrc interface using the C++ BLAST databases API.
NOTE: This file contains work in progress and the APIs are likely to change, please do not rely on th...
vector< SSeqLoc > TSeqLocVector
Vector of sequence locations.
Definition: sseqloc.hpp:129
Uint1 * sequence_start
Start of sequence, usually one byte before sequence as that byte is a NULL sentinel byte.
Definition: blast_def.h:244
Parameters for setting up effective lengths and search spaces.
Options used for gapped extension These include: a.
Int4 compositionBasedStats
mode of compositional adjustment to use; if zero then compositional adjustment is not used
Computed values used as parameters for gapped alignments.
Structure supporting the gapped alignment.
Int4 query_stop
query end offseet of current alignment
Int4 subject_start
subject start offset current alignment
Int4 query_start
query start offset of current alignment
Int4 subject_stop
subject end offset of current alignment
Int4 score
Return value: alignment score.
GapEditScript * edit_script
The traceback (gap) information.
The structure to hold all HSPs for a given sequence after the gapped alignment.
Definition: blast_hits.h:153
Boolean do_not_reallocate
Is reallocation of the hsp_array allowed?
Definition: blast_hits.h:161
Int4 oid
The ordinal id of the subject sequence this HSP list is for.
Definition: blast_hits.h:154
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
Int4 hsp_max
The maximal number of HSPs allowed to be saved.
Definition: blast_hits.h:160
Int4 allocated
The allocated size of the hsp_array.
Definition: blast_hits.h:159
The structure to contain all BLAST results, for multiple queries.
Definition: blast_hits.h:183
BlastHitList ** hitlist_array
Array of results for individual query sequences.
Definition: blast_hits.h:185
Default implementation of BlastHSPStream.
A wrap of data structure used to create a writer.
ADT definition of BlastHSPWriter.
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
BlastSeg query
Query sequence info.
Definition: blast_hits.h:131
Int4 context
Context number of query.
Definition: blast_hits.h:133
Int4 num
How many HSP's are linked together for sum statistics evaluation? If unset (0), this HSP is not part ...
Definition: blast_hits.h:135
BlastSeg subject
Subject sequence info.
Definition: blast_hits.h:132
GapEditScript * gap_info
ALL gapped alignment is here.
Definition: blast_hits.h:134
Int4 score
This HSP's raw score.
Definition: blast_hits.h:127
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
Options used when evaluating and saving hits These include: a.
Parameter block that contains a pointer to BlastHitSavingOptions and the values derived from it.
Structure used for scoring calculations.
Definition: blast_stat.h:177
Scoring options block Used to produce the BlastScoreBlk structure This structure may be needed for lo...
Boolean gapped_calculation
gap-free search if FALSE
Scoring parameters block Contains scoring-related information that is actually used for the blast sea...
Int4 end
End of hsp.
Definition: blast_hits.h:99
Int4 gapped_start
Where the gapped extension started.
Definition: blast_hits.h:100
Int2 frame
Translation frame.
Definition: blast_hits.h:97
Int4 offset
Start of hsp.
Definition: blast_hits.h:98
Used to hold a set of positions, mostly used for filtering.
Definition: blast_def.h:204
Structure used as the second argument to functions satisfying the GetSeqBlkFnPtr signature,...
Definition: blast_seqsrc.h:257
Int4 oid
Oid in BLAST database, index in an array of sequences, etc [in].
Definition: blast_seqsrc.h:259
Boolean reset_ranges
This option allows the BLAST engine to communicate with the BlastSeqSrc that the offset ranges for a ...
Definition: blast_seqsrc.h:273
EBlastEncoding encoding
Encoding of sequence, i.e.
Definition: blast_seqsrc.h:263
Boolean check_oid_exclusion
Check whether an OID is excluded due to overlapping filtering.
Definition: blast_seqsrc.h:279
BlastSeqSrcSetRangesArg * ranges
Definition: blast_seqsrc.h:286
BLAST_SequenceBlk * seq
Sequence to return, if NULL, it should allocated by GetSeqBlkFnPtr (using BlastSeqBlkNew or BlastSetU...
Definition: blast_seqsrc.h:284
Structure used as the argument to function SetRanges.
Definition: blast_seqsrc.h:208
Int4 * ranges
Ranges in sorted order [in].
Definition: blast_seqsrc.h:219
Int4 num_ranges
Number of actual ranges contained.
Definition: blast_seqsrc.h:216
Complete type definition of Blast Sequence Source ADT.
Definition: blast_seqsrc.c:43
Structure to hold the a message from the core of the BLAST engine.
Definition: blast_message.h:70
Edit script: linked list of correspondencies between two sequences.
Definition: gapinfo.h:57
static string query
Declares the CTBlastnOptionsHandle class.
Utility stuff for more convenient using of Boost.Test library.
Defines location of test data folder at NCBI.
static const char * NCBI_GetTestDataPath(void)
Get the directory where test data is stored at NCBI.
NOTE: This file contains work in progress and the APIs are likely to change, please do not rely on th...
BlastHSPList * s_GetHSPList(int num_hsps, int oid, int len, int range_len)
BOOST_AUTO_TEST_CASE(testHSPUpdateWithTraceback)
static CS_CONTEXT * context
Definition: will_convert.c:21
void free(voidpf ptr)
voidp malloc(uInt size)
voidp calloc(uInt items, uInt size)
Modified on Wed Jul 17 13:19:12 2024 by modify_doxy.py rev. 669887