-
Notifications
You must be signed in to change notification settings - Fork 10
/
main.cpp
764 lines (666 loc) · 25.6 KB
/
main.cpp
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
// ***********************************************************************
//
// Vite: A C++ library for distributed-memory graph clustering
// using MPI+OpenMP
//
// Daniel Chavarria (daniel.chavarria@pnnl.gov)
// Antonino Tumeo (antonino.tumeo@pnnl.gov)
// Mahantesh Halappanavar (hala@pnnl.gov)
// Pacific Northwest National Laboratory
//
// Hao Lu (luhowardmark@wsu.edu)
// Sayan Ghosh (sayan.ghosh@wsu.edu)
// Ananth Kalyanaraman (ananth@eecs.wsu.edu)
// Washington State University
//
// ***********************************************************************
//
// Copyright (2017) 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.
//
// ************************************************************************
#include <sys/resource.h>
#include <sys/time.h>
#include <unistd.h>
#include <cassert>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <sstream>
#include <string>
#include <cctype>
#include <omp.h>
#include <mpi.h>
#include "rebuild.hpp"
#include "louvain.hpp"
#include "compare.hpp"
#include "utils.hpp"
std::ofstream ofs;
std::ofstream ofcks;
static std::string inputFileName, outputFileName, colorArgs;
static std::string groundTruthFileName;
static int me, nprocs;
// coloring related
static bool coloring = false;
static int maxColors = 8;
static bool vertexOrdering = false;
static bool runOnePhase = false;
static bool singleColorIter = false;
static bool generateGraph = false;
static bool justProcessGraph = false;
static GraphElem numVerticesGenGraph = 0;
static GraphWeight randomEdgePercent = 0.0;
static bool readBalanced = false;
static int ranksPerNode = 1;
static bool outputFiles = false;
static bool thresholdScaling = false;
// early termination related
static bool earlyTerm = false;
static int ETType = 0;
static GraphWeight ETDelta = 1.0;
// community comparison with ground-truth
static bool compareCommunities = false;
static bool isGroundTruthZeroBased = true;
// parse command line parameters
static void parseCommandLine(const int argc, char * const argv[]);
int main(int argc, char *argv[])
{
double t0, t1, t2, t3;
MPI_Init(&argc, &argv);
MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
MPI_Comm_rank(MPI_COMM_WORLD, &me);
parseCommandLine(argc, argv);
// read ground truth info if compareCommunities
// is turned ON, serial
// Ground Truth file is expected to be of the format generated
// by LFR-gen by Fortunato, et al.
// https://sites.google.com/site/santofortunato/inthepress2
std::vector<GraphElem> commGroundTruth;
if (me == 0 /* root */ && compareCommunities) {
t1 = MPI_Wtime();
loadGroundTruthFile(commGroundTruth, groundTruthFileName, isGroundTruthZeroBased);
t2 = MPI_Wtime();
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "Time taken to load ground truth file ("
<< groundTruthFileName << "): " << (t2-t1) << std::endl;
#else
ofs << "Time taken to load ground truth file ("
<< groundTruthFileName << "): " << (t2-t1) << std::endl;
#endif
}
MPI_Barrier(MPI_COMM_WORLD); // wait till GT file is loaded
#if defined(DONT_CREATE_DIAG_FILES)
#else
std::ostringstream oss;
std::ostringstream oss2;
oss2<< "check.out." <<me;
oss << "dat.out." << me;
ofs.open(oss.str());
ofcks.open(oss2.str());
#endif
rusage rus;
createCommunityMPIType();
createEdgeMPIType();
double td0, td1;
td0 = MPI_Wtime();
DistGraph *dg = nullptr;
GraphElem teps = 0;
// load the input data file and distribute data
if (generateGraph) {
generateInMemGraph(me, nprocs, dg, numVerticesGenGraph, randomEdgePercent, outputFileName);
}
else {
if (readBalanced)
loadDistGraphMPIIOBalanced(me, nprocs, ranksPerNode, dg, inputFileName);
else
loadDistGraphMPIIO(me, nprocs, ranksPerNode, dg, inputFileName);
}
MPI_Barrier(MPI_COMM_WORLD);
assert(dg);
#ifdef PRINT_DIST_STATS
dg->printStats();
#endif
td1 = MPI_Wtime();
getrusage(RUSAGE_SELF, &rus);
#ifdef DETAILOP
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "Time to create distributed graph: " << (td1 - td0) << std::endl;
std::cout << "Memory used: " << (static_cast<double>(rus.ru_maxrss) * 1024.0) / 1048576.0 <<std::endl;
#else
ofs << "Time to create distributed graph: " << (td1 - td0) << std::endl;
ofs << "Memory used: " << (static_cast<double>(rus.ru_maxrss) * 1024.0) / 1048576.0 <<std::endl;
#endif
#else
if(me == 0)
ofs << "Time to create distributed graph: " << (td1 - td0) << std::endl;
#endif
#ifdef USE_OMP
if(me == 0) {
int nt = 0;
#pragma omp parallel
{
nt = omp_get_max_threads();
}
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "Threads Enabled : " << nt << std::endl;
#else
ofs << "Threads Enabled : " << nt << std::endl;
#endif
}
#endif
if (justProcessGraph) {
if (me == 0)
std::cout << "Option to just process the graph is turned on, so exiting without invoking community detection." << std::endl;
return 0;
}
GraphWeight currMod = -1.0, prevMod = -1.0;
double total=(td1 - td0), ctime=0.0;
int phase = 0, short_phase = 0;
int iters = 0, tot_iters = 0;
ColorElem numColors;
ColorVector colors;
CommunityVector cvect;
GraphWeight threshold;
std::vector<GraphElem> ssizes, rsizes, svdata, rvdata;
size_t ssz = 0U, rsz = 0U;
const GraphElem nv = dg->getTotalNumVertices();
// gather all communities in root
std::vector<GraphElem> commAll, cvectAll;
// initialize map
if (me == 0 && (outputFiles || compareCommunities))
commAll.resize(nv, -1);
MPI_Barrier(MPI_COMM_WORLD);
// outermost loop
while(1) {
double ptotal = 0.0;
// if threshold-scaling is ON, then
// use larger threshold towards the
// beginning, and smaller one towards
// the end
if (thresholdScaling && !runOnePhase) {
if (short_phase > 12)
short_phase = 0;
if (short_phase >= 0 && short_phase <= 2)
threshold = 1.0E-3;
else if (short_phase >= 3 && short_phase <= 6)
threshold = 1.0E-4;
else if (short_phase >= 7 && short_phase <= 9)
threshold = 1.0E-5;
else if (short_phase >= 10 && short_phase <= 12)
threshold = 1.0E-6;
}
else
threshold = 1.0E-6;
t1 = MPI_Wtime();
// only invoke coloring for first phase when the graph is the largest
if (coloring && (phase == 0)) {
t1 = MPI_Wtime();
numColors = distColoringMultiHashMinMax(me, nprocs, *dg, colors, (maxColors/2), MAX_COVG, singleColorIter);
#if defined(DONT_CREATE_DIAG_FILES)
if (me == 0) std::cout << "Number of colors (2*nHash): " << numColors << std::endl;
#else
if (me == 0) ofs << "Number of colors (2*nHash): " << numColors << std::endl;
#endif
t0 = MPI_Wtime();
if(me == 0)
#if defined(DONT_CREATE_DIAG_FILES)
std::cout<< "Coloring Time: "<<t0-t1<<std::endl;
#else
ofs<< "Coloring Time: "<<t0-t1<<std::endl;
#endif
if (earlyTerm && ETType == 1) {
currMod = distLouvainMethodWithColoring(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, true);
}
else if (earlyTerm && ETType == 2) {
currMod = distLouvainMethodWithColoring(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, ETDelta, true);
}
else if (earlyTerm && ETType == 3) {
currMod = distLouvainMethodWithColoring(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, false);
}
else if (earlyTerm && ETType == 4) {
currMod = distLouvainMethodWithColoring(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, ETDelta, false);
}
else {
currMod = distLouvainMethodWithColoring(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters);
}
}
else if (vertexOrdering && (phase == 0)) {
t1 = MPI_Wtime();
numColors = distColoringMultiHashMinMax(me, nprocs, *dg, colors, (maxColors/2), MAX_COVG, singleColorIter);
#if defined(DONT_CREATE_DIAG_FILES)
if (me == 0) std::cout << "Number of colors (2*nHash): " << numColors << std::endl;
#else
if (me == 0) ofs << "Number of colors (2*nHash): " << numColors << std::endl;
#endif
t0 = MPI_Wtime();
if(me == 0)
#if defined(DONT_CREATE_DIAG_FILES)
std::cout<< "Coloring Time: "<<t0-t1<<std::endl;
#else
ofs<< "Coloring Time: "<<t0-t1<<std::endl;
#endif
if (earlyTerm && ETType == 1) {
currMod = distLouvainMethodVertexOrder(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, true);
}
else if (earlyTerm && ETType == 2) {
currMod = distLouvainMethodVertexOrder(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, ETDelta, true);
}
else if (earlyTerm && ETType == 3) {
currMod = distLouvainMethodVertexOrder(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, false);
}
else if (earlyTerm && ETType == 4) {
currMod = distLouvainMethodVertexOrder(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters, ETDelta, false);
}
else {
currMod = distLouvainMethodVertexOrder(me, nprocs, *dg, numColors+1, colors,
ssz, rsz, ssizes, rsizes, svdata, rvdata, cvect, currMod, threshold,
iters);
}
}
else {
if (earlyTerm && ETType == 1) {
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, threshold, iters, true);
}
else if (earlyTerm && ETType == 2) {
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, threshold, iters, ETDelta, true);
}
else if (earlyTerm && ETType == 3) {
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, threshold, iters, false);
}
else if (earlyTerm && ETType == 4) {
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, threshold, iters, ETDelta, false);
}
else {
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, threshold, iters);
}
}
t0 = MPI_Wtime();
tot_iters += iters;
ptotal += (t0-t1);
ctime += (t0-t1);
if((currMod - prevMod) > threshold) {
/// Store communities in every phase
if (outputFiles || compareCommunities) {
// gather cvect into root
gatherAllComm(0 /*root*/, me, nprocs, cvectAll, cvect);
// update community info
if (me == 0) {
#ifdef DEBUG_PRINTF
std::cout << "Updated community list per level into global array at root..." << std::endl;
#endif
// ensure clusters are contiguous
const GraphElem Nc = cvectAll.size();
std::vector<GraphElem> cvectRen(Nc, 0);
std::map<GraphElem, GraphElem> kval;
GraphElem nidx = 0;
std::copy(cvectAll.begin(), cvectAll.end(), cvectRen.begin());
std::sort(cvectRen.begin(), cvectRen.end()); // default < cmp
// create a dictionary and map to new indices
for (GraphElem n = 0; n < Nc; n++)
{
if (kval.find(cvectRen[n]) == kval.end())
{
kval.insert(std::pair<GraphElem,GraphElem>(cvectRen[n], nidx));
nidx++;
}
}
// update cvectAll
for (GraphElem c = 0; c < Nc; c++)
cvectAll[c] = kval[cvectAll[c]];
cvectRen.clear();
kval.clear();
if (phase == 0)
std::copy(cvectAll.begin(), cvectAll.end(), commAll.begin());
else {
for (GraphElem p = 0; p < commAll.size(); p++)
commAll[p] = cvectAll[commAll[p]];
}
}
MPI_Barrier(MPI_COMM_WORLD);
}
/// Create new graph and rebuild
if (!runOnePhase) {
t3 = MPI_Wtime();
distbuildNextLevelGraph(nprocs, me, dg, ssz, rsz,
ssizes, rsizes, svdata, rvdata, cvect);
t2 = MPI_Wtime();
if(me == 0) {
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "**************************" << std::endl;
std::cout<< "Rebuild Time: "<<t2-t3<<std::endl;
#else
ofs << " **************************" << std::endl;
ofs<< "Rebuild Time: "<<t2-t3<<std::endl;
#endif
}
ptotal+=(t2-t3);
}
}
else { // modularity gain is not enough, exit phase.
// if thresholdScaling is ON, and the program exits
// in between the cycle, make sure it at least runs
// ONCE with 1E-6 as threshold before exiting
GraphWeight tex_1 = 0.0, tex_2 = 0.0;
if (thresholdScaling && !runOnePhase && phase < 10) {
tex_1 = MPI_Wtime();
currMod = distLouvainMethod(me, nprocs, *dg, ssz, rsz, ssizes, rsizes,
svdata, rvdata, cvect, currMod, 1.0E-6, iters);
tex_2 = MPI_Wtime();
ptotal += (tex_2 - tex_1);
}
break;
}
if(me == 0 ) {
teps += dg->getTotalNumEdges() * tot_iters;
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "Level "<< phase << ", Modularity: " << currMod <<", Clustering time: "<<t0-t1<< ", Iterations: "
<< tot_iters << std::endl;
#else
ofs << "Level "<< phase << ", Modularity: " << currMod <<", Clustering time: "<<t0-t1<< ", Iterations: "
<< tot_iters << std::endl;
#endif
}
MPI_Barrier(MPI_COMM_WORLD);
prevMod = currMod;
currMod = -1.0;
if(me == 0) {
#if defined(DONT_CREATE_DIAG_FILES)
std::cout<< "Phase Total time: " << ptotal << std::endl;
#else
ofs<< "Phase Total time: " << ptotal << std::endl;
ofcks << "**************************" << std::endl;
#endif
}
total += ptotal;
phase++;
// Exit if only running for a single phase
if (runOnePhase) {
if (me == 0)
std::cout << "Exiting at the end of first phase." << std::endl;
break;
}
if (thresholdScaling && !runOnePhase)
short_phase++;
// high phase count, break anyway
if (phase == TERMINATION_PHASE_COUNT || (tot_iters > 10000)) {
if(me == 0) {
std::cout << "Program ran for " << TERMINATION_PHASE_COUNT
<< " phases or 10K iterations, exiting..." << std::endl;
std::cout << "HINT: Check input graph weights..." << std::endl;
std::cout << "HINT: Compile with -DDEBUG_PRINTF..." << std::endl;
}
break;
}
} // end of phases
MPI_Barrier(MPI_COMM_WORLD);
double tot_time = 0.0;
MPI_Reduce(&total, &tot_time, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
double ctot_time = 0.0;
MPI_Reduce(&ctime, &ctot_time, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
if(me == 0) {
#if defined(DONT_CREATE_DIAG_FILES)
std::cout << "-------------------------------------------------------" << std::endl;
std::cout << "Average total time (secs.): " << (tot_time/nprocs) << std::endl;
std::cout << "Average time for clustering (secs.): " << (ctot_time/nprocs) << std::endl;
std::cout << "MODS (final modularity * average clustering time): " << (currMod * (ctot_time/nprocs)) << std::endl;
std::cout << "TEPS (traversed edges across clustering iterations): " << (double)teps/(double)(ctot_time/nprocs) << std::endl;
std::cout << "-------------------------------------------------------" << std::endl;
#else
ofs << "--------------------------------------------------------------" << std::endl;
ofs<< "Average (over processes) total time (secs.): "<< (tot_time/nprocs)<<std::endl;
ofs << "Average time for clustering (secs.): " << (ctot_time/nprocs) << std::endl;
ofs << "MODS (final modularity * average clustering time): " << (currMod * (ctot_time/nprocs)) << std::endl;
ofs<< "TEPS (traversed edges acrss clustering iterations): " << (double)teps/(double)(ctot_time/nprocs) <<std::endl;
ofs << "--------------------------------------------------------------" << std::endl;
#endif
}
// dump community information in a file
if (outputFiles) {
if(me == 0) {
std::string outFileName = inputFileName;
outFileName += ".communities";
std::ofstream ofs;
ofs.open(outFileName.c_str(), std::ios::out);
if (!ofs) {
std::cerr << "Error creating community file: " << outFileName << std::endl;
exit(EXIT_FAILURE);
}
// write the data
for (GraphElem i = 0; i < nv; i++) {
ofs << commAll[i] << "\n";
}
ofs.close();
std::cout << std::endl;
std::cout << "-----------------------------------------------------" << std::endl;
std::cout << "Saved community information (" << commAll.size()
<< " vertices) on file: " << outFileName << std::endl;
std::cout << "-----------------------------------------------------" << std::endl;
}
MPI_Barrier(MPI_COMM_WORLD);
}
// compare computed communities with ground truth
if(compareCommunities) {
if(me == 0)
compare_communities(commGroundTruth, commAll);
// rest of the processes wait on barrier
MPI_Barrier(MPI_COMM_WORLD);
}
commAll.clear();
cvectAll.clear();
commGroundTruth.clear();
destroyCommunityMPIType();
destroyEdgeMPIType();
colors.clear();
delete dg;
#if defined(DONT_CREATE_DIAG_FILES)
#else
#if defined(DETAILOP)
ofs.close();
ofcks.close();
#else
if(me == 0)
ofs.close();
#endif
#endif
MPI_Finalize();
return 0;
} // main
void parseCommandLine(const int argc, char * const argv[])
{
int ret;
char *temp; // check empty values
while ((ret = getopt(argc, argv, "f:bc:od:r:t:a:ig:zpn:e:s:j")) != -1) {
switch (ret) {
case 'f':
inputFileName.assign(optarg);
break;
case 'b':
readBalanced = true;
break;
case 'c':
{
colorArgs.assign(optarg);
std::stringstream ss(colorArgs);
std::string s;
std::vector<std::string> args;
while (std::getline(ss, s, ' ')) {
args.push_back(s);
}
if (args.size() > 1) {
maxColors = std::stol(args[0]);
std::transform(args[1].begin(), args[1].end(), args[1].begin(), ::tolower);
std::istringstream is(args[1]);
is >> std::boolalpha >> singleColorIter;
}
else
maxColors = std::stol(args[0]);
}
coloring = true;
break;
case 'o':
outputFiles = true;
break;
case 'd':
{
colorArgs.assign(optarg);
std::stringstream ss(colorArgs);
std::string s;
std::vector<std::string> args;
while (std::getline(ss, s, ' ')) {
args.push_back(s);
}
if (args.size() > 1) {
maxColors = std::stol(args[0]);
std::transform(args[1].begin(), args[1].end(), args[1].begin(), ::tolower);
std::istringstream is(args[1]);
is >> std::boolalpha >> singleColorIter;
}
else
maxColors = std::stol(args[0]);
}
vertexOrdering = true;
break;
case 'r':
ranksPerNode = atoi(optarg);
break;
case 't':
earlyTerm = true;
ETType = atoi(optarg);
break;
case 'a':
ETDelta = atof(optarg);
break;
case 'i':
thresholdScaling = true;
break;
case 'g':
groundTruthFileName.assign(optarg);
compareCommunities = true;
break;
case 'z':
isGroundTruthZeroBased = false;
break;
case 'p':
runOnePhase = true;
break;
case 'n':
generateGraph = true;
numVerticesGenGraph = atol(optarg);
break;
case 'e':
randomEdgePercent = atof(optarg);
break;
case 's':
outputFileName.assign(optarg);
break;
case 'j':
justProcessGraph = true;
break;
default:
assert(0 && "Should not reach here!!");
break;
}
}
if (me == 0 && !generateGraph && !outputFileName.empty()) {
std::cout << "Passing an output file has no effect for real-world graphs." << std::endl;
}
if (me == 0 && earlyTerm && (ETType < 1 || ETType > 4)) {
std::cerr << "early-term-type parameter is between 1-4 : -t 1 OR -t 2 OR -t 3 OR -t 4" << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && earlyTerm && (ETType == 2 || ETType == 4) && (ETDelta > 1.0 || ETDelta < 0.0)) {
std::cerr << "early-term-alpha must be between 0 and 1: -t 2 -a 0.5" << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && !generateGraph && inputFileName.empty()) {
std::cerr << "Must specify a binary file name with -f" << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && compareCommunities && groundTruthFileName.empty()) {
std::cerr << "Must specify a ground-truth file name with -g" << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && (coloring && vertexOrdering)) {
std::cerr << "Cannot turn on both coloring (-c) and vertex ordering (-d). Re-run using one of the options." << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && (runOnePhase && thresholdScaling)) {
std::cerr << "Cannot turn on threshold-cycling (-i) with a single phase run (-p). Re-run using one of the options." << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (me == 0 && !generateGraph && (randomEdgePercent > 0.0)) {
std::cerr << "Must specify -n <...> for graph generation first and then -p <...> to add random edges to it." << std::endl;
MPI_Abort(MPI_COMM_WORLD, -99);
}
if (generateGraph) {
// check if number of nodes is divisible by #processes
if ((numVerticesGenGraph % nprocs) != 0) {
if (me == 0) {
std::cout << "[ERROR] For in-memory graph generation, number of vertices must be perfectly divisible by number of processes." << std::endl;
std::cout << "Exiting..." << std::endl;
}
MPI_Abort(MPI_COMM_WORLD, -99);
}
// check if processes are power of 2
if (((nprocs != 0) && !(nprocs & (nprocs - 1)))) {
}
else {
if (me == 0) {
std::cout << "[ERROR] For in-memory graph generation, number of processes must be a power of 2." << std::endl;
std::cout << "Exiting..." << std::endl;
}
MPI_Abort(MPI_COMM_WORLD, -99);
}
}
} // parseCommandLine