forked from ECP-ExaGraph/miniVite
-
Notifications
You must be signed in to change notification settings - Fork 0
/
graph.hpp
1215 lines (1012 loc) · 45.1 KB
/
graph.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
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// ***********************************************************************
//
// miniVite
//
// ***********************************************************************
//
// Copyright (2018) 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 GRAPH_HPP
#define GRAPH_HPP
#include <iostream>
#include <algorithm>
#include <vector>
#include <string>
#include <fstream>
#include <sstream>
#include <climits>
#include <array>
#include <unordered_map>
#include <mpi.h>
#include "utils.hpp"
unsigned seed;
struct Edge
{
GraphElem tail_;
GraphWeight weight_;
Edge(): tail_(-1), weight_(0.0) {}
};
struct EdgeTuple
{
GraphElem ij_[2];
GraphWeight w_;
EdgeTuple(GraphElem i, GraphElem j, GraphWeight w):
ij_{i, j}, w_(w)
{}
EdgeTuple(GraphElem i, GraphElem j):
ij_{i, j}, w_(1.0)
{}
EdgeTuple():
ij_{-1, -1}, w_(0.0)
{}
};
// per process graph instance
class Graph
{
public:
Graph():
lnv_(-1), lne_(-1), nv_(-1),
ne_(-1), comm_(MPI_COMM_WORLD)
{
MPI_Comm_size(comm_, &size_);
MPI_Comm_rank(comm_, &rank_);
}
Graph(GraphElem lnv, GraphElem lne,
GraphElem nv, GraphElem ne,
MPI_Comm comm=MPI_COMM_WORLD):
lnv_(lnv), lne_(lne),
nv_(nv), ne_(ne),
comm_(comm)
{
MPI_Comm_size(comm_, &size_);
MPI_Comm_rank(comm_, &rank_);
edge_indices_.resize(lnv_+1, 0);
edge_list_.resize(lne_); // this is usually populated later
parts_.resize(size_+1);
parts_[0] = 0;
for (GraphElem i = 1; i < size_+1; i++)
parts_[i] = ((nv_ * i) / size_);
}
~Graph()
{
edge_list_.clear();
edge_indices_.clear();
parts_.clear();
}
// update vertex partition information
void repart(std::vector<GraphElem> const& parts)
{ memcpy(parts_.data(), parts.data(), sizeof(GraphElem)*(size_+1)); }
// TODO FIXME put asserts like the following
// everywhere function member of Graph class
void set_edge_index(GraphElem const vertex, GraphElem const e0)
{
#if defined(DEBUG_BUILD)
assert((vertex >= 0) && (vertex <= lnv_));
assert((e0 >= 0) && (e0 <= lne_));
edge_indices_.at(vertex) = e0;
#else
edge_indices_[vertex] = e0;
#endif
}
void edge_range(GraphElem const vertex, GraphElem& e0,
GraphElem& e1) const
{
e0 = edge_indices_[vertex];
e1 = edge_indices_[vertex+1];
}
// collective
void set_nedges(GraphElem lne)
{
lne_ = lne;
edge_list_.resize(lne_);
// compute total number of edges
ne_ = 0;
MPI_Allreduce(&lne_, &ne_, 1, MPI_GRAPH_TYPE, MPI_SUM, comm_);
}
GraphElem get_base(const int rank) const
{ return parts_[rank]; }
GraphElem get_bound(const int rank) const
{ return parts_[rank+1]; }
GraphElem get_range(const int rank) const
{ return (parts_[rank+1] - parts_[rank] + 1); }
int get_owner(const GraphElem vertex) const
{
const std::vector<GraphElem>::const_iterator iter =
std::upper_bound(parts_.begin(), parts_.end(), vertex);
return (iter - parts_.begin() - 1);
}
GraphElem get_lnv() const { return lnv_; }
GraphElem get_lne() const { return lne_; }
GraphElem get_nv() const { return nv_; }
GraphElem get_ne() const { return ne_; }
MPI_Comm get_comm() const { return comm_; }
// return edge and active info
// ----------------------------
Edge const& get_edge(GraphElem const index) const
{ return edge_list_[index]; }
Edge& set_edge(GraphElem const index)
{ return edge_list_[index]; }
// local <--> global index translation
// -----------------------------------
GraphElem local_to_global(GraphElem idx)
{ return (idx + get_base(rank_)); }
GraphElem global_to_local(GraphElem idx)
{ return (idx - get_base(rank_)); }
// w.r.t passed rank
GraphElem local_to_global(GraphElem idx, int rank)
{ return (idx + get_base(rank)); }
GraphElem global_to_local(GraphElem idx, int rank)
{ return (idx - get_base(rank)); }
// print edge list (with weights)
void print(bool print_weight = true) const
{
if (lne_ < MAX_PRINT_NEDGE)
{
for (int p = 0; p < size_; p++)
{
MPI_Barrier(comm_);
if (p == rank_)
{
std::cout << "###############" << std::endl;
std::cout << "Process #" << p << ": " << std::endl;
std::cout << "###############" << std::endl;
GraphElem base = get_base(p);
for (GraphElem i = 0; i < lnv_; i++)
{
GraphElem e0, e1;
edge_range(i, e0, e1);
if (print_weight) { // print weights (default)
for (GraphElem e = e0; e < e1; e++)
{
Edge const& edge = get_edge(e);
std::cout << i+base << " " << edge.tail_ << " " << edge.weight_ << std::endl;
}
}
else { // don't print weights
for (GraphElem e = e0; e < e1; e++)
{
Edge const& edge = get_edge(e);
std::cout << i+base << " " << edge.tail_ << std::endl;
}
}
}
MPI_Barrier(comm_);
}
}
}
else
{
if (rank_ == 0)
std::cout << "Graph size per process is {" << lnv_ << ", " << lne_ <<
"}, which will overwhelm STDOUT." << std::endl;
}
}
// print statistics about edge distribution
void print_dist_stats()
{
long sumdeg = 0, maxdeg = 0;
long lne = (long) lne_;
MPI_Reduce(&lne, &sumdeg, 1, MPI_LONG, MPI_SUM, 0, comm_);
MPI_Reduce(&lne, &maxdeg, 1, MPI_LONG, MPI_MAX, 0, comm_);
long my_sq = lne*lne;
long sum_sq = 0;
MPI_Reduce(&my_sq, &sum_sq, 1, MPI_LONG, MPI_SUM, 0, comm_);
double average = (double) sumdeg / size_;
double avg_sq = (double) sum_sq / size_;
double var = avg_sq - (average*average);
double stddev = sqrt(var);
MPI_Barrier(comm_);
if (rank_ == 0)
{
std::cout << std::endl;
std::cout << "-------------------------------------------------------" << std::endl;
std::cout << "Graph edge distribution characteristics" << std::endl;
std::cout << "-------------------------------------------------------" << std::endl;
std::cout << "Number of vertices: " << nv_ << std::endl;
std::cout << "Number of edges: " << ne_ << std::endl;
std::cout << "Maximum number of edges: " << maxdeg << std::endl;
std::cout << "Average number of edges: " << average << std::endl;
std::cout << "Expected value of X^2: " << avg_sq << std::endl;
std::cout << "Variance: " << var << std::endl;
std::cout << "Standard deviation: " << stddev << std::endl;
std::cout << "-------------------------------------------------------" << std::endl;
}
}
// public variables
std::vector<GraphElem> edge_indices_;
std::vector<Edge> edge_list_;
private:
GraphElem lnv_, lne_, nv_, ne_;
std::vector<GraphElem> parts_;
MPI_Comm comm_;
int rank_, size_;
};
// read in binary edge list files
// using MPI I/O
class BinaryEdgeList
{
public:
BinaryEdgeList() :
M_(-1), N_(-1),
M_local_(-1), N_local_(-1),
comm_(MPI_COMM_WORLD)
{}
BinaryEdgeList(MPI_Comm comm) :
M_(-1), N_(-1),
M_local_(-1), N_local_(-1),
comm_(comm)
{}
// read a file and return a graph
Graph* read(int me, int nprocs, int ranks_per_node, std::string file)
{
int file_open_error;
MPI_File fh;
MPI_Status status;
// specify the number of aggregates
MPI_Info info;
MPI_Info_create(&info);
int naggr = (ranks_per_node > 1) ? (nprocs/ranks_per_node) : ranks_per_node;
if (naggr >= nprocs)
naggr = 1;
std::stringstream tmp_str;
tmp_str << naggr;
std::string str = tmp_str.str();
MPI_Info_set(info, "cb_nodes", str.c_str());
file_open_error = MPI_File_open(comm_, file.c_str(), MPI_MODE_RDONLY, info, &fh);
MPI_Info_free(&info);
if (file_open_error != MPI_SUCCESS)
{
std::cout << " Error opening file! " << std::endl;
MPI_Abort(comm_, -99);
}
// read the dimensions
MPI_File_read_all(fh, &M_, sizeof(GraphElem), MPI_BYTE, &status);
MPI_File_read_all(fh, &N_, sizeof(GraphElem), MPI_BYTE, &status);
M_local_ = ((M_*(me + 1)) / nprocs) - ((M_*me) / nprocs);
// create local graph
Graph *g = new Graph(M_local_, 0, M_, N_);
// Let N = array length and P = number of processors.
// From j = 0 to P-1,
// Starting point of array on processor j = floor(N * j / P)
// Length of array on processor j = floor(N * (j + 1) / P) - floor(N * j / P)
uint64_t tot_bytes=(M_local_+1)*sizeof(GraphElem);
MPI_Offset offset = 2*sizeof(GraphElem) + ((M_*me) / nprocs)*sizeof(GraphElem);
// read in INT_MAX increments if total byte size is > INT_MAX
if (tot_bytes < INT_MAX)
MPI_File_read_at(fh, offset, &g->edge_indices_[0], tot_bytes, MPI_BYTE, &status);
else
{
int chunk_bytes=INT_MAX;
uint8_t *curr_pointer = (uint8_t*) &g->edge_indices_[0];
uint64_t transf_bytes = 0;
while (transf_bytes < tot_bytes)
{
MPI_File_read_at(fh, offset, curr_pointer, chunk_bytes, MPI_BYTE, &status);
transf_bytes += chunk_bytes;
offset += chunk_bytes;
curr_pointer += chunk_bytes;
if ((tot_bytes - transf_bytes) < INT_MAX)
chunk_bytes = tot_bytes - transf_bytes;
}
}
N_local_ = g->edge_indices_[M_local_] - g->edge_indices_[0];
g->set_nedges(N_local_);
tot_bytes = N_local_*(sizeof(Edge));
offset = 2*sizeof(GraphElem) + (M_+1)*sizeof(GraphElem) + g->edge_indices_[0]*(sizeof(Edge));
if (tot_bytes < INT_MAX)
MPI_File_read_at(fh, offset, &g->edge_list_[0], tot_bytes, MPI_BYTE, &status);
else
{
int chunk_bytes=INT_MAX;
uint8_t *curr_pointer = (uint8_t*)&g->edge_list_[0];
uint64_t transf_bytes = 0;
while (transf_bytes < tot_bytes)
{
MPI_File_read_at(fh, offset, curr_pointer, chunk_bytes, MPI_BYTE, &status);
transf_bytes += chunk_bytes;
offset += chunk_bytes;
curr_pointer += chunk_bytes;
if ((tot_bytes - transf_bytes) < INT_MAX)
chunk_bytes = (tot_bytes - transf_bytes);
}
}
MPI_File_close(&fh);
for(GraphElem i=1; i < M_local_+1; i++)
g->edge_indices_[i] -= g->edge_indices_[0];
g->edge_indices_[0] = 0;
return g;
}
// find a distribution such that every
// process own equal number of edges (serial)
void find_balanced_num_edges(int nprocs, std::string file, std::vector<GraphElem>& mbins)
{
FILE *fp;
GraphElem nv, ne; // #vertices, #edges
std::vector<GraphElem> nbins(nprocs,0);
fp = fopen(file.c_str(), "rb");
if (fp == NULL)
{
std::cout<< " Error opening file! " << std::endl;
return;
}
// read nv and ne
fread(&nv, sizeof(GraphElem), 1, fp);
fread(&ne, sizeof(GraphElem), 1, fp);
// bin capacity
GraphElem nbcap = (ne / nprocs), ecount_idx, past_ecount_idx = 0;
int p = 0;
for (GraphElem m = 0; m < nv; m++)
{
fread(&ecount_idx, sizeof(GraphElem), 1, fp);
// bins[p] >= capacity only for the last process
if ((nbins[p] < nbcap) || (p == (nprocs - 1)))
nbins[p] += (ecount_idx - past_ecount_idx);
// increment p as long as p is not the last process
// worst case: excess edges piled up on (p-1)
if ((nbins[p] >= nbcap) && (p < (nprocs - 1)))
p++;
mbins[p+1]++;
past_ecount_idx = ecount_idx;
}
fclose(fp);
// prefix sum to store indices
for (int k = 1; k < nprocs+1; k++)
mbins[k] += mbins[k-1];
nbins.clear();
}
// read a file and return a graph
// uses a balanced distribution
// (approximately equal #edges per process)
Graph* read_balanced(int me, int nprocs, int ranks_per_node, std::string file)
{
int file_open_error;
MPI_File fh;
MPI_Status status;
std::vector<GraphElem> mbins(nprocs+1,0);
// find #vertices per process such that
// each process roughly owns equal #edges
if (me == 0)
{
find_balanced_num_edges(nprocs, file, mbins);
std::cout << "Trying to achieve equal edge distribution across processes." << std::endl;
}
MPI_Barrier(comm_);
MPI_Bcast(mbins.data(), nprocs+1, MPI_GRAPH_TYPE, 0, comm_);
// specify the number of aggregates
MPI_Info info;
MPI_Info_create(&info);
int naggr = (ranks_per_node > 1) ? (nprocs/ranks_per_node) : ranks_per_node;
if (naggr >= nprocs)
naggr = 1;
std::stringstream tmp_str;
tmp_str << naggr;
std::string str = tmp_str.str();
MPI_Info_set(info, "cb_nodes", str.c_str());
file_open_error = MPI_File_open(comm_, file.c_str(), MPI_MODE_RDONLY, info, &fh);
MPI_Info_free(&info);
if (file_open_error != MPI_SUCCESS)
{
std::cout << " Error opening file! " << std::endl;
MPI_Abort(comm_, -99);
}
// read the dimensions
MPI_File_read_all(fh, &M_, sizeof(GraphElem), MPI_BYTE, &status);
MPI_File_read_all(fh, &N_, sizeof(GraphElem), MPI_BYTE, &status);
M_local_ = mbins[me+1] - mbins[me];
// create local graph
Graph *g = new Graph(M_local_, 0, M_, N_);
// readjust parts with new vertex partition
g->repart(mbins);
uint64_t tot_bytes=(M_local_+1)*sizeof(GraphElem);
MPI_Offset offset = 2*sizeof(GraphElem) + mbins[me]*sizeof(GraphElem);
// read in INT_MAX increments if total byte size is > INT_MAX
if (tot_bytes < INT_MAX)
MPI_File_read_at(fh, offset, &g->edge_indices_[0], tot_bytes, MPI_BYTE, &status);
else
{
int chunk_bytes=INT_MAX;
uint8_t *curr_pointer = (uint8_t*) &g->edge_indices_[0];
uint64_t transf_bytes = 0;
while (transf_bytes < tot_bytes)
{
MPI_File_read_at(fh, offset, curr_pointer, chunk_bytes, MPI_BYTE, &status);
transf_bytes += chunk_bytes;
offset += chunk_bytes;
curr_pointer += chunk_bytes;
if ((tot_bytes - transf_bytes) < INT_MAX)
chunk_bytes = tot_bytes - transf_bytes;
}
}
N_local_ = g->edge_indices_[M_local_] - g->edge_indices_[0];
g->set_nedges(N_local_);
tot_bytes = N_local_*(sizeof(Edge));
offset = 2*sizeof(GraphElem) + (M_+1)*sizeof(GraphElem) + g->edge_indices_[0]*(sizeof(Edge));
if (tot_bytes < INT_MAX)
MPI_File_read_at(fh, offset, &g->edge_list_[0], tot_bytes, MPI_BYTE, &status);
else
{
int chunk_bytes=INT_MAX;
uint8_t *curr_pointer = (uint8_t*)&g->edge_list_[0];
uint64_t transf_bytes = 0;
while (transf_bytes < tot_bytes)
{
MPI_File_read_at(fh, offset, curr_pointer, chunk_bytes, MPI_BYTE, &status);
transf_bytes += chunk_bytes;
offset += chunk_bytes;
curr_pointer += chunk_bytes;
if ((tot_bytes - transf_bytes) < INT_MAX)
chunk_bytes = (tot_bytes - transf_bytes);
}
}
MPI_File_close(&fh);
for(GraphElem i=1; i < M_local_+1; i++)
g->edge_indices_[i] -= g->edge_indices_[0];
g->edge_indices_[0] = 0;
mbins.clear();
return g;
}
private:
GraphElem M_;
GraphElem N_;
GraphElem M_local_;
GraphElem N_local_;
MPI_Comm comm_;
};
// RGG graph
// 1D vertex distribution
class GenerateRGG
{
public:
GenerateRGG(GraphElem nv, MPI_Comm comm = MPI_COMM_WORLD)
{
nv_ = nv;
comm_ = comm;
MPI_Comm_rank(comm_, &rank_);
MPI_Comm_size(comm_, &nprocs_);
// neighbors
up_ = down_ = MPI_PROC_NULL;
if (nprocs_ > 1) {
if (rank_ > 0 && rank_ < (nprocs_ - 1)) {
up_ = rank_ - 1;
down_ = rank_ + 1;
}
if (rank_ == 0)
down_ = 1;
if (rank_ == (nprocs_ - 1))
up_ = rank_ - 1;
}
n_ = nv_ / nprocs_;
// check if number of nodes is divisible by #processes
if ((nv_ % nprocs_) != 0) {
if (rank_ == 0) {
std::cout << "[ERROR] Number of vertices must be perfectly divisible by number of processes." << std::endl;
std::cout << "Exiting..." << std::endl;
}
MPI_Abort(comm_, -99);
}
// check if processes are power of 2
if (!is_pwr2(nprocs_)) {
if (rank_ == 0) {
std::cout << "[ERROR] Number of processes must be a power of 2." << std::endl;
std::cout << "Exiting..." << std::endl;
}
MPI_Abort(comm_, -99);
}
// calculate r(n)
GraphWeight rc = sqrt((GraphWeight)log(nv)/(GraphWeight)(PI*nv));
GraphWeight rt = sqrt((GraphWeight)2.0736/(GraphWeight)nv);
rn_ = (rc + rt)/(GraphWeight)2.0;
assert(((GraphWeight)1.0/(GraphWeight)nprocs_) > rn_);
MPI_Barrier(comm_);
}
// create RGG and returns Graph
// TODO FIXME use OpenMP wherever possible
// use Euclidean distance as edge weight
// for random edges, choose from (0,1)
// otherwise, use unit weight throughout
Graph* generate(bool isLCG, bool unitEdgeWeight = true, GraphWeight randomEdgePercent = 0.0)
{
// Generate random coordinate points
std::vector<GraphWeight> X, Y, X_up, Y_up, X_down, Y_down;
if (isLCG)
X.resize(2*n_);
else
X.resize(n_);
Y.resize(n_);
if (up_ != MPI_PROC_NULL) {
X_up.resize(n_);
Y_up.resize(n_);
}
if (down_ != MPI_PROC_NULL) {
X_down.resize(n_);
Y_down.resize(n_);
}
// create local graph
Graph *g = new Graph(n_, 0, nv_, nv_);
// generate random number within range
// X: 0, 1
// Y: rank_*1/p, (rank_+1)*1/p,
GraphWeight rec_np = (GraphWeight)(1.0/(GraphWeight)nprocs_);
GraphWeight lo = rank_* rec_np;
GraphWeight hi = lo + rec_np;
assert(hi > lo);
// measure the time to generate random numbers
MPI_Barrier(MPI_COMM_WORLD);
double st = MPI_Wtime();
if (!isLCG) {
// set seed (declared an extern in utils)
seed = (unsigned)reseeder(1);
#if defined(PRINT_RANDOM_XY_COORD)
for (int k = 0; k < nprocs_; k++) {
if (k == rank_) {
std::cout << "Random number generated on Process#" << k << " :" << std::endl;
for (GraphElem i = 0; i < n_; i++) {
X[i] = genRandom<GraphWeight>(0.0, 1.0);
Y[i] = genRandom<GraphWeight>(lo, hi);
std::cout << "X, Y: " << X[i] << ", " << Y[i] << std::endl;
}
}
MPI_Barrier(comm_);
}
#else
for (GraphElem i = 0; i < n_; i++) {
X[i] = genRandom<GraphWeight>(0.0, 1.0);
Y[i] = genRandom<GraphWeight>(lo, hi);
}
#endif
}
else { // LCG
// X | Y
// e.g seeds: 1741, 3821
// create LCG object
// seed to generate x0
LCG xr(/*seed*/1, X.data(), 2*n_, comm_);
// generate random numbers between 0-1
xr.generate();
// rescale xr further between lo-hi
// and put the numbers in Y taking
// from X[n]
xr.rescale(Y.data(), n_, lo);
#if defined(PRINT_RANDOM_XY_COORD)
for (int k = 0; k < nprocs_; k++) {
if (k == rank_) {
std::cout << "Random number generated on Process#" << k << " :" << std::endl;
for (GraphElem i = 0; i < n_; i++) {
std::cout << "X, Y: " << X[i] << ", " << Y[i] << std::endl;
}
}
MPI_Barrier(comm_);
}
#endif
}
double et = MPI_Wtime();
double tt = et - st;
double tot_tt = 0.0;
MPI_Reduce(&tt, &tot_tt, 1, MPI_DOUBLE, MPI_SUM, 0, comm_);
if (rank_ == 0) {
double tot_avg = (tot_tt/nprocs_);
std::cout << "Average time to generate " << 2*n_
<< " random numbers using LCG (in s): "
<< tot_avg << std::endl;
}
// ghost(s)
// cross edges, each processor
// communicates with up or/and down
// neighbor only
std::vector<EdgeTuple> sendup_edges, senddn_edges;
std::vector<EdgeTuple> recvup_edges, recvdn_edges;
std::vector<EdgeTuple> edgeList;
// counts, indexing: [2] = {up - 0, down - 1}
// TODO can't we use MPI_INT
std::array<GraphElem, 2> send_sizes = {0, 0}, recv_sizes = {0, 0};
#if defined(CHECK_NUM_EDGES)
GraphElem numEdges = 0;
#endif
// local
for (GraphElem i = 0; i < n_; i++) {
for (GraphElem j = i + 1; j < n_; j++) {
// euclidean distance:
// 2D: sqrt((px-qx)^2 + (py-qy)^2)
GraphWeight dx = X[i] - X[j];
GraphWeight dy = Y[i] - Y[j];
GraphWeight ed = sqrt(dx*dx + dy*dy);
// are the two vertices within the range?
if (ed <= rn_) {
// local to global index
const GraphElem g_i = g->local_to_global(i);
const GraphElem g_j = g->local_to_global(j);
if (!unitEdgeWeight) {
edgeList.emplace_back(i, g_j, ed);
edgeList.emplace_back(j, g_i, ed);
}
else {
edgeList.emplace_back(i, g_j);
edgeList.emplace_back(j, g_i);
}
#if defined(CHECK_NUM_EDGES)
numEdges += 2;
#endif
g->edge_indices_[i+1]++;
g->edge_indices_[j+1]++;
}
}
}
MPI_Barrier(comm_);
// communicate ghost coordinates with neighbors
const int x_ndown = X_down.empty() ? 0 : n_;
const int y_ndown = Y_down.empty() ? 0 : n_;
const int x_nup = X_up.empty() ? 0 : n_;
const int y_nup = Y_up.empty() ? 0 : n_;
MPI_Sendrecv(X.data(), n_, MPI_WEIGHT_TYPE, up_, SR_X_UP_TAG,
X_down.data(), x_ndown, MPI_WEIGHT_TYPE, down_, SR_X_UP_TAG,
comm_, MPI_STATUS_IGNORE);
MPI_Sendrecv(X.data(), n_, MPI_WEIGHT_TYPE, down_, SR_X_DOWN_TAG,
X_up.data(), x_nup, MPI_WEIGHT_TYPE, up_, SR_X_DOWN_TAG,
comm_, MPI_STATUS_IGNORE);
MPI_Sendrecv(Y.data(), n_, MPI_WEIGHT_TYPE, up_, SR_Y_UP_TAG,
Y_down.data(), y_ndown, MPI_WEIGHT_TYPE, down_, SR_Y_UP_TAG,
comm_, MPI_STATUS_IGNORE);
MPI_Sendrecv(Y.data(), n_, MPI_WEIGHT_TYPE, down_, SR_Y_DOWN_TAG,
Y_up.data(), y_nup, MPI_WEIGHT_TYPE, up_, SR_Y_DOWN_TAG,
comm_, MPI_STATUS_IGNORE);
// exchange ghost vertices / cross edges
if (nprocs_ > 1) {
if (up_ != MPI_PROC_NULL) {
for (GraphElem i = 0; i < n_; i++) {
for (GraphElem j = i + 1; j < n_; j++) {
GraphWeight dx = X[i] - X_up[j];
GraphWeight dy = Y[i] - Y_up[j];
GraphWeight ed = sqrt(dx*dx + dy*dy);
if (ed <= rn_) {
const GraphElem g_i = g->local_to_global(i);
const GraphElem g_j = j + up_*n_;
if (!unitEdgeWeight) {
sendup_edges.emplace_back(j, g_i, ed);
edgeList.emplace_back(i, g_j, ed);
}
else {
sendup_edges.emplace_back(j, g_i);
edgeList.emplace_back(i, g_j);
}
#if defined(CHECK_NUM_EDGES)
numEdges++;
#endif
g->edge_indices_[i+1]++;
}
}
}
// send up sizes
send_sizes[0] = sendup_edges.size();
}
if (down_ != MPI_PROC_NULL) {
for (GraphElem i = 0; i < n_; i++) {
for (GraphElem j = i + 1; j < n_; j++) {
GraphWeight dx = X[i] - X_down[j];
GraphWeight dy = Y[i] - Y_down[j];
GraphWeight ed = sqrt(dx*dx + dy*dy);
if (ed <= rn_) {
const GraphElem g_i = g->local_to_global(i);
const GraphElem g_j = j + down_*n_;
if (!unitEdgeWeight) {
senddn_edges.emplace_back(j, g_i, ed);
edgeList.emplace_back(i, g_j, ed);
}
else {
senddn_edges.emplace_back(j, g_i);
edgeList.emplace_back(i, g_j);
}
#if defined(CHECK_NUM_EDGES)
numEdges++;
#endif
g->edge_indices_[i+1]++;
}
}
}
// send down sizes
send_sizes[1] = senddn_edges.size();
}
}
MPI_Barrier(comm_);
// communicate ghost vertices with neighbors
// send/recv buffer sizes
MPI_Sendrecv(&send_sizes[0], 1, MPI_GRAPH_TYPE, up_, SR_SIZES_UP_TAG,
&recv_sizes[1], 1, MPI_GRAPH_TYPE, down_, SR_SIZES_UP_TAG,
comm_, MPI_STATUS_IGNORE);
MPI_Sendrecv(&send_sizes[1], 1, MPI_GRAPH_TYPE, down_, SR_SIZES_DOWN_TAG,
&recv_sizes[0], 1, MPI_GRAPH_TYPE, up_, SR_SIZES_DOWN_TAG,
comm_, MPI_STATUS_IGNORE);
// resize recv buffers
if (recv_sizes[0] > 0)
recvup_edges.resize(recv_sizes[0]);
if (recv_sizes[1] > 0)
recvdn_edges.resize(recv_sizes[1]);
// send/recv both up and down
MPI_Sendrecv(sendup_edges.data(), send_sizes[0]*sizeof(struct EdgeTuple), MPI_BYTE,
up_, SR_UP_TAG, recvdn_edges.data(), recv_sizes[1]*sizeof(struct EdgeTuple),
MPI_BYTE, down_, SR_UP_TAG, comm_, MPI_STATUS_IGNORE);
MPI_Sendrecv(senddn_edges.data(), send_sizes[1]*sizeof(struct EdgeTuple), MPI_BYTE,
down_, SR_DOWN_TAG, recvup_edges.data(), recv_sizes[0]*sizeof(struct EdgeTuple),
MPI_BYTE, up_, SR_DOWN_TAG, comm_, MPI_STATUS_IGNORE);
// update local #edges
// down
if (down_ != MPI_PROC_NULL) {
for (GraphElem i = 0; i < recv_sizes[1]; i++) {
#if defined(CHECK_NUM_EDGES)
numEdges++;
#endif
if (!unitEdgeWeight)
edgeList.emplace_back(recvdn_edges[i].ij_[0], recvdn_edges[i].ij_[1], recvdn_edges[i].w_);
else
edgeList.emplace_back(recvdn_edges[i].ij_[0], recvdn_edges[i].ij_[1]);
g->edge_indices_[recvdn_edges[i].ij_[0]+1]++;
}
}
// up
if (up_ != MPI_PROC_NULL) {
for (GraphElem i = 0; i < recv_sizes[0]; i++) {
#if defined(CHECK_NUM_EDGES)
numEdges++;
#endif
if (!unitEdgeWeight)
edgeList.emplace_back(recvup_edges[i].ij_[0], recvup_edges[i].ij_[1], recvup_edges[i].w_);
else
edgeList.emplace_back(recvup_edges[i].ij_[0], recvup_edges[i].ij_[1]);
g->edge_indices_[recvup_edges[i].ij_[0]+1]++;
}
}
// add random edges based on
// randomEdgePercent
if (randomEdgePercent > 0.0) {
const GraphElem pnedges = (edgeList.size()/2);
GraphElem tot_pnedges = 0;
MPI_Allreduce(&pnedges, &tot_pnedges, 1, MPI_GRAPH_TYPE, MPI_SUM, comm_);
// extra #edges per process
const GraphElem nrande = (((GraphElem)(randomEdgePercent * (GraphWeight)tot_pnedges))/100);
GraphElem pnrande = 0.0;
// TODO FIXME try to ensure a fair edge distibution
if (nrande < nprocs_) {
if (rank_ == (nprocs_ - 1))
pnrande += nrande;
}
else {
pnrande = nrande / nprocs_;
const GraphElem pnrem = nrande % nprocs_;
if (pnrem != 0) {
if (rank_ == (nprocs_ - 1))
pnrande += pnrem;
}
}
// add pnrande edges
// send/recv buffers
std::vector<std::vector<EdgeTuple>> rand_edges(nprocs_);
std::vector<EdgeTuple> sendrand_edges, recvrand_edges;
// outgoing/incoming send/recv sizes
// TODO FIXME if number of randomly added edges are above
// INT_MAX, weird things will happen, fix it
std::vector<int> sendrand_sizes(nprocs_), recvrand_sizes(nprocs_);
#if defined(PRINT_EXTRA_NEDGES)
int extraEdges = 0;
#endif
#if defined(DEBUG_PRINTF)
for (int i = 0; i < nprocs_; i++) {
if (i == rank_) {
std::cout << "[" << i << "]Target process for random edge insertion between "
<< lo << " and " << hi << std::endl;
}
MPI_Barrier(comm_);
}
#endif
// make sure each process has a
// different seed this time since
// we want random edges
unsigned rande_seed = (unsigned)(time(0)^getpid());
GraphWeight weight = 1.0;
std::hash<GraphElem> reh;
// cannot use genRandom if it's already been seeded
std::default_random_engine re(rande_seed);
std::uniform_int_distribution<GraphElem> IR, JR;
std::uniform_real_distribution<GraphWeight> IJW;
for (GraphElem k = 0; k < pnrande; k++) {