BitMagic-C++
xsample07.cpp
Go to the documentation of this file.
1/*
2Copyright(c) 2020 Anatoliy Kuznetsov(anatoliy_kuznetsov at yahoo.com)
3
4Licensed under the Apache License, Version 2.0 (the "License");
5you may not use this file except in compliance with the License.
6You may obtain a copy of the License at
7
8 http://www.apache.org/licenses/LICENSE-2.0
9
10Unless required by applicable law or agreed to in writing, software
11distributed under the License is distributed on an "AS IS" BASIS,
12WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13See the License for the specific language governing permissions and
14limitations under the License.
15
16For more information please visit: http://bitmagic.io
17*/
18
19/** \example xsample07.cpp
20
21 Use of bvector<> for k-mer fingerprint K should be short,
22 no minimizers are used here (the approach can be used to create a
23 scheme with minimizers).
24
25 This example:
26 - loads a FASTA file (single DNA molecule is expected)
27 - generates K-mers for the specified K
28 builds a fingerprint bit-vector (presence-absence vector)
29 - runs k-mer counting (fastest method uses sort based algorithm)
30 builds a k-mer frequency compressed sparse vector (TF-vector)
31 - shows how to use TF vector to pick top N% of most frequent k-mers
32 and exclude them (as over-represented)
33
34 \sa bm::sparse_vector
35 \sa bm::rsc_sparse_vector
36 \sa bm::rank_range_split
37 \sa bm::rsc_sparse_vector<>::merge_not_null
38
39*/
40
41/*! \file xsample07.cpp
42 \brief Example: Use of bvector<> for k-mer fingerprint
43 K should be short, no minimizers here
44*/
45
46#include <assert.h>
47#include <stdlib.h>
48
49#include <iostream>
50#include <vector>
51#include <map>
52#include <algorithm>
53#include <utility>
54
55#include <future>
56#include <thread>
57#include <mutex>
58#include <atomic>
59
60#include "bm64.h" // use 48-bit vectors
61#include "bmalgo.h"
62#include "bmserial.h"
63#include "bmaggregator.h"
64#include "bmsparsevec_compr.h"
65#include "bmsparsevec_algo.h"
66#include "bmundef.h" /* clear the pre-proc defines from BM */
67
68// BitMagic utilities for debug and timings
69#include "bmdbg.h"
70#include "bmtimer.h"
71
72#include "dna_finger.h"
73
74using namespace std;
75
76
77
78// Arguments
79//
80std::string ifa_name;
81std::string ikd_name;
82std::string ikd_counts_name;
83std::string kh_name;
84std::string ikd_rep_name;
85std::string ikd_freq_name;
86bool is_diag = false;
87bool is_timing = false;
88bool is_bench = false;
89unsigned ik_size = 8;
90unsigned parallel_jobs = 4;
91unsigned f_percent = 5; // percent of k-mers we try to clear as over-represented
92
93#include "cmd_args.h"
94
95
96
97// Global types
98//
99typedef std::vector<char> vector_char_type;
103typedef std::map<unsigned, unsigned> histogram_map_u32;
104
105
106// Global vars
107//
110std::atomic_ullong k_mer_progress_count(0);
111
112
113/// really simple FASTA parser (one entry per file)
114///
115static
116int load_FASTA(const std::string& fname, vector_char_type& seq_vect)
117{
118 bm::chrono_taker tt1(cout, "1. Parse FASTA", 1, &timing_map);
119
120 seq_vect.resize(0);
121 std::ifstream fin(fname.c_str(), std::ios::in);
122 if (!fin.good())
123 return -1;
124
125 std::string line;
126 for (unsigned i = 0; std::getline(fin, line); ++i)
127 {
128 if (line.empty() ||
129 line.front() == '>')
130 continue;
131 for (std::string::iterator it = line.begin(); it != line.end(); ++it)
132 seq_vect.push_back(*it);
133 } // for
134 return 0;
135}
136
137inline
138bool get_DNA_code(char bp, bm::id64_t& dna_code)
139{
140 switch (bp)
141 {
142 case 'A':
143 dna_code = 0; // 00
144 break;
145 case 'T':
146 dna_code = 1; // 01
147 break;
148 case 'G':
149 dna_code = 2; // 10
150 break;
151 case 'C':
152 dna_code = 3; // 11
153 break;
154 default: // ambiguity codes are ignored (for simplicity)
155 return false;
156 }
157 return true;
158}
159
160/// Calculate k-mer as an unsigned long integer
161///
162/// @return true - if k-mer is "true" (not 'NNNNNN')
163///
164inline
165bool get_kmer_code(const char* dna,
166 size_t pos, unsigned k_size,
167 bm::id64_t& k_mer)
168{
169 // generate k-mer
170 //
171 bm::id64_t k_acc = 0;
172 unsigned shift = 0;
173 dna += pos;
174 for (size_t i = 0; i < k_size; ++i)
175 {
176 char bp = dna[i];
177 bm::id64_t dna_code;
178 bool valid = get_DNA_code(bp, dna_code);
179 if (!valid)
180 return false;
181 k_acc |= (dna_code << shift); // accumulate new code within 64-bit accum
182 shift += 2; // each DNA base pair needs 2-bits to store
183 } // for i
184 k_mer = k_acc;
185 return true;
186}
187
188
189/// Translate integer code to DNA letter
190///
191inline
192char int2DNA(unsigned code)
193{
194 static char lut[] = { 'A', 'T', 'G', 'C', 'N', '$' };
195 if (code < 5)
196 return lut[code];
197 assert(0);
198 return 'N';
199}
200
201/// Translate k-mer code into ATGC DNA string
202///
203/// @param dna - target string
204/// @param k_mer - k-mer code
205/// @param k_size -
206inline
207void translate_kmer(std::string& dna, bm::id64_t kmer_code, unsigned k_size)
208{
209 dna.resize(k_size);
210 for (size_t i = 0; i < k_size; ++i)
211 {
212 unsigned dna_code = unsigned(kmer_code & 3);
213 char bp = int2DNA(dna_code);
214 dna[i] = bp;
215 kmer_code >>= 2;
216 } // for i
217 assert(!kmer_code);
218}
219
220
221/// QA function to validate if reverse k-mer decode gives the same string
222///
223inline
224void validate_k_mer(const char* dna,
225 size_t pos, unsigned k_size,
226 bm::id64_t k_mer)
227{
228 for (size_t i = 0; i < k_size; ++i)
229 {
230 char bp = dna[pos+i];
231 unsigned dna_code = unsigned(k_mer & 3ull);
232 char bp_c = int2DNA(dna_code);
233 if (bp != bp_c)
234 {
235 if (bp == 'N' && bp_c != 'A')
236 {
237 cerr << bp << " " << bp_c << endl;
238 cerr << "Error! N code mismatch at pos = " << pos+i
239 << endl;
240 exit(1);
241 }
242 }
243 k_mer >>= 2;
244 } // for i
245 if (k_mer)
246 {
247 cerr << "Error! non-zero k-mer remainder at " << pos << endl;
248 exit(1);
249 }
250}
251
252/// Auxiliary function to do sort+unique on a vactor of ints
253/// removes duplicate elements
254///
255template<typename VECT>
256void sort_unique(VECT& vect)
257{
258 std::sort(vect.begin(), vect.end());
259 auto last = std::unique(vect.begin(), vect.end());
260 vect.erase(last, vect.end());
261}
262
263
264/// Auxiliary function to do sort+unique on a vactor of ints and save results
265/// in a counts vector
266///
267template<typename VECT, typename COUNT_VECT>
268void sort_count(VECT& vect, COUNT_VECT& cvect)
269{
270 if (!vect.size())
271 return;
272 std::sort(vect.begin(), vect.end());
273 typename VECT::value_type prev = vect[0];
274 typename COUNT_VECT::value_type cnt = 1;
275 auto vsize = vect.size();
276 size_t i = 1;
277 for (; i < vsize; ++i)
278 {
279 auto v = vect[i];
280 if (v == prev)
281 {
282 ++cnt;
283 continue;
284 }
285 cvect.inc_not_null(prev, cnt);
286 prev = v; cnt = 1;
287 } // for i
288
289 cvect.inc_not_null(prev, cnt);
290 assert(cvect.in_sync());
291}
292
293/**
294 This function turns each k-mer into an integer number and encodes it
295 in a bit-vector (presense vector)
296 The natural limitation here is that integer has to be less tha 48-bits
297 (limitations of bm::bvector<>)
298 This method build a presense k-mer fingerprint vector which can be
299 used for Jaccard distance comparison.
300
301 @param bv - [out] - target bit-vector
302 @param seq_vect - [out] DNA sequence vector
303 @param k-size - dimention for k-mer generation
304 */
305template<typename BV>
307 const vector_char_type& seq_vect,
308 unsigned k_size,
309 bool check)
310{
311 const bm::id64_t chunk_size = 400000000;
312 bm::chrono_taker tt1(cout, "2. Generate k-mers", 1, &timing_map);
313
314 bv.clear();
315 bv.init(); // need to explicitly init to use bvector<>::set_bit_no_check()
316 if (seq_vect.empty())
317 return;
318 const char* dna_str = &seq_vect[0];
319
320 std::vector<bm::id64_t> k_buf;
321 k_buf.reserve(chunk_size);
322
323 {
324 bm::id64_t k_mer_code;
325 vector_char_type::size_type dna_sz = seq_vect.size()-(k_size-1);
326 vector_char_type::size_type pos = 0;
327 bool valid = false;
328 for (; pos < dna_sz; ++pos)
329 {
330 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
331 if (valid)
332 {
333 k_buf.push_back(k_mer_code);
334 break;
335 }
336 } // for pos
337
338 const unsigned k_shift = (k_size-1) * 2;
339 if (valid)
340 {
341 for (++pos; pos < dna_sz; ++pos)
342 {
343 bm::id64_t bp_code;
344 valid = get_DNA_code(dna_str[pos + (k_size - 1)], bp_code);
345 if (!valid)
346 {
347 pos += k_size; // wind fwrd to the next BP char
348 for (; pos < dna_sz; ++pos) // search for the next valid k-mer
349 {
350 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
351 if (valid)
352 {
353 k_buf.push_back(k_mer_code);
354 break;
355 }
356 }
357 continue;
358 }
359 // shift out the previous base pair code, OR the new arrival
360 k_mer_code = ((k_mer_code >> 2) | (bp_code << k_shift));
361
362 // generated k-mer codes are accumulated in buffer for sorting
363 k_buf.push_back(k_mer_code);
364
365 if (check)
366 {
367 validate_k_mer(dna_str, pos, k_size, k_mer_code);
368 bm::id64_t k_check;
369 valid = get_kmer_code(dna_str, pos, k_size, k_check);
370 assert(valid);
371 assert(k_check == k_mer_code);
372 }
373
374 if (k_buf.size() >= chunk_size) // sorting check.point
375 {
376 sort_unique(k_buf);
377 if (k_buf.size())
378 {
379 bv.set(&k_buf[0], k_buf.size(), bm::BM_SORTED); // fast bulk set
380 k_buf.resize(0);
381 bv.optimize(); // periodically re-optimize to save memory
382 }
383
384 float pcnt = float(pos) / float(dna_sz);
385 pcnt *= 100;
386 cout << "\r" << unsigned(pcnt) << "% of " << dna_sz
387 << " (" << (pos+1) <<") "
388 << flush;
389 }
390 } // for pos
391 }
392
393 if (k_buf.size()) // add last incomplete chunk here
394 {
395 sort_unique(k_buf);
396 bv.set(&k_buf[0], k_buf.size(), bm::BM_SORTED); // fast bulk set
397
398 cout << "Unique k-mers: " << k_buf.size() << endl;
399 }
400 }
401 bv.optimize();
402}
403
404/// k-mer counting algorithm using reference sequence,
405/// regenerates k-mer codes, sorts them and counts
406///
407inline
408void count_kmers(const vector_char_type& seq_vect,
409 unsigned k_size,
410 rsc_sparse_vector_u32& kmer_counts)
411{
412 const bm::id64_t chunk_size = 400000000;
413 if (seq_vect.empty())
414 return;
415 const char* dna_str = &seq_vect[0];
416 std::vector<bm::id64_t> k_buf;
417 k_buf.reserve(chunk_size);
418
419 bm::id64_t k_mer_code;
420 vector_char_type::size_type dna_sz = seq_vect.size()-(k_size-1);
421 vector_char_type::size_type pos = 0;
422 bool valid = false;
423 for (; pos < dna_sz; ++pos)
424 {
425 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
426 if (valid)
427 {
428 k_buf.push_back(k_mer_code);
429 break;
430 }
431 } // for pos
432 const unsigned k_shift = (k_size-1) * 2;
433 if (valid)
434 {
435 for (++pos; pos < dna_sz; ++pos)
436 {
437 bm::id64_t bp_code;
438 valid = get_DNA_code(dna_str[pos + (k_size - 1)], bp_code);
439 if (!valid)
440 {
441 pos += k_size; // wind fwrd to the next BP char
442 for (; pos < dna_sz; ++pos) // search for the next valid k-mer
443 {
444 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
445 if (valid)
446 {
447 k_buf.push_back(k_mer_code);
448 break;
449 }
450 }
451 continue;
452 }
453 // shift out the previous base pair code, OR the new arrival
454 k_mer_code = ((k_mer_code >> 2) | (bp_code << k_shift));
455 // generated k-mer codes are accumulated in buffer for sorting
456 k_buf.push_back(k_mer_code);
457
458 if (k_buf.size() == chunk_size) // sorting point
459 {
460 sort_count(k_buf, kmer_counts);
461 k_buf.resize(0);
462
463 float pcnt = float(pos) / float(dna_sz);
464 pcnt *= 100;
465 cout << "\r" << unsigned(pcnt) << "% of " << dna_sz
466 << " (" << (pos+1) <<") "
467 << flush;
468 }
469
470 } // for pos
471 }
472 sort_count(k_buf, kmer_counts);
473}
474
475/// Functor to process job batch (task)
476///
477template<typename BV>
479{
480public:
481 typedef BV bvector_type;
483
484 /// constructor
485 ///
487 unsigned k_size,
488 size_type from,
489 size_type to,
490 rsc_sparse_vector_u32& kmer_counts)
491 : m_seq_vect(seq_vect), m_k_size(k_size), m_from(from), m_to(to),
492 m_kmer_counts(kmer_counts)
493 {}
494
496 : m_seq_vect(func.m_seq_vect), m_k_size(func.m_k_size),
497 m_from(func.m_from), m_to(func.m_to),
498 m_kmer_counts(func.m_kmer_counts)
499 {}
500
501 /// Main logic (functor)
503 {
504 const bvector_type* bv_null = m_kmer_counts.get_null_bvector();
505 rsc_sparse_vector_u32 kmer_counts_part(*bv_null);
506 kmer_counts_part.sync();
507
508 const bm::id64_t chunk_size = 2000000;
509 if (m_seq_vect.empty())
510 return;
511
512 const char* dna_str = &m_seq_vect[0];
513 std::vector<bm::id64_t> k_buf;
514 k_buf.reserve(chunk_size);
515
516 const auto k_size = m_k_size;
517
518 bm::id64_t k_mer_code(0);
519 vector_char_type::size_type dna_sz = m_seq_vect.size()-(m_k_size-1);
520 vector_char_type::size_type pos = 0;
521 bool valid = false;
522 for (; pos < dna_sz; ++pos)
523 {
524 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
525 if (valid)
526 {
527 if (k_mer_code >= m_from && k_mer_code <= m_to)
528 k_buf.push_back(k_mer_code);
529 break;
530 }
531 } // for pos
532
533 const unsigned k_shift = (k_size-1) * 2;
534 if (valid)
535 {
536 for (++pos; pos < dna_sz; ++pos)
537 {
538 bm::id64_t bp_code;
539 valid = get_DNA_code(dna_str[pos + (k_size - 1)], bp_code);
540 if (!valid)
541 {
542 pos += k_size; // wind fwrd to the next BP char
543 for (; pos < dna_sz; ++pos) // search for the next valid k-mer
544 {
545 valid = get_kmer_code(dna_str, pos, k_size, k_mer_code);
546 if (valid)
547 {
548 if (k_mer_code >= m_from && k_mer_code <= m_to)
549 k_buf.push_back(k_mer_code);
550 break;
551 }
552 }
553 continue;
554 }
555 // shift out the previous base pair code, OR the new arrival
556 k_mer_code = ((k_mer_code >> 2) | (bp_code << k_shift));
557 // generated k-mer codes are accumulated in buffer for sorting
558 if (k_mer_code >= m_from && k_mer_code <= m_to)
559 {
560 k_buf.push_back(k_mer_code);
561 if (k_buf.size() == chunk_size) // sorting point
562 {
563 sort_count(k_buf, kmer_counts_part);
564 k_buf.resize(0);
565 }
566 }
567 } // for pos
568 }
569
570 sort_count(k_buf, kmer_counts_part);
571
572 // merge results
573 {
574 static std::mutex mtx_counts_lock;
575 std::lock_guard<std::mutex> guard(mtx_counts_lock);
576
577 // merge_not_null fast merge for non-overlapping ranges of
578 // rsc_sparse_vector
579 m_kmer_counts.merge_not_null(kmer_counts_part);
580 }
581 }
582
583private:
584 const vector_char_type& m_seq_vect;
585 unsigned m_k_size;
586 size_type m_from;
587 size_type m_to;
588 rsc_sparse_vector_u32& m_kmer_counts;
589};
590
591/// MT k-mer counting
592///
593template<typename BV>
594void count_kmers_parallel(const BV& bv_kmers,
595 const vector_char_type& seq_vect,
596 rsc_sparse_vector_u32& kmer_counts,
597 unsigned k_size,
598 unsigned concurrency)
599{
600 typedef typename BV::size_type bv_size_type;
601 typedef std::vector<std::pair<bv_size_type, bv_size_type> > bv_ranges_vector;
602
603 bv_ranges_vector pair_vect;
604
605 bv_size_type cnt = bv_kmers.count();
606 bv_size_type split_rank = cnt / concurrency; // target population count per job
607
608 if (split_rank < concurrency || concurrency < 2)
609 {
610 count_kmers(seq_vect, k_size, kmer_counts); // run single threaded
611 return;
612 }
613 // run split algorithm to determine equal weight ranges for parallel
614 // processing
615 bm::rank_range_split(bv_kmers, split_rank, pair_vect);
616
617 // Create parallel async tasks running on a range of source sequence
618 //
619 std::vector<std::future<void> > futures;
620 futures.reserve(concurrency);
621
622 for (size_t k = 0; k < pair_vect.size(); ++k)
623 {
624 futures.emplace_back(std::async(std::launch::async,
626 pair_vect[k].first,
627 pair_vect[k].second,
628 kmer_counts)));
629 } // for k
630
631 // wait for all jobs to finish, print progress report
632 //
633 for (auto& e : futures)
634 {
635 unsigned m = 0;
636 while(1)
637 {
638 std::future_status status = e.wait_for(std::chrono::seconds(60));
639 if (status == std::future_status::ready)
640 break;
641 cout << "\r" << ++m << " min" << flush;
642 } // while
643 } // for
644 cout << endl;
645}
646
647
648/// k-mer counting method using Bitap algorithm for occurence search
649/// this method is significantly slower than direct regeneration of k-mer
650/// codes and sorting count
651///
652template<typename BV>
653void count_kmers(const BV& bv_kmers, rsc_sparse_vector_u32& kmer_counts)
654{
655 std::string kmer_str;
656 typename BV::size_type cnt = 0; // progress report counter
657 typename BV::enumerator en = bv_kmers.first();
658 for ( ;en.valid(); ++en, ++cnt)
659 {
660 auto k_mer_code = *en;
661 translate_kmer(kmer_str, k_mer_code, ik_size); // translate k-mer code to string
662
663 // find list of sequence positions where k-mer is found
664 // (uncomment if you need a full search)
665 /*
666 typename BV::size_type bv_count;
667 {
668 std::vector<typename BV::size_type> km_search;
669 dna_scanner.Find(kmer_str, km_search);
670 bv_count = km_search.size();
671 }
672 */
673
674 typename BV::size_type count = dna_scanner.FindCount(kmer_str);
675 assert(count);
676 kmer_counts.set(k_mer_code, (unsigned)count);
677
678 if ((cnt % 1000) == 0)
679 cout << "\r" << cnt << flush;
680
681 } // for en
682}
683
684/**
685 k-mer counting job functor class using bm::aggregator<>
686
687 Functor received its range of k-mers in the presence-absense
688 bit-vector then follows it to run the search-counting algorithm
689 using DNA fingerprints common for all job functors.
690
691 bm::aggregator<> cannot be shared across threads,
692 so functor creates its own
693 */
694template<typename DNA_Scan>
696{
697public:
700
701 /// constructor
702 ///
703 Counting_JobFunctor(const DNA_Scan& parent_scanner,
704 const bvector_type& bv_kmers,
705 size_type from,
706 size_type to,
707 rsc_sparse_vector_u32& kmer_counts)
708 : m_parent_scanner(parent_scanner), m_bv_kmers(bv_kmers),
709 m_from(from), m_to(to), m_kmer_counts(kmer_counts)
710 {}
711
712 /// copy-ctor
714 : m_parent_scanner(func.m_parent_scanner), m_bv_kmers(func.m_bv_kmers),
715 m_from(func.m_from), m_to(func.m_to), m_kmer_counts(func.m_kmer_counts)
716 {}
717
718 /// Main logic (functor)
720 {
721 std::string kmer_str;
722 size_type cnt = 0; // progress report counter
723
724 static std::mutex mtx_counts_lock;
725 bvector_type bv_search_res;
726
727 typename bvector_type::enumerator en = m_bv_kmers.get_enumerator(m_from);
728 for ( ;en.valid(); ++en, ++cnt)
729 {
730 auto k_mer_code = *en;
731 if (k_mer_code > m_to)
732 break;
733 translate_kmer(kmer_str, k_mer_code, ik_size); // translate k-mer code to string
734
735 // setup the aggregator to perform search
736 //
737 m_Agg.reset();
738 m_Agg.set_compute_count(true); // disable full search, only count
739 for (size_t i = 0; i < kmer_str.size(); ++i)
740 {
741 const bvector_type& bv_mask = m_parent_scanner.GetVector(kmer_str[i]);
742 m_Agg.add(&bv_mask);
743 }
744
745 m_Agg.combine_shift_right_and(bv_search_res);
746
747 // Note we get search count from the Aggregator, not from search
748 // result vector, which will be empty,
749 // because we set_compute_count(true)
750 //
751 size_type search_count = m_Agg.count();
752
753 // counts are shared across threads, use locked access
754 // to save the results
755 // TODO: implement results buffering to avoid mutex overhead
756 {
757 std::lock_guard<std::mutex> guard(mtx_counts_lock);
758 m_kmer_counts.set(k_mer_code, unsigned(search_count));
759 assert(m_kmer_counts.in_sync());
760 }
761
762 k_mer_progress_count.fetch_add(1);
763
764 } // for en
765 }
766
767private:
768 const DNA_Scan& m_parent_scanner;
769 const bvector_type& m_bv_kmers;
770 size_type m_from;
771 size_type m_to;
772 typename DNA_Scan::aggregator_type m_Agg;
773
774 rsc_sparse_vector_u32& m_kmer_counts;
775};
776
777/**
778 Runs k-mer counting in parallel
779*/
780template<typename BV>
781void count_kmers_parallel(const BV& bv_kmers,
782 rsc_sparse_vector_u32& kmer_counts,
783 unsigned concurrency)
784{
785 typedef typename BV::size_type bv_size_type;
786 typedef std::vector<std::pair<bv_size_type, bv_size_type> > bv_ranges_vector;
787
788 bv_ranges_vector pair_vect;
789
790 bv_size_type cnt = bv_kmers.count();
791 bv_size_type split_rank = cnt / concurrency; // target population count per job
792
793 if (split_rank < concurrency || concurrency < 2)
794 {
795 count_kmers(bv_kmers, kmer_counts); // run single threaded
796 return;
797 }
798
799 // run split algorithm to determine equal weight ranges for parallel
800 // processing
801 bm::rank_range_split(bv_kmers, split_rank, pair_vect);
802
803 // Create parallel async tasks running on a range of source sequence
804 //
805 std::vector<std::future<void> > futures;
806 futures.reserve(concurrency);
807
808 for (size_t k = 0; k < pair_vect.size(); ++k)
809 {
810 futures.emplace_back(std::async(std::launch::async,
812 pair_vect[k].first,
813 pair_vect[k].second,
814 kmer_counts)));
815 } // for k
816
817 // wait for all jobs to finish, print progress report
818 //
819 for (auto& e : futures)
820 {
821 unsigned long long c_prev = 0;
822 while(1)
823 {
824 std::future_status status = e.wait_for(std::chrono::seconds(60));
825 if (status == std::future_status::ready)
826 break;
827
828 // progress report (entertainment)
829 //
830 unsigned long long c = k_mer_progress_count;
831 auto delta = c - c_prev;
832 c_prev = c;
833
834 auto remain_cnt = cnt - c;
835 auto remain_min = remain_cnt / delta;
836 cout << "\r" << c << ": progress per minute=" << delta;
837 if (remain_min < 120)
838 {
839 cout << " wait for " << remain_min << "m " << flush;
840 }
841 else
842 {
843 auto remain_h = remain_min / 60;
844 cout << " wait for " << remain_h << "h " << flush;
845 }
846 } // while
847 } // for
848 cout << endl;
849
850}
851
852
853/// Compute a map of how often each k-mer frequency is observed in the
854/// k-mer counts vector
855/// @param hmap - [out] histogram map
856/// @param kmer_counts - [in] kmer counts vector
857///
858static
860 const rsc_sparse_vector_u32& kmer_counts)
861{
863 kmer_counts.get_null_bvector();
864 auto en = bv_null->first();
865 for (; en.valid(); ++en)
866 {
867 auto kmer_code = *en;
868 auto kmer_count = kmer_counts.get(kmer_code);
869 auto mit = hmap.find(kmer_count);
870 if (mit == hmap.end())
871 hmap[kmer_count] = 1;
872 else
873 mit->second++;
874 } // for
875}
876
877/// Save TSV report of k-mer frequences
878/// (reverse sorted, most frequent k-mers first)
879///
880static
881void report_hmap(const string& fname, const histogram_map_u32& hmap)
882{
883 ofstream outf;
884 outf.open(fname, ios::out | ios::trunc );
885
886 outf << "kmer count \t number of kmers\n";
887
888 auto it = hmap.rbegin(); auto it_end = hmap.rend();
889 for (; it != it_end; ++it)
890 {
891 outf << it->first << "\t" << it->second << endl;
892 }
893}
894
895/// Create vector, representing subset of k-mers of high frequency
896///
897/// @param frequent_bv[out] - bit-vector of frequent k-mers (subset of all k-mers)
898/// @param hmap - histogram map of all k-mers
899/// @param kmer_counts - kmer frequency(counts) vector
900/// @param percent - percent of frequent k-mers to build a subset (5%)
901/// percent here is of total number of k-mers (not percent of all occurences)
902/// @param k_size - K mer size
903///
904template<typename BV>
905void compute_frequent_kmers(BV& frequent_bv,
906 const histogram_map_u32& hmap,
907 const rsc_sparse_vector_u32& kmer_counts,
908 unsigned percent,
909 unsigned k_size)
910{
911 (void)k_size;
912 frequent_bv.clear();
913
914 if (!percent)
915 return;
916
917 // scanner class for fast search for values in sparse vector
919 BV bv_found; // search results vector
920
922 kmer_counts.get_null_bvector();
923
924 auto total_kmers = bv_null->count();
925 bm::id64_t target_f_count = (total_kmers * percent) / 100; // how many frequent k-mers we need to pick
926
927 auto it = hmap.rbegin();
928 auto it_end = hmap.rend();
929 for (; it != it_end; ++it)
930 {
931 auto kmer_count = it->first;
932
933 scanner.find_eq(kmer_counts, kmer_count, bv_found); // seach for all values == 25
934 auto found_cnt = bv_found.count();
935 (void)found_cnt;
936 assert(found_cnt);
937 {
938 bm::bvector<>::enumerator en = bv_found.first();
939 for (; en.valid(); ++en)
940 {
941 auto kmer_code = *en;
942 unsigned k_count = kmer_counts.get(kmer_code);
943
944 if (k_count == 1) // unique k-mer, ignore
945 continue;
946
947 if (it->second == 1)
948 {
949 assert(k_count == kmer_count);
950 }
951 frequent_bv.set(kmer_code);
952 if (kmer_count >= target_f_count)
953 {
954 frequent_bv.optimize();
955 return;
956 }
957 target_f_count -= 1;
958
959 } // for en
960 }
961 } // for it
962 frequent_bv.optimize();
963}
964
965
966int main(int argc, char *argv[])
967{
968 vector_char_type seq_vect; // read FASTA sequence
969 bm::bvector<> bv_kmers; // k-mer presense(-absence) vector
970
971 try
972 {
973 auto ret = parse_args(argc, argv);
974 if (ret != 0)
975 {
976 cerr << "cmd-line parse error. " << endl;
977 return ret;
978 }
979
980 cout << "concurrency=" << parallel_jobs << endl;
981
982 if (!ifa_name.empty()) // FASTA file load
983 {
984 // limitation: loads a single molecule only
985 //
986 auto res = load_FASTA(ifa_name, seq_vect);
987 if (res != 0)
988 return res;
989 std::cout << "FASTA sequence size=" << seq_vect.size() << std::endl;
990 }
991
992 if (seq_vect.size())
993 {
994 cout << "k-mer generation for k=" << ik_size << endl;
995
996 generate_k_mer_bvector(bv_kmers, seq_vect, ik_size, is_diag);
997
998 cout << "Found " << bv_kmers.count() << " k-mers." << endl;
999
1000 if (is_diag)
1001 {
1002 bm::print_bvector_stat(cout,bv_kmers);
1003 size_t blob_size = bm::compute_serialization_size(bv_kmers);
1004 cout << "DNA 2-bit coded FASTA size=" << seq_vect.size()/4 << endl;
1005 cout << " Compressed size=" << blob_size << endl;
1006 }
1007 }
1008
1009 if (!ikd_name.empty())
1010 {
1011 bm::chrono_taker tt1(cout, "3. k-mer serialization and save", 1, &timing_map);
1012
1013 bm::SaveBVector(ikd_name.c_str(), bv_kmers);
1014 }
1015
1016 if (seq_vect.size())
1017 {
1018 bm::chrono_taker tt1(cout, "4. Build DNA fingerprints (bulk, parallel)", 1, &timing_map);
1020 }
1021
1022 if (seq_vect.size() &&
1023 (!ikd_counts_name.empty() || !ikd_rep_name.empty()))
1024 {
1025 rsc_sparse_vector_u32 rsc_kmer_counts(bv_kmers); // rank-select sparse vector for counts
1026 rsc_kmer_counts.sync();
1027 rsc_sparse_vector_u32 rsc_kmer_counts2(bv_kmers);
1028 rsc_kmer_counts2.sync();
1029
1030
1031 cout << " Searching for k-mer counts..." << endl;
1032
1033 if (is_diag) // compute reference counts using slower algorithm for verification
1034 {
1035 cout << " ... using bm::aggregator<>" << endl;
1036 bm::chrono_taker tt1(cout, "5a. k-mer counting (bm::aggregator<>)", 1, &timing_map);
1037
1038 count_kmers_parallel(bv_kmers, rsc_kmer_counts, parallel_jobs);
1039 //count_kmers(bv_kmers, rsc_kmer_counts);
1040
1041 rsc_kmer_counts.optimize();
1042 }
1043
1044 {
1045 cout << " ... using std::sort() and count" << endl;
1046 bm::chrono_taker tt1(cout, "5. k-mer counting std::sort()", 1, &timing_map);
1047
1048 count_kmers_parallel(bv_kmers, seq_vect,
1049 rsc_kmer_counts2, ik_size, parallel_jobs);
1050
1051 rsc_kmer_counts2.optimize();
1052
1053 if (!ikd_counts_name.empty())
1054 {
1055 int res = bm::file_save_svector(rsc_kmer_counts2, ikd_counts_name);
1056 if (res)
1057 {
1058 cerr << "Error: Count vector save failed!" << endl;
1059 exit(1);
1060 }
1061 }
1062
1063 if (is_diag) // verification
1064 {
1065 bool eq = rsc_kmer_counts.equal(rsc_kmer_counts2);
1066 if(!eq)
1067 {
1069 bool found = bm::sparse_vector_find_first_mismatch(rsc_kmer_counts,
1070 rsc_kmer_counts2,
1071 idx);
1072 auto v1 = rsc_kmer_counts.get(idx);
1073 auto v2 = rsc_kmer_counts2.get(idx);
1074 cerr << "Mismatch at idx=" << idx << " v1=" << v1 << " v2=" << v2 << endl;
1075 assert(found); (void)found;
1076 assert(eq);
1077 cerr << "Integrity check failed!" << endl;
1078 exit(1);
1079 }
1080 }
1081 }
1082
1083 // build histogram of k-mer counts
1084 // as a map of kmer_count -> number of occurences of k-mer
1085 //
1086
1087 histogram_map_u32 hmap;
1088 {
1089 bm::chrono_taker tt1(cout, "6. build histogram of k-mer frequencies", 1, &timing_map);
1090 compute_kmer_histogram(hmap, rsc_kmer_counts2);
1091 }
1092 if (!kh_name.empty())
1093 {
1094 report_hmap(kh_name, hmap);
1095 }
1096
1097 // here we build a build-vector of frequent (say top 5%) k-mers
1098 // (if needed we can exclude it because they likely represent repeats
1099 //
1100 bm::bvector<> bv_freq(bm::BM_GAP);
1101 {
1102 bm::chrono_taker tt1(cout, "7. Build vector of frequent k-mers", 1, &timing_map);
1103 compute_frequent_kmers(bv_freq, hmap,
1104 rsc_kmer_counts2, f_percent, ik_size);
1105 }
1106 if (!ikd_freq_name.empty())
1107 {
1108 bm::SaveBVector(ikd_freq_name.c_str(), bv_freq);
1109 }
1110
1111 cout << "Found frequent k-mers: " << bv_freq.count() << endl;
1112
1113 if (!ikd_rep_name.empty())
1114 {
1115 // exclude frequent k-mers (logical SUBtraction)
1116 //
1117 bv_kmers.bit_sub(bv_freq);
1118 bv_kmers.optimize();
1119
1120 bm::SaveBVector(ikd_rep_name.c_str(), bv_kmers);
1121 }
1122 }
1123
1124 if (is_timing)
1125 {
1126 std::cout << std::endl << "Performance:" << std::endl;
1128 }
1129
1130 }
1131 catch(std::exception& ex)
1132 {
1133 std::cerr << ex.what() << std::endl;
1134 return 1;
1135 }
1136
1137
1138
1139 return 0;
1140}
1141
Algorithms for fast aggregation of N bvectors.
Algorithms for bvector<> (main include)
Serialization / compression of bvector<>. Set theoretical operations on compressed BLOBs.
Algorithms for bm::sparse_vector.
Compressed sparse container rsc_sparse_vector<> for integer types.
Timing utilities for benchmarking (internal)
pre-processor un-defines to avoid global space pollution (internal)
k-mer counting job functor class using bm::aggregator<>
Definition: xsample07.cpp:696
DNA_Scan::bvector_type bvector_type
Definition: xsample07.cpp:698
Counting_JobFunctor(const DNA_Scan &parent_scanner, const bvector_type &bv_kmers, size_type from, size_type to, rsc_sparse_vector_u32 &kmer_counts)
constructor
Definition: xsample07.cpp:703
bvector_type::size_type size_type
Definition: xsample07.cpp:699
Counting_JobFunctor(const Counting_JobFunctor &func)
copy-ctor
Definition: xsample07.cpp:713
void operator()()
Main logic (functor)
Definition: xsample07.cpp:719
Utility for keeping all DNA finger print vectors and search using various techniques.
Definition: xsample04.cpp:181
void BuildParallel(const vector< char > &sequence, unsigned threads)
Build fingerprint bit-vectors using bulk insert iterator and parallel processing.
Definition: xsample04a.cpp:234
Functor to process job batch (task)
Definition: xsample07.cpp:479
void operator()()
Main logic (functor)
Definition: xsample07.cpp:502
SortCounting_JobFunctor(const vector_char_type &seq_vect, unsigned k_size, size_type from, size_type to, rsc_sparse_vector_u32 &kmer_counts)
constructor
Definition: xsample07.cpp:486
SortCounting_JobFunctor(const SortCounting_JobFunctor &func)
Definition: xsample07.cpp:495
bvector_type::size_type size_type
Definition: xsample07.cpp:482
Constant iterator designed to enumerate "ON" bits.
Definition: bm.h:603
bool valid() const BMNOEXCEPT
Checks if iterator is still valid.
Definition: bm.h:283
Bitvector Bit-vector container with runtime compression of bits.
Definition: bm.h:115
size_type count() const BMNOEXCEPT
population count (count of ON bits)
Definition: bm.h:2366
bm::bvector< Alloc > & bit_sub(const bm::bvector< Alloc > &bv1, const bm::bvector< Alloc > &bv2, typename bm::bvector< Alloc >::optmode opt_mode=opt_none)
3-operand SUB : this := bv1 MINUS bv2 SUBtraction is also known as AND NOT
Definition: bm.h:6092
void optimize(bm::word_t *temp_block=0, optmode opt_mode=opt_compress, statistics *stat=0)
Optimize memory bitvector's memory allocation.
Definition: bm.h:3600
bvector_size_type size_type
Definition: bm.h:121
Utility class to collect performance measurements and statistics.
Definition: bmtimer.h:41
std::map< std::string, statistics > duration_map_type
test name to duration map
Definition: bmtimer.h:66
static void print_duration_map(TOut &tout, const duration_map_type &dmap, format fmt=ct_time)
Definition: bmtimer.h:150
void merge_not_null(rsc_sparse_vector< Val, SV > &csv)
merge two vectors (argument gets destroyed) It is important that both vectors have the same NULL vect...
void set(size_type idx, value_type v)
set specified element with bounds checking and automatic resize
void optimize(bm::word_t *temp_block=0, typename bvector_type::optmode opt_mode=bvector_type::opt_compress, statistics *stat=0)
run memory optimization for all vector slices
void sync(bool force)
Re-calculate rank-select index for faster access to vector.
bool equal(const rsc_sparse_vector< Val, SV > &csv) const BMNOEXCEPT
check if another vector has the same content
value_type get(size_type idx) const BMNOEXCEPT
get specified element without bounds checking
const bvector_type * get_null_bvector() const BMNOEXCEPT
Get bit-vector of assigned values (or NULL)
bool in_sync() const BMNOEXCEPT
returns true if prefix sum table is in sync with the vector
algorithms for sparse_vector scan/search
void find_eq(const SV &sv, value_type value, bvector_type &bv_out)
find all sparse vector elements EQ to search value
succinct sparse vector with runtime compression using bit-slicing / transposition method
Definition: bmsparsevec.h:87
@ BM_SORTED
input set is sorted (ascending order)
Definition: bmconst.h:205
@ BM_GAP
GAP compression is ON.
Definition: bmconst.h:147
void rank_range_split(const BV &bv, typename BV::size_type rank, PairVect &target_v)
Algorithm to identify bit-vector ranges (splits) for the rank.
Definition: bmalgo.h:394
bool sparse_vector_find_first_mismatch(const SV &sv1, const SV &sv2, typename SV::size_type &midx, bm::null_support null_proc=bm::use_null)
Find first mismatch (element which is different) between two sparse vectors (uses linear scan in bit-...
unsigned long long int id64_t
Definition: bmconst.h:35
bm::bvector ::size_type bv_size_type
Definition: sample24.cpp:54
std::vector< std::pair< bv_size_type, bv_size_type > > bv_ranges_vector
Definition: sample24.cpp:55
static int parse_args(int argc, char *argv[])
Definition: xsample03.cpp:110
bm::aggregator< bm::bvector<> > aggregator_type
Definition: xsample04.cpp:144
std::vector< char > vector_char_type
Definition: xsample06.cpp:133
DNA_FingerprintScanner< bm::bvector<> > dna_scanner_type
Definition: xsample07.cpp:100
int main(int argc, char *argv[])
Definition: xsample07.cpp:966
void validate_k_mer(const char *dna, size_t pos, unsigned k_size, bm::id64_t k_mer)
QA function to validate if reverse k-mer decode gives the same string.
Definition: xsample07.cpp:224
void compute_frequent_kmers(BV &frequent_bv, const histogram_map_u32 &hmap, const rsc_sparse_vector_u32 &kmer_counts, unsigned percent, unsigned k_size)
Create vector, representing subset of k-mers of high frequency.
Definition: xsample07.cpp:905
void generate_k_mer_bvector(BV &bv, const vector_char_type &seq_vect, unsigned k_size, bool check)
This function turns each k-mer into an integer number and encodes it in a bit-vector (presense vector...
Definition: xsample07.cpp:306
std::string ikd_rep_name
Definition: xsample07.cpp:84
unsigned ik_size
Definition: xsample07.cpp:89
void translate_kmer(std::string &dna, bm::id64_t kmer_code, unsigned k_size)
Translate k-mer code into ATGC DNA string.
Definition: xsample07.cpp:207
static void report_hmap(const string &fname, const histogram_map_u32 &hmap)
Save TSV report of k-mer frequences (reverse sorted, most frequent k-mers first)
Definition: xsample07.cpp:881
bool is_bench
Definition: xsample07.cpp:88
std::string kh_name
Definition: xsample07.cpp:83
std::map< unsigned, unsigned > histogram_map_u32
Definition: xsample07.cpp:103
static void compute_kmer_histogram(histogram_map_u32 &hmap, const rsc_sparse_vector_u32 &kmer_counts)
Compute a map of how often each k-mer frequency is observed in the k-mer counts vector.
Definition: xsample07.cpp:859
bm::sparse_vector< unsigned, bm::bvector<> > sparse_vector_u32
Definition: xsample07.cpp:101
void count_kmers(const vector_char_type &seq_vect, unsigned k_size, rsc_sparse_vector_u32 &kmer_counts)
k-mer counting algorithm using reference sequence, regenerates k-mer codes, sorts them and counts
Definition: xsample07.cpp:408
dna_scanner_type dna_scanner
Definition: xsample07.cpp:109
std::string ikd_freq_name
Definition: xsample07.cpp:85
std::string ikd_counts_name
Definition: xsample07.cpp:82
std::string ifa_name
Definition: xsample07.cpp:80
std::string ikd_name
Definition: xsample07.cpp:81
bm::chrono_taker ::duration_map_type timing_map
Definition: xsample07.cpp:108
bool get_kmer_code(const char *dna, size_t pos, unsigned k_size, bm::id64_t &k_mer)
Calculate k-mer as an unsigned long integer.
Definition: xsample07.cpp:165
std::vector< char > vector_char_type
Definition: xsample07.cpp:99
bool is_diag
Definition: xsample07.cpp:86
bm::rsc_sparse_vector< unsigned, sparse_vector_u32 > rsc_sparse_vector_u32
Definition: xsample07.cpp:102
void sort_unique(VECT &vect)
Auxiliary function to do sort+unique on a vactor of ints removes duplicate elements.
Definition: xsample07.cpp:256
char int2DNA(unsigned code)
Translate integer code to DNA letter.
Definition: xsample07.cpp:192
void count_kmers_parallel(const BV &bv_kmers, const vector_char_type &seq_vect, rsc_sparse_vector_u32 &kmer_counts, unsigned k_size, unsigned concurrency)
MT k-mer counting.
Definition: xsample07.cpp:594
static int load_FASTA(const std::string &fname, vector_char_type &seq_vect)
really simple FASTA parser (one entry per file)
Definition: xsample07.cpp:116
unsigned f_percent
Definition: xsample07.cpp:91
void sort_count(VECT &vect, COUNT_VECT &cvect)
Auxiliary function to do sort+unique on a vactor of ints and save results in a counts vector.
Definition: xsample07.cpp:268
unsigned parallel_jobs
Definition: xsample07.cpp:90
bool get_DNA_code(char bp, bm::id64_t &dna_code)
Definition: xsample07.cpp:138
bool is_timing
Definition: xsample07.cpp:87
std::atomic_ullong k_mer_progress_count(0)
bm::bvector bvector_type
Definition: xsample07a.cpp:94