1#ifndef BP_MAXIMAL_MATCHING_H
2#define BP_MAXIMAL_MATCHING_H
4#include "CombBLAS/CombBLAS.h"
22template <
typename Par_DCSC_Bool,
typename IT>
24 FullyDistVec<IT, IT>& mateCol2Row, FullyDistVec<IT, IT>& degColRecv,
int type,
bool rand=
true)
30 MPI_Comm comm =
A.getcommgrid()->GetWorld();
31 MPI_Comm_size(comm,&
nprocs);
32 MPI_Comm_rank(comm,&myrank);
37 nthreads = omp_get_num_threads();
41 FullyDistVec<IT, IT> degCol = degColRecv;
44 FullyDistSpVec<IT, IT> unmatchedRow(mateRow2Col, [](
IT mate){
return mate==-1;});
45 FullyDistSpVec<IT, IT> degColSG(
A.getcommgrid(),
A.getncol());
50 FullyDistSpVec<IT, VertexType> unmatchedCol(
A.getcommgrid(),
A.getncol());
52 unmatchedCol = EWiseApply<VertexType>(unmatchedCol, degCol,
58 FullyDistSpVec<IT, VertexType> fringeRow(
A.getcommgrid(),
A.getnrow());
59 FullyDistSpVec<IT, IT> fringeRow2(
A.getcommgrid(),
A.getnrow());
60 FullyDistSpVec<IT, VertexType> deg1Col(
A.getcommgrid(),
A.getncol());
63 IT curUnmatchedCol = unmatchedCol.getnnz();
64 IT curUnmatchedRow = unmatchedRow.getnnz();
67 double tStart = MPI_Wtime();
68 std::vector<std::vector<double> > timing;
73 cout <<
"=======================================================\n";
74 cout <<
"@@@@@@ Number of processes: " <<
nprocs << endl;
75 cout <<
"=======================================================\n";
76 cout <<
"It | UMRow | UMCol | newlyMatched | Time "<< endl;
77 cout <<
"=======================================================\n";
83 while(curUnmatchedCol !=0 && curUnmatchedRow!=0 && newlyMatched != 0 )
88 unmatchedCol = EWiseApply<VertexType>(unmatchedCol, degCol,
99 std::vector<double> times;
100 double t1 = MPI_Wtime();
103 SpMV<Select2ndMinSR<bool, VertexType>>(
A, unmatchedCol, fringeRow,
false);
107 SpMV<Select2ndMinSR<bool, VertexType>>(
A, unmatchedCol, fringeRow,
false);
111 deg1Col = EWiseApply<VertexType>(unmatchedCol, degCol,
116 if(deg1Col.getnnz()>9)
117 SpMV<Select2ndMinSR<bool, VertexType>>(
A, deg1Col, fringeRow,
false);
119 SpMV<Select2ndMinSR<bool, VertexType>>(
A, unmatchedCol, fringeRow,
false);
122 fringeRow = EWiseApply<VertexType>(fringeRow, mateRow2Col,
127 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
133 fringeRow2 = EWiseApply<IT>(fringeRow, mateRow2Col,
138 FullyDistSpVec<IT, IT> newMatchedCols = fringeRow2.Invert(
A.getncol());
139 FullyDistSpVec<IT, IT> newMatchedRows = newMatchedCols.Invert(
A.getnrow());
140 mateCol2Row.Set(newMatchedCols);
141 mateRow2Col.Set(newMatchedRows);
142 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
147 unmatchedRow.Select(mateRow2Col, [](
IT mate){
return mate==-1;});
148 unmatchedCol.Select(mateCol2Row, [](
IT mate){
return mate==-1;});
153 newMatchedRows.Apply([](
IT val){
return 1;});
154 SpMV< SelectPlusSR<bool, IT>>(AT, newMatchedRows, degColSG,
false);
156 degCol.EWiseApply(degColSG,
157 [](
IT old_deg,
IT new_deg){
return old_deg-new_deg;},
158 [](
IT old_deg,
IT new_deg){
return true;},
159 false,
static_cast<IT>(0));
161 unmatchedCol = EWiseApply<VertexType>(unmatchedCol, degCol,
166 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
171 newlyMatched = newMatchedCols.getnnz();
172 times.push_back(std::accumulate(times.begin(), times.end(), 0.0));
173 timing.push_back(times);
177 printf(
"%3d %10lld %10lld %10lld %18lf\n", iteration , curUnmatchedRow, curUnmatchedCol, newlyMatched, times.back());
180 curUnmatchedCol = unmatchedCol.getnnz();
181 curUnmatchedRow = unmatchedRow.getnnz();
186 IT cardinality = mateRow2Col.Count([](
IT mate){
return mate!=-1;});
187 std::vector<double> totalTimes(timing[0].
size(),0);
188 for(
int i=0; i<timing.size(); i++)
190 for(
int j=0; j<timing[i].size(); j++)
192 totalTimes[j] += timing[i][j];
200 cout <<
"==========================================================\n";
201 cout <<
"\n================individual timings =======================\n";
202 cout <<
" SpMV Update-Match Update-UMC Total "<< endl;
203 cout <<
"==========================================================\n";
204 for(
int i=0; i<timing.size(); i++)
206 for(
int j=0; j<timing[i].size(); j++)
208 printf(
"%12.5lf ", timing[i][j]);
213 cout <<
"-------------------------------------------------------\n";
214 for(
int i=0; i<totalTimes.size(); i++)
215 printf(
"%12.5lf ", totalTimes[i]);
218 std::cout <<
"****** maximal matching runtime ********\n";
219 std::cout <<
"nprocesses nthreads ncores algorithm Unmatched-Rows Cardinality Total Time***\n";
220 std::cout <<
nprocs <<
" " << nthreads <<
" " <<
nprocs * nthreads <<
" ";
221 if(type ==
DMD) std::cout <<
"DMD";
222 else if(type ==
GREEDY) std::cout <<
"Greedy";
223 else if(type ==
KARP_SIPSER) std::cout <<
"Karp-Sipser";
226 printf(
"%lld %lld %lf\n", curUnmatchedRow, cardinality, totalTimes.back());
227 std::cout <<
"-------------------------------------------------------\n\n";
239template <
typename Par_MAT_Double,
typename IT>
240void WeightedGreedy(Par_MAT_Double &
A, FullyDistVec<IT, IT>& mateRow2Col,
241 FullyDistVec<IT, IT>& mateCol2Row, FullyDistVec<IT, IT>& degCol)
244 typedef VertexTypeML < IT, double>
VertexType;
249 nthreads = omp_get_num_threads();
252 PreAllocatedSPA<VertexType> SPA(
A.seq(), nthreads*4);
254 MPI_Comm comm =
A.getcommgrid()->GetWorld();
255 MPI_Comm_size(comm,&
nprocs);
256 MPI_Comm_rank(comm,&myrank);
260 FullyDistSpVec<IT, IT> unmatchedRow(mateRow2Col, [](
IT mate){
return mate==-1;});
264 FullyDistSpVec<IT, VertexType> unmatchedCol(
A.getcommgrid(),
A.getncol());
266 unmatchedCol = EWiseApply<VertexType>(unmatchedCol, degCol,
272 FullyDistSpVec<IT, VertexType> fringeRow(
A.getcommgrid(),
A.getnrow());
273 FullyDistSpVec<IT, IT> fringeRow2(
A.getcommgrid(),
A.getnrow());
274 FullyDistSpVec<IT, VertexType> fringeRow3(
A.getcommgrid(),
A.getnrow());
276 IT curUnmatchedCol = unmatchedCol.getnnz();
277 IT curUnmatchedRow = unmatchedRow.getnnz();
281 std::vector<std::vector<double> > timing;
286 cout <<
"=======================================================\n";
287 cout <<
"@@@@@@ Number of processes: " <<
nprocs << endl;
288 cout <<
"=======================================================\n";
289 cout <<
"It | UMRow | UMCol | newlyMatched | Time "<< endl;
290 cout <<
"=======================================================\n";
296 while(curUnmatchedCol !=0 && curUnmatchedRow!=0 && newlyMatched != 0 )
303 std::vector<double> times;
304 double t1 = MPI_Wtime();
306 SpMV<WeightMaxMLSR<double, VertexType>>(
A, unmatchedCol, fringeRow,
false, SPA);
309 fringeRow = EWiseApply<VertexType>(fringeRow, mateRow2Col,
314 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
320 fringeRow2 = EWiseApply<IT>(fringeRow, mateRow2Col,
325 FullyDistSpVec<IT, IT> newMatchedCols = fringeRow2.Invert(
A.getncol());
326 FullyDistSpVec<IT, IT> newMatchedRows = newMatchedCols.Invert(
A.getnrow());
327 mateCol2Row.Set(newMatchedCols);
328 mateRow2Col.Set(newMatchedRows);
329 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
335 unmatchedRow.Select(mateRow2Col, [](
IT mate){
return mate==-1;});
336 unmatchedCol.Select(mateCol2Row, [](
IT mate){
return mate==-1;});
337 if(myrank == 0){times.push_back(MPI_Wtime()-t1); t1 = MPI_Wtime();}
342 newlyMatched = newMatchedCols.getnnz();
343 times.push_back(std::accumulate(times.begin(), times.end(), 0.0));
344 timing.push_back(times);
348 printf(
"%3d %10lld %10lld %10lld %18lf\n", iteration , curUnmatchedRow, curUnmatchedCol, newlyMatched, times.back());
351 curUnmatchedCol = unmatchedCol.getnnz();
352 curUnmatchedRow = unmatchedRow.getnnz();
357 IT cardinality = mateRow2Col.Count([](
IT mate){
return mate!=-1;});
358 std::vector<double> totalTimes(timing[0].
size(),0);
359 for(
int i=0; i<timing.size(); i++)
361 for(
int j=0; j<timing[i].size(); j++)
363 totalTimes[j] += timing[i][j];
371 cout <<
"==========================================================\n";
372 cout <<
"\n================individual timings =======================\n";
373 cout <<
" SpMV Update-Match Update-UMC Total "<< endl;
374 cout <<
"==========================================================\n";
375 for(
int i=0; i<timing.size(); i++)
377 for(
int j=0; j<timing[i].size(); j++)
379 printf(
"%12.5lf ", timing[i][j]);
384 cout <<
"-------------------------------------------------------\n";
385 for(
int i=0; i<totalTimes.size(); i++)
386 printf(
"%12.5lf ", totalTimes[i]);
390 std::cout <<
"****** maximal matching runtime ********\n";
391 std::cout <<
"Unmatched-Rows Cardinality Total Time***\n";
392 printf(
"%lld %lld %lf\n", curUnmatchedRow, cardinality, totalTimes.back());
393 std::cout <<
"-------------------------------------------------------\n\n";
402template <
class Par_DCSC_Bool,
class IT,
class NT>
407 MPI_Comm comm =
A.getcommgrid()->GetWorld();
408 MPI_Comm_rank(comm,&myrank);
409 FullyDistSpVec<IT, IT> fringeRow(
A.getcommgrid(),
A.getnrow());
410 FullyDistSpVec<IT, IT> fringeCol(
A.getcommgrid(),
A.getncol());
411 FullyDistSpVec<IT, IT> unmatchedRow(mateRow2Col, [](
IT mate){
return mate==-1;});
412 FullyDistSpVec<IT, IT> unmatchedCol(mateCol2Row, [](
IT mate){
return mate==-1;});
413 unmatchedRow.setNumToInd();
414 unmatchedCol.setNumToInd();
417 SpMV<Select2ndMinSR<bool, VertexType>>(
A, unmatchedCol, fringeRow,
false);
418 fringeRow =
EWiseMult(fringeRow, mateRow2Col,
true, (
IT) -1);
419 if(fringeRow.getnnz() != 0)
422 std::cout <<
"Not maximal matching!!\n";
428 SpMV<Select2ndMinSR<bool, VertexType>>(tA, unmatchedRow, fringeCol,
false);
429 fringeCol =
EWiseMult(fringeCol, mateCol2Row,
true, (
IT) -1);
430 if(fringeCol.getnnz() != 0)
433 std::cout <<
"Not maximal matching**!!\n";
void MaximalMatching(Par_DCSC_Bool &A, Par_DCSC_Bool &AT, FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row, FullyDistVec< IT, IT > °ColRecv, int type, bool rand=true)
bool isMaximalmatching(Par_DCSC_Bool &A, FullyDistVec< IT, NT > &mateRow2Col, FullyDistVec< IT, NT > &mateCol2Row)
void WeightedGreedy(Par_MAT_Double &A, FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row, FullyDistVec< IT, IT > °Col)
Dcsc< IU, typename promote_trait< NU1, NU2 >::T_promote > EWiseMult(const Dcsc< IU, NU1 > &A, const Dcsc< IU, NU2 > *B, bool exclude)