123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191 |
- #include <nvBowtie/bowtie2/cuda/mapping.h>
- #include <nvBowtie/bowtie2/cuda/mapping_impl.h>
- namespace nvbio {
- namespace bowtie2 {
- namespace cuda {
- __global__
- void gather_ranges_kernel(
- const uint32 count,
- const uint32 n_reads,
- const SeedHitDequeArrayDeviceView hits,
- const uint32* hit_counts_scan,
- uint64* out_ranges)
- {
- const uint32 thread_id = threadIdx.x + BLOCKDIM*blockIdx.x;
- if (thread_id >= count) return;
-
-
- const uint32 read_id = upper_bound_index( thread_id, hit_counts_scan, n_reads );
-
-
- const uint32 count_offset = read_id ? hit_counts_scan[read_id-1] : 0u;
- const uint32 range_id = thread_id - count_offset;
- const SeedHit* hits_data = hits.get_data( read_id );
- const uint2 range = hits_data[ range_id ].get_range();
-
-
- out_ranges[ thread_id ] = range.y - range.x;
- }
- void gather_ranges(
- const uint32 count,
- const uint32 n_reads,
- const SeedHitDequeArrayDeviceView hits,
- const uint32* hit_counts_scan,
- uint64* out_ranges)
- {
- const int blocks = (count + BLOCKDIM-1) / BLOCKDIM;
- gather_ranges_kernel<<<blocks, BLOCKDIM>>>( count, n_reads, hits, hit_counts_scan, out_ranges );
- }
- void map_whole_read(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- const nvbio::cuda::PingPongQueuesView<uint32> queues,
- uint8* reseed,
- SeedHitDequeArrayDeviceView hits,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_whole_read_t( read_batch, fmi, rfmi, queues, reseed, hits, params, fw, rc );
- }
- void map_exact(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- const uint32 retry,
- const nvbio::cuda::PingPongQueuesView<uint32> queues,
- uint8* reseed,
- SeedHitDequeArrayDeviceView hits,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_exact_t( read_batch, fmi, rfmi, retry, queues, reseed, hits, params, fw, rc );
- }
- void map_exact(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- SeedHitDequeArrayDeviceView hits,
- const uint2 seed_range,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_exact_t( read_batch, fmi, rfmi, hits, seed_range, params, fw, rc );
- }
- void map_approx(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- const uint32 retry,
- const nvbio::cuda::PingPongQueuesView<uint32> queues,
- uint8* reseed,
- SeedHitDequeArrayDeviceView hits,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_approx_t( read_batch, fmi, rfmi, retry, queues, reseed, hits, params, fw, rc );
- }
- void map_approx(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- SeedHitDequeArrayDeviceView hits,
- const uint2 seed_range,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_approx_t( read_batch, fmi, rfmi, hits, seed_range, params, fw, rc );
- }
- void map(
- const ReadsDef::type& read_batch,
- const FMIndexDef::type fmi,
- const FMIndexDef::type rfmi,
- const uint32 retry,
- const nvbio::cuda::PingPongQueuesView<uint32> queues,
- uint8* reseed,
- SeedHitDequeArrayDeviceView hits,
- const ParamsPOD params,
- const bool fw,
- const bool rc)
- {
- map_t( read_batch, fmi, rfmi, retry, queues, reseed, hits, params, fw, rc );
- }
- }
- }
- }
|