-
Notifications
You must be signed in to change notification settings - Fork 3
/
hashfastric.hpp
854 lines (707 loc) · 21.8 KB
/
hashfastric.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
// ***********************************************************************
//
// TRIC
//
// ***********************************************************************
//
// Copyright (2019) Battelle Memorial Institute
// All rights reserved.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions
// are met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the copyright holder nor the names of its
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
// "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
// LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
// FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
// COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
// INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
// BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
// LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
// CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
// LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
// ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
// POSSIBILITY OF SUCH DAMAGE.
//
// ************************************************************************
#pragma once
#ifndef HASH_TFASTRIC_HPP
#define HASH_TFASTRIC_HPP
#include "graph.hpp"
#include <numeric>
#include <utility>
#include <cstring>
#include <iomanip>
#include <limits>
#ifndef TAG_DATA
#define TAG_DATA 100
#endif
#ifndef BLOOMFILTER_TOL
#define BLOOMFILTER_TOL 1E-09
#endif
#include "murmurhash/MurmurHash3.h"
class Bloomfilter
{
public:
Bloomfilter(GraphElem n, GraphWeight p=BLOOMFILTER_TOL)
: n_(pow(2, std::ceil(log(n)/log(2)))), p_(p)
{
m_ = std::ceil((n_ * log(p_)) / log(1 / pow(2, log(2))));
k_ = std::round((m_ / n_) * log(2));
if (k_%2 != 0)
k_ += 1;
hashes_.resize(k_);
bits_.resize(m_);
std::fill(bits_.begin(), bits_.end(), '0');
if (k_ == 0)
throw std::invalid_argument("Bloomfilter could not be initialized: k must be larger than 0");
}
Bloomfilter(GraphElem n, GraphElem k, GraphWeight p)
: n_(pow(2, std::ceil(log(n)/log(2)))), k_(k), p_(p)
{
m_ = std::ceil((n_ * log(p_)) / log(1 / pow(2, log(2))));
if (k_%2 != 0)
k_ += 1;
hashes_.resize(k_);
bits_.resize(m_);
std::fill(bits_.begin(), bits_.end(), '0');
if (k_ == 0)
throw std::invalid_argument("Bloomfilter could not be initialized: k must be larger than 0");
}
void insert(GraphElem const& i, GraphElem const& j)
{
hash(i, j);
for (GraphElem k = 0; k < k_; k++)
bits_[hashes_[k]] = '1';
}
void print() const
{
std::cout << "-------------Bloom filter statistics-------------" << std::endl;
std::cout << "Maximum number of Items (n): " << n_ << std::endl;
std::cout << "Probability of False Positives (p): " << p_ << std::endl;
std::cout << "Number of bits in filter (m): " << m_ << std::endl;
std::cout << "Number of hash functions (k): " << k_ << std::endl;
std::cout << "-------------------------------------------------" << std::endl;
}
void clear()
{
bits_.clear();
hashes_.clear();
}
bool contains(GraphElem i, GraphElem j)
{
hash(i, j);
for (GraphElem k = 0; k < k_; k++)
{
if (bits_[hashes_[k]] == '0')
return false;
}
return true;
}
GraphElem nbits() const
{ return m_; }
const char* data() const
{ return bits_.data(); }
char* data()
{ return bits_.data(); }
// "nucular" options, use iff
// you know what you're doing
void copy_from(char* dest)
{ std::memcpy(dest, bits_.data(), m_); }
void copy_to(char* source)
{ std::memcpy(bits_.data(), source, m_); }
void zfill()
{ std::fill(bits_.begin(), bits_.end(), '0'); }
private:
GraphElem n_, m_, k_;
GraphWeight p_;
void hash( uint64_t lhs, uint64_t rhs )
{
uint64_t key[2] = {lhs, rhs};
for (uint64_t n = 0; n < k_; n+=2)
{
MurmurHash3_x64_128 ( &key, 2*sizeof(uint64_t), 0, &hashes_[n] );
hashes_[n] = hashes_[n] % m_;
hashes_[n+1] = hashes_[n+1] % m_;
}
}
std::vector<char> bits_;
std::vector<uint64_t> hashes_;
};
class MapVec
{
public:
MapVec(): data_() {}
~MapVec()
{ data_.clear(); }
inline void insert(int key, int value)
{
if (data_.count(key) > 0)
data_[key].emplace_back(value);
else
{
data_.emplace(key, std::vector<int>());
data_[key].emplace_back(value);
}
}
inline bool contains(int key, int value)
{
if (std::find_if(data_[key].begin(), data_[key].end(),
[value](int const& element){ return element == value;}) == data_[key].end())
return false;
return true;
}
inline void clear()
{
for (auto it = data_.begin(); it != data_.end(); ++it)
it->second.clear();
data_.clear();
}
GraphElem size() const
{ return data_.size(); }
private:
std::unordered_map<int,std::vector<int>> data_;
};
class TriangulateAggrBufferedHash
{
public:
TriangulateAggrBufferedHash(Graph* g, const GraphElem bufsize):
g_(g), sbuf_ctr_(nullptr), sbuf_(nullptr), rbuf_(nullptr), pdegree_(0),
sreq_(nullptr), erange_(nullptr), vcount_(nullptr), ntriangles_(0), nghosts_(0),
out_nghosts_(0), in_nghosts_(0), pindex_(0), prev_m_(nullptr), prev_k_(nullptr),
stat_(nullptr), ebf_(nullptr), targets_(0), bufsize_(0)
#if defined(STORE_PG_INFO)
, pbf_(nullptr)
#endif
{
comm_ = g_->get_comm();
MPI_Comm_size(comm_, &size_);
MPI_Comm_rank(comm_, &rank_);
const GraphElem lnv = g_->get_lnv();
const GraphElem nv = g_->get_nv();
vcount_ = new GraphElem[lnv]();
erange_ = new GraphElem[nv*2]();
std::vector<int> rdispls(size_, 0), sdispls(size_, 0), scounts(size_, 0), rcounts(size_, 0);
std::vector<int> source_counts(size_, 0);
GraphElem *send_count = new GraphElem[size_]();
GraphElem *recv_count = new GraphElem[size_]();
double t0 = MPI_Wtime();
// store edge ranges
GraphElem base = g_->get_base(rank_);
for (GraphElem i = 0; i < lnv; i++)
{
GraphElem e0, e1;
g_->edge_range(i, e0, e1);
if ((e0 + 1) == e1)
continue;
Edge const& edge_s = g_->get_edge(e0);
Edge const& edge_t = g_->get_edge(e1-1);
erange_[(i + base)*2] = edge_s.tail_;
erange_[(i + base)*2+1] = edge_t.tail_;
}
MPI_Barrier(comm_);
MPI_Allreduce(MPI_IN_PLACE, erange_, nv*2, MPI_GRAPH_TYPE,
MPI_SUM, comm_);
GraphElem nedges = 0;
// setup bloom filter and perform local counting
for (GraphElem i = 0; i < lnv; i++)
{
GraphElem e0, e1, tup[2];
g_->edge_range(i, e0, e1);
if ((e0 + 1) == e1)
continue;
for (GraphElem m = e0; m < e1; m++)
{
Edge const& edge_m = g_->get_edge(m);
const int owner = g_->get_owner(edge_m.tail_);
if (owner != rank_)
{
if (std::find(targets_.begin(), targets_.end(), owner)
== targets_.end())
targets_.push_back(owner);
nedges += 1;
}
else
{
if (m < (e1 - 1))
{
tup[0] = edge_m.tail_;
for (GraphElem n = m + 1; n < e1; n++)
{
Edge const& edge_n = g_->get_edge(n);
tup[1] = edge_n.tail_;
if (check_edgelist(tup))
ntriangles_ += 1;
}
}
GraphElem l0, l1;
bool is_remote = false;
const GraphElem lv = g_->global_to_local(edge_m.tail_);
g_->edge_range(lv, l0, l1);
for (GraphElem l = l0; l < l1; l++)
{
Edge const& edge = g_->get_edge(l);
if (g_->get_owner(edge.tail_) != rank_)
{
is_remote = true;
break;
}
}
if (is_remote)
nedges += 1;
}
}
}
MPI_Barrier(comm_);
double t1 = MPI_Wtime();
double p_tot = t1 - t0, t_tot = 0.0;
MPI_Reduce(&p_tot, &t_tot, 1, MPI_DOUBLE, MPI_SUM, 0, comm_);
if (rank_ == 0)
{
std::cout << "Average time for local counting and misc. during instantiation (secs.): "
<< ((double)(t_tot / (double)size_)) << std::endl;
}
pdegree_ = targets_.size();
for (int i = 0; i < pdegree_; i++)
pindex_.insert({targets_[i], i});
t0 = MPI_Wtime();
MPI_Barrier(comm_);
const int targets_size = targets_.size();
MPI_Allgather(&targets_size, 1, MPI_INT, source_counts.data(), 1, MPI_INT, comm_);
int sdisp = 0, rdisp = 0;
for (int p = 0; p < size_; p++)
{
rdispls[p] = rdisp;
rdisp += source_counts[p];
}
std::vector<int> source_data(rdisp, 0);
MPI_Allgatherv(targets_.data(), targets_size, MPI_INT, source_data.data(),
source_counts.data(), rdispls.data(), MPI_INT, comm_);
if (nedges)
ebf_ = new Bloomfilter(nedges*2);
//ebf_ = new Bloomfilter(nedges*2, 10, 1.0E-5);
#if defined(STORE_PG_INFO)
#if defined(USE_BLOOMF_PG)
pbf_ = static_cast<Bloomfilter*>(new Bloomfilter(rdisp*2, 8, 1.0E-8));
#else
pbf_ = static_cast<MapVec*>(new MapVec());
#endif
#if 0
for (int p = 0; p < size_; p++)
{
if (source_counts[p] > 0)
{
for (int n = 0; n < source_counts[p]; n++)
pbf_->insert(p, source_data[n + rdispls[p]]);
}
}
#endif
for (int p = 0; p < pdegree_; p++)
for (int n = 0; n < source_counts[targets_[p]]; n++)
pbf_->insert(targets_[p], source_data[n + rdispls[targets_[p]]]);
#endif
for (GraphElem i = 0; i < lnv; i++)
{
GraphElem e0, e1;
g_->edge_range(i, e0, e1);
if ((e0 + 1) == e1)
continue;
for (GraphElem m = e0; m < e1; m++)
{
Edge const& edge_m = g_->get_edge(m);
const int owner = g_->get_owner(edge_m.tail_);
if (owner != rank_)
{
ebf_->insert(g_->local_to_global(i), edge_m.tail_);
if (m < (e1 - 1))
{
for (GraphElem n = m + 1; n < e1; n++)
{
Edge const& edge_n = g_->get_edge(n);
if (!edge_within_max(edge_m.tail_, edge_n.tail_))
break;
if (!edge_above_min(edge_m.tail_, edge_n.tail_) || !edge_above_min(edge_n.tail_, edge_m.tail_))
continue;
#if defined(STORE_PG_INFO)
if (!is_connected_pes(owner, g_->get_owner(edge_n.tail_)))
continue;
#endif
send_count[owner] += 1;
vcount_[i] += 1;
}
}
}
else
{
GraphElem l0, l1;
bool is_remote = false;
const GraphElem lv = g_->global_to_local(edge_m.tail_);
g_->edge_range(lv, l0, l1);
for (GraphElem l = l0; l < l1; l++)
{
Edge const& edge = g_->get_edge(l);
if (g_->get_owner(edge.tail_) != rank_)
{
is_remote = true;
break;
}
}
if (is_remote)
ebf_->insert(g_->local_to_global(i), edge_m.tail_);
}
}
}
if (rank_ == 0)
ebf_->print();
// outgoing/incoming data and buffer size
MPI_Alltoall(send_count, 1, MPI_GRAPH_TYPE, recv_count, 1, MPI_GRAPH_TYPE, comm_);
for (GraphElem p = 0; p < size_; p++)
{
out_nghosts_ += send_count[p];
in_nghosts_ += recv_count[p];
}
nghosts_ = out_nghosts_ + in_nghosts_;
bufsize_ = ((nghosts_*3) < bufsize) ? (nghosts_*3) : bufsize;
MPI_Allreduce(MPI_IN_PLACE, &bufsize_, 1, MPI_GRAPH_TYPE, MPI_MAX, comm_);
if (rank_ == 0)
std::cout << "Adjusted Per-PE buffer count: " << bufsize_ << std::endl;
rbuf_ = new GraphElem[bufsize_];
sbuf_ = new GraphElem[pdegree_*bufsize_];
sbuf_ctr_ = new GraphElem[pdegree_]();
prev_k_ = new GraphElem[pdegree_];
prev_m_ = new GraphElem[pdegree_];
stat_ = new char[pdegree_];
sreq_ = new MPI_Request[pdegree_];
std::fill(sreq_, sreq_ + pdegree_, MPI_REQUEST_NULL);
std::fill(prev_k_, prev_k_ + pdegree_, -1);
std::fill(prev_m_, prev_m_ + pdegree_, -1);
std::fill(stat_, stat_ + pdegree_, '0');
MPI_Barrier(comm_);
#if defined(DEBUG_PRINTF)
if (rank_ == 0)
{
std::cout << "Edge range per vertex (#ID: <range>): " << std::endl;
for (int i = 0, j = 0; i < nv*2; i+=2, j++)
std::cout << j << ": " << erange_[i] << ", " << erange_[i+1] << std::endl;
}
#endif
delete []send_count;
delete []recv_count;
sdispls.clear();
rdispls.clear();
scounts.clear();
rcounts.clear();
sdispls.clear();
rdispls.clear();
source_counts.clear();
source_data.clear();
}
~TriangulateAggrBufferedHash() {}
void clear()
{
delete []sbuf_;
delete []rbuf_;
delete []sbuf_ctr_;
delete []sreq_;
delete []prev_k_;
delete []prev_m_;
delete []stat_;
delete []vcount_;
delete []erange_;
#if defined(STORE_PG_INFO)
if (pbf_)
{
pbf_->clear();
delete pbf_;
}
#endif
if (ebf_)
{
ebf_->clear();
delete ebf_;
}
pindex_.clear();
targets_.clear();
}
// x/y are processes
#if defined(STORE_PG_INFO)
inline bool is_connected_pes(GraphElem const& x, GraphElem const& y)
{
if ((x == y) || pbf_->contains(x,y) || pbf_->contains(y,x))
return true;
return false;
}
#endif
void nbsend(GraphElem owner)
{
if (sbuf_ctr_[pindex_[owner]] > 0)
{
MPI_Isend(&sbuf_[pindex_[owner]*bufsize_], sbuf_ctr_[pindex_[owner]],
MPI_GRAPH_TYPE, owner, TAG_DATA, comm_, &sreq_[pindex_[owner]]);
}
}
void nbsend()
{
for (int const& p : targets_)
nbsend(p);
}
inline void lookup_edges()
{
const GraphElem lnv = g_->get_lnv();
for (GraphElem i = 0; i < lnv; i++)
{
if (vcount_[i] == 0) // all edges processed, move on
continue;
GraphElem e0, e1;
g_->edge_range(i, e0, e1);
if ((e0 + 1) == e1)
continue;
for (GraphElem m = e0; m < e1-1; m++)
{
EdgeStat& edge = g_->get_edge_stat(m);
const int owner = g_->get_owner(edge.edge_->tail_);
const GraphElem pidx = pindex_[owner];
const GraphElem disp = pidx*bufsize_;
if (owner != rank_ && edge.active_)
{
if (stat_[pidx] == '1')
continue;
if (m >= prev_m_[pidx])
{
if (sbuf_ctr_[pidx] == (bufsize_-1))
{
prev_m_[pidx] = m;
prev_k_[pidx] = -1;
stat_[pidx] = '1'; // messages in-flight
nbsend(owner);
continue;
}
sbuf_[disp+sbuf_ctr_[pidx]] = edge.edge_->tail_;
sbuf_ctr_[pidx] += 1;
for (GraphElem n = ((prev_k_[pidx] == -1) ? (m + 1) : prev_k_[pidx]); n < e1; n++)
{
Edge const& edge_n = g_->get_edge(n);
if (!edge_within_max(edge.edge_->tail_, edge_n.tail_))
break;
if (!edge_above_min(edge.edge_->tail_, edge_n.tail_) || !edge_above_min(edge_n.tail_, edge.edge_->tail_))
continue;
#if defined(STORE_PG_INFO)
if (!is_connected_pes(owner, g_->get_owner(edge_n.tail_)))
continue;
#endif
if (sbuf_ctr_[pidx] == (bufsize_-1))
{
prev_m_[pidx] = m;
prev_k_[pidx] = n;
sbuf_[disp+sbuf_ctr_[pidx]] = -1; // demarcate vertex boundary
sbuf_ctr_[pidx] += 1;
stat_[pidx] = '1';
nbsend(owner);
break;
}
sbuf_[disp+sbuf_ctr_[pidx]] = edge_n.tail_;
sbuf_ctr_[pidx] += 1;
out_nghosts_ -= 1;
vcount_[i] -= 1;
}
if (stat_[pidx] == '0')
{
prev_m_[pidx] = m;
prev_k_[pidx] = -1;
edge.active_ = false;
if (sbuf_ctr_[pidx] == (bufsize_-1))
{
sbuf_[disp+sbuf_ctr_[pidx]] = -1;
sbuf_ctr_[pidx] += 1;
stat_[pidx] = '1';
nbsend(owner);
}
else
{
sbuf_[disp+sbuf_ctr_[pidx]] = -1;
sbuf_ctr_[pidx] += 1;
}
}
}
}
}
}
}
inline bool check_edgelist(GraphElem tup[2])
{
GraphElem e0, e1;
const GraphElem lv = g_->global_to_local(tup[0]);
g_->edge_range(lv, e0, e1);
for (GraphElem e = e0; e < e1; e++)
{
Edge const& edge = g_->get_edge(e);
if (tup[1] == edge.tail_)
return true;
if (edge.tail_ > tup[1])
break;
}
return false;
}
inline bool edge_between_range(GraphElem x, GraphElem y) const
{
if ((y >= erange_[x*2]) && (y <= erange_[x*2+1]))
return true;
return false;
}
inline bool edge_above_min(GraphElem x, GraphElem y) const
{
if (y >= erange_[x*2])
return true;
return false;
}
inline bool edge_within_max(GraphElem x, GraphElem y) const
{
if (y <= erange_[x*2+1])
return true;
return false;
}
inline void process_messages()
{
MPI_Status status;
int flag = -1;
GraphElem tup[2] = {-1,-1}, k = 0, prev = 0;
int count = 0, source = -1;
MPI_Iprobe(MPI_ANY_SOURCE, TAG_DATA, comm_, &flag, &status);
if (flag)
{
source = status.MPI_SOURCE;
MPI_Get_count(&status, MPI_GRAPH_TYPE, &count);
MPI_Recv(rbuf_, count, MPI_GRAPH_TYPE, source,
TAG_DATA, comm_, MPI_STATUS_IGNORE);
}
else
return;
while(1)
{
if (k == count)
break;
if (rbuf_[k] == -1)
{
k += 1;
prev = k;
continue;
}
GraphElem curr_count = 0;
tup[0] = rbuf_[k];
for (GraphElem m = k + 1; m < count; m++)
{
if (rbuf_[m] == -1)
{
curr_count = m + 1;
break;
}
tup[1] = rbuf_[m];
if (ebf_->contains(tup[0], tup[1]))
ntriangles_ += 1;
in_nghosts_ -= 1;
}
k += (curr_count - prev);
prev = k;
}
}
inline GraphElem count()
{
#if defined(USE_ALLREDUCE_FOR_EXIT)
GraphElem count;
#else
bool done = false, nbar_active = false;
MPI_Request nbar_req = MPI_REQUEST_NULL;
#endif
bool sends_done = false;
int *inds = new int[pdegree_];
int over = -1;
#if defined(USE_ALLREDUCE_FOR_EXIT)
while(1)
#else
while(!done)
#endif
{
if (out_nghosts_ == 0)
{
if (!sends_done)
{
nbsend();
sends_done = true;
}
}
else
lookup_edges();
process_messages();
MPI_Testsome(pdegree_, sreq_, &over, inds, MPI_STATUSES_IGNORE);
if (over != MPI_UNDEFINED)
{
for (int i = 0; i < over; i++)
{
GraphElem idx = static_cast<GraphElem>(inds[i]);
sbuf_ctr_[idx] = 0;
stat_[idx] = '0';
}
}
#if defined(USE_ALLREDUCE_FOR_EXIT)
count = in_nghosts_;
MPI_Allreduce(MPI_IN_PLACE, &count, 1, MPI_GRAPH_TYPE, MPI_SUM, comm_);
if (count == 0)
break;
#else
if (nbar_active)
{
int test_nbar = -1;
MPI_Test(&nbar_req, &test_nbar, MPI_STATUS_IGNORE);
done = !test_nbar ? false : true;
}
else
{
if (in_nghosts_ == 0)
{
MPI_Ibarrier(comm_, &nbar_req);
nbar_active = true;
}
}
#endif
#if defined(DEBUG_PRINTF)
std::cout << "in/out: " << in_nghosts_ << ", " << out_nghosts_ << std::endl;
#endif
}
GraphElem ttc = 0, ltc = ntriangles_;
MPI_Barrier(comm_);
MPI_Reduce(<c, &ttc, 1, MPI_GRAPH_TYPE, MPI_SUM, 0, comm_);
free(inds);
return (ttc / 3);
}
private:
Graph* g_;
GraphElem ntriangles_, bufsize_, nghosts_, out_nghosts_, in_nghosts_, pdegree_;
GraphElem *sbuf_, *rbuf_, *prev_k_, *prev_m_, *sbuf_ctr_, *vcount_, *erange_;
MPI_Request *sreq_;
char *stat_;
Bloomfilter *ebf_;
#if defined(STORE_PG_INFO)
#if defined(USE_BLOOMF_PG)
Bloomfilter *pbf_;
#else
MapVec *pbf_;
#endif
#endif
std::vector<int> targets_;
int rank_, size_;
std::unordered_map<int, int> pindex_;
MPI_Comm comm_;
};
#endif