COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
BPMaximumMatching.cpp
Go to the documentation of this file.
1
2#ifdef THREADED
3#ifndef _OPENMP
4#define _OPENMP
5#endif
6
7#include <omp.h>
8int cblas_splits = 1;
9#endif
10
11#include "../CombBLAS.h"
12#include <mpi.h>
13#include <sys/time.h>
14#include <iostream>
15#include <functional>
16#include <algorithm>
17#include <vector>
18#include <string>
19#include <sstream>
20
21
22#include "BPMaximalMatching.h"
23#include "BPMaximumMatching.h"
24
25using namespace std;
26using namespace combblas;
27
28
33
34// algorithmic options
36int init;
38bool fewexp;
41string ofname;
42
43
44
45/*
46 Remove isolated vertices and purmute
47 */
49{
50
51 int nprocs, myrank;
54
55
58 FullyDistVec<int64_t, int64_t> nonisoRowV; // id's of non-isolated (connected) Row vertices
59 FullyDistVec<int64_t, int64_t> nonisoColV; // id's of non-isolated (connected) Col vertices
60 FullyDistVec<int64_t, int64_t> nonisov; // id's of non-isolated (connected) vertices
61
62 A.Reduce(*ColSums, Column, plus<int64_t>(), static_cast<int64_t>(0));
63 A.Reduce(*RowSums, Row, plus<int64_t>(), static_cast<int64_t>(0));
64
65 // this steps for general graph
66 /*
67 ColSums->EWiseApply(*RowSums, plus<int64_t>()); not needed for bipartite graph
68 nonisov = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
69 nonisov.RandPerm(); // so that A(v,v) is load-balanced (both memory and time wise)
70 A.operator()(nonisov, nonisov, true); // in-place permute to save memory
71 */
72
73 // this steps for bipartite graph
74 nonisoColV = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
75 nonisoRowV = RowSums->FindInds(bind2nd(greater<int64_t>(), 0));
76 delete ColSums;
77 delete RowSums;
78
79
80 {
81 nonisoColV.RandPerm();
82 nonisoRowV.RandPerm();
83 }
84
85
86 int64_t nrows1=A.getnrow(), ncols1=A.getncol(), nnz1 = A.getnnz();
87 double avgDeg1 = (double) nnz1/(nrows1+ncols1);
88
89
90 A.operator()(nonisoRowV, nonisoColV, true);
91
92 int64_t nrows2=A.getnrow(), ncols2=A.getncol(), nnz2 = A.getnnz();
93 double avgDeg2 = (double) nnz2/(nrows2+ncols2);
94
95
96 if(myrank == 0)
97 {
98 cout << "ncol nrows nedges deg \n";
99 cout << nrows1 << " " << ncols1 << " " << nnz1 << " " << avgDeg1 << " \n";
100 cout << nrows2 << " " << ncols2 << " " << nnz2 << " " << avgDeg2 << " \n";
101 }
102
104
105
106}
107
108
110{
111 int myrank;
113 if(myrank == 0)
114 {
115 cout << "\n-------------- usage --------------\n";
116 cout << "Usage (random matrix): ./bpmm <er|g500|ssca> <Scale> <EDGEFACTOR> <init><diropt><prune><graft>\n";
117 cout << "Usage (input matrix): ./bpmm <input> <matrix> <init><diropt><prune><graft>\n\n";
118
119 cout << " \n-------------- meaning of arguments ----------\n";
120 cout << "** er: Erdos-Renyi, g500: Graph500 benchmark, ssca: SSCA benchmark\n";
121 cout << "** scale: matrix dimention is 2^scale\n";
122 cout << "** edgefactor: average degree of vertices\n";
123 cout << "** (optional) init : maximal matching algorithm used to initialize\n ";
124 cout << " none: noinit, greedy: greedy init , ks: Karp-Sipser, dmd: dynamic mindegree\n";
125 cout << " default: none\n";
126 cout << "** (optional) randMaximal: random parent selection in greedy/Karp-Sipser\n" ;
127 //cout << "** (optional) diropt: employ direction-optimized BFS\n" ;
128 cout << "** (optional) prune: discard trees as soon as an augmenting path is found\n" ;
129 //cout << "** (optional) graft: employ tree grafting\n" ;
130 cout << "** (optional) moreSplit: more splitting of Matrix.\n" ;
131 cout << "** (optional) randPerm: Randomly permute the matrix for load balance.\n" ;
132 cout << "** (optional) saveMatching: Save the matching vector in a file (filename: inputfile_matching.txt).\n" ;
133 cout << "(order of optional arguments does not matter)\n";
134
135
136 cout << " \n-------------- examples ----------\n";
137 cout << "Example: mpirun -np 4 ./bpmm g500 18 16" << endl;
138 cout << "Example: mpirun -np 4 ./bpmm g500 18 16 ks diropt graft" << endl;
139 cout << "Example: mpirun -np 4 ./bpmm input cage12.mtx randPerm ks diropt graft\n" << endl;
140 }
141}
142
143void GetOptions(char* argv[], int argc)
144{
145 string allArg="";
146 for(int i=0; i<argc; i++)
147 {
148 allArg += string(argv[i]);
149 }
150
151 if(allArg.find("prune")!=string::npos)
152 prune = true;
153 if(allArg.find("fewexp")!=string::npos)
154 fewexp = true;
155 if(allArg.find("moreSplit")!=string::npos)
156 moreSplit = true;
157 if(allArg.find("saveMatching")!=string::npos)
158 saveMatching=true;
159 if(allArg.find("randMM")!=string::npos)
160 randMM = true;
161 if(allArg.find("randMaximal")!=string::npos)
162 randMaximal = true;
163 if(allArg.find("randPerm")!=string::npos)
164 randPerm = true;
165 if(allArg.find("greedy")!=string::npos)
166 init = GREEDY;
167 else if(allArg.find("ks")!=string::npos)
169 else if(allArg.find("dmd")!=string::npos)
170 init = DMD;
171 else
172 init = NO_INIT;
173
174}
175
177{
179 tinfo.str("");
180 tinfo << "\n---------------------------------\n";
181 tinfo << "Calling maximum-cardinality matching with options: " << endl;
182 tinfo << " init: ";
183 if(init == NO_INIT) tinfo << " no-init ";
184 if(init == KARP_SIPSER) tinfo << " Karp-Sipser, ";
185 if(init == DMD) tinfo << " dynamic mindegree, ";
186 if(init == GREEDY) tinfo << " greedy, ";
187 if(randMaximal) tinfo << " random parent selection in greedy/Karp-Sipser, ";
188 if(prune) tinfo << " tree pruning, ";
189 if(moreSplit) tinfo << " moreSplit ";
190 if(randPerm) tinfo << " Randomly permute the matrix for load balance ";
191 if(saveMatching) tinfo << " Write the matcing in a file";
192 tinfo << "\n---------------------------------\n\n";
194
195}
196
198{
199 FullyDistVec<int64_t, int64_t> mateRow2Col ( A.getcommgrid(), A.getnrow(), (int64_t) -1);
200 FullyDistVec<int64_t, int64_t> mateCol2Row ( A.getcommgrid(), A.getncol(), (int64_t) -1);
201
202 // best option
203 init = DMD; randMaximal = false; randMM = true; prune = true;
207 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
208 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
209
210 // best option + KS
211 init = KARP_SIPSER; randMaximal = true; randMM = true; prune = true;
215 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
216 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
217
218
219 // best option + Greedy
220 init = GREEDY; randMaximal = true; randMM = true; prune = true;
224 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
225 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
226
227 // best option + No init
228 init = NO_INIT; randMaximal = false; randMM = true; prune = true;
231 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
232 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
233
234
235 // best option - randMM
236 init = DMD; randMaximal = false; randMM = false; prune = true;
240 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
241 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
242
243
244 // best option - prune
245 init = DMD; randMaximal = false; randMM = true; prune = false;
249 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
250 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
251
252}
253
254
255// default experiment
257{
258 FullyDistVec<int64_t, int64_t> mateRow2Col ( A.getcommgrid(), A.getnrow(), (int64_t) -1);
259 FullyDistVec<int64_t, int64_t> mateCol2Row ( A.getcommgrid(), A.getncol(), (int64_t) -1);
260
261 // best option
262 init = DMD; randMaximal = false; randMM = true; prune = true;
266 if(saveMatching && ofname!="")
267 {
268 mateRow2Col.ParallelWrite(ofname,false,false);
269 }
270 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
271 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
272
273}
274
275
276
277
278
279
280
282{
283
284 int nprocs, myrank;
287
288 FullyDistVec<int64_t, int64_t> mateRow2Col ( A.getcommgrid(), A.getnrow(), (int64_t) -1);
289 FullyDistVec<int64_t, int64_t> mateCol2Row ( A.getcommgrid(), A.getncol(), (int64_t) -1);
290
291 double time_start;
292
293 // best option
294 init = DMD; randMaximal = false; randMM = true; prune = true;
297 double time_dmd = MPI_Wtime()-time_start;
298 int64_t cardDMD = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
299
303 int64_t mmcardDMD = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
304 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
305 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
306
307 // best option + KS
308 init = KARP_SIPSER; randMaximal = true; randMM = true; prune = true;
311 double time_ks = MPI_Wtime()-time_start;
312 int64_t cardKS = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
313
317 int64_t mmcardKS = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
318 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
319 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
320
321
322 // best option + Greedy
323 init = GREEDY; randMaximal = true; randMM = true; prune = true;
327 int64_t cardGreedy = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
328
332 int64_t mmcardGreedy = mateRow2Col.Count([](int64_t mate){return mate!=-1;});
333 mateRow2Col.Apply([](int64_t val){return (int64_t) -1;});
334 mateCol2Row.Apply([](int64_t val){return (int64_t) -1;});
335
336 if(myrank == 0)
337 {
338 cout << "\n maximal matching experiment \n";
339 cout << cardGreedy << " " << mmcardGreedy << " " << time_greedy << " " << time_mm_greedy << " " << cardKS << " " << mmcardKS << " " << time_ks << " " << time_mm_ks << " " << cardDMD << " " << mmcardDMD << " " << time_dmd << " " << time_mm_dmd << " \n";
340 }
341}
342
343
344
345int main(int argc, char* argv[])
346{
347
348 // ------------ initialize MPI ---------------
349 int provided;
352 {
353 printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
355 }
356 int nprocs, myrank;
359 if(argc < 3)
360 {
361 ShowUsage();
362 MPI_Finalize();
363 return -1;
364 }
365
366 init = DMD;
367 randMaximal = false;
368 prune = false;
369 randMM = true;
370 moreSplit = false;
371 fewexp=false;
372 saveMatching = false;
373 ofname = "";
374 randPerm = false;
375
376 SpParHelper::Print("***** I/O and other preprocessing steps *****\n");
377 // ------------ Process input arguments and build matrix ---------------
378 {
379
382 double t01, t02;
383 if(string(argv[1]) == string("input")) // input option
384 {
385 ABool = new Par_DCSC_Bool();
386
387 string filename(argv[2]);
388 tinfo.str("");
389 tinfo << "\n**** Reading input matrix: " << filename << " ******* " << endl;
391 t01 = MPI_Wtime();
392 ABool->ParallelReadMM(filename, true, maximum<bool>()); // one-based matrix market file
393 t02 = MPI_Wtime();
394 ABool->PrintInfo();
395 tinfo.str("");
396 tinfo << "Reader took " << t02-t01 << " seconds" << endl;
398 GetOptions(argv+3, argc-3);
399 if(saveMatching)
400 {
401 ofname = filename + ".matching.out";
402 }
403
404
405 }
406 else if(argc < 4)
407 {
408 ShowUsage();
409 MPI_Finalize();
410 return -1;
411 }
412 else
413 {
414 unsigned scale = (unsigned) atoi(argv[2]);
415 unsigned EDGEFACTOR = (unsigned) atoi(argv[3]);
416 double initiator[4];
417 if(string(argv[1]) == string("er"))
418 {
419 initiator[0] = .25;
420 initiator[1] = .25;
421 initiator[2] = .25;
422 initiator[3] = .25;
423 if(myrank==0)
424 cout << "Randomly generated ER matric\n";
425 }
426 else if(string(argv[1]) == string("g500"))
427 {
428 initiator[0] = .57;
429 initiator[1] = .19;
430 initiator[2] = .19;
431 initiator[3] = .05;
432 if(myrank==0)
433 cout << "Randomly generated G500 matric\n";
434 }
435 else if(string(argv[1]) == string("ssca"))
436 {
437 initiator[0] = .6;
438 initiator[1] = .4/3;
439 initiator[2] = .4/3;
440 initiator[3] = .4/3;
441 if(myrank==0)
442 cout << "Randomly generated SSCA matric\n";
443 }
444 else
445 {
446 if(myrank == 0)
447 printf("The input type - %s - is not recognized.\n", argv[2]);
449 }
450
451 SpParHelper::Print("Generating input matrix....\n");
452 t01 = MPI_Wtime();
455 ABool = new Par_DCSC_Bool(*DEL, false);
456 delete DEL;
457 t02 = MPI_Wtime();
458 ABool->PrintInfo();
459 tinfo.str("");
460 tinfo << "Generator took " << t02-t01 << " seconds" << endl;
462
464 //removeIsolated(*ABool);
465 SpParHelper::Print("Generated matrix symmetricized....\n");
466 ABool->PrintInfo();
467
468 GetOptions(argv+4, argc-4);
469
470 }
471
472
473 if(randPerm)
474 {
475 // randomly permute for load balance
476 SpParHelper::Print("Performing random permutation of matrix.\n");
479 prow.iota(ABool->getnrow(), 0);
480 pcol.iota(ABool->getncol(), 0);
481 prow.RandPerm();
482 pcol.RandPerm();
483 (*ABool)(prow, pcol, true);
484 SpParHelper::Print("Performed random permutation of matrix.\n");
485 }
486
487
490 AT.Transpose();
491
492 // Reduce is not multithreaded, so I am doing it here
494 A.Reduce(degCol, Column, plus<int64_t>(), static_cast<int64_t>(0));
495
496 int nthreads;
497#ifdef THREADED
498#pragma omp parallel
499 {
500 int splitPerThread = 1;
502 nthreads = omp_get_num_threads();
504 }
505 tinfo.str("");
506 tinfo << "Threading activated with " << nthreads << " threads, and matrix split into "<< cblas_splits << " parts" << endl;
508 A.ActivateThreading(cblas_splits); // note: crash on empty matrix
509 AT.ActivateThreading(cblas_splits);
510#endif
511
512
513 SpParHelper::Print("**************************************************\n\n");
515 }
516 MPI_Finalize();
517 return 0;
518}
519
520
bool moreSplit
int main(int argc, char *argv[])
void experiment_maximal(Par_DCSC_Bool &A, Par_DCSC_Bool &AT, FullyDistVec< int64_t, int64_t > degCol)
bool randMaximal
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > Par_DCSC_Bool
bool fewexp
bool prune
SpParMat< int64_t, bool, SpCCols< int64_t, bool > > Par_CSC_Bool
int init
void GetOptions(char *argv[], int argc)
SpParMat< int64_t, double, SpDCCols< int64_t, double > > Par_DCSC_Double
void defaultExp(Par_DCSC_Bool &A, Par_DCSC_Bool &AT, FullyDistVec< int64_t, int64_t > degCol)
void ShowUsage()
bool randMM
string ofname
void showCurOptions()
bool saveMatching
bool randPerm
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > Par_DCSC_int64_t
void removeIsolated(Par_DCSC_Bool &A)
void experiment(Par_DCSC_Bool &A, Par_DCSC_Bool &AT, FullyDistVec< int64_t, int64_t > degCol)
#define GREEDY
#define KARP_SIPSER
#define DMD
#define NO_INIT
int cblas_splits
Definition DirOptBFS.cpp:72
#define EDGEFACTOR
Definition DirOptBFS.cpp:81
void GenGraph500Data(double initiator[4], int log_numverts, int edgefactor, bool scramble=false, bool packed=false)
static void Print(const std::string &s)
long int64_t
Definition compat.h:21
@ Column
Definition SpDefs.h:114
void Symmetricize(PARMAT &A)
Definition ReadMatDist.h:29
void maximumMatching(SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row, bool prune=true, bool randMM=false, bool maximizeWeight=false)
void MaximalMatching(Par_DCSC_Bool &A, Par_DCSC_Bool &AT, FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row, FullyDistVec< IT, IT > &degColRecv, int type, bool rand=true)
double A