COMBINATORIAL_BLAS  1.6
MD.cpp
Go to the documentation of this file.
1 #define DETERMINISTIC
2 #include "CombBLAS/CombBLAS.h"
3 #include <mpi.h>
4 #include <sys/time.h>
5 #include <iostream>
6 #include <functional>
7 #include <algorithm>
8 #include <vector>
9 #include <string>
10 #include <sstream>
11 
12 
13 #define EDGEFACTOR 16
14 using namespace std;
15 using namespace combblas;
16 
17 
18 
19 template <typename PARMAT>
20 void Symmetricize(PARMAT & A)
21 {
22  // boolean addition is practically a "logical or"
23  // therefore this doesn't destruct any links
24  PARMAT AT = A;
25  AT.Transpose();
26  AT.RemoveLoops(); // not needed for boolean matrices, but no harm in keeping it
27  A += AT;
28 }
29 
30 
31 
32 struct SelectMinSR
33 {
34  typedef int64_t T_promote;
35  static T_promote id(){ return -1; };
36  static bool returnedSAID() { return false; }
37  static MPI_Op mpi_op() { return MPI_MIN; };
38 
39  static T_promote add(const T_promote & arg1, const T_promote & arg2)
40  {
41  return std::min(arg1, arg2);
42  }
43 
44  static T_promote multiply(const bool & arg1, const T_promote & arg2)
45  {
46  return arg2;
47  }
48 
49  static void axpy(bool a, const T_promote & x, T_promote & y)
50  {
51  y = std::min(y, x);
52  }
53 };
54 
55 
58 void MD(PSpMat_Int64 & A);
59 
60 
61 int main(int argc, char* argv[])
62 {
63  int nprocs, myrank;
64  MPI_Init(&argc, &argv);
65  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
66  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
67  if(argc < 3)
68  {
69  if(myrank == 0)
70  {
71  cout << "Usage: ./md <rmat|er|input> <scale|filename>" << endl;
72  cout << "Example: mpirun -np 4 ./md rmat 20" << endl;
73  cout << "Example: mpirun -np 4 ./md er 20" << endl;
74  cout << "Example: mpirun -np 4 ./md input a.mtx" << endl;
75 
76  }
77  MPI_Finalize();
78  return -1;
79  }
80  {
81  PSpMat_Bool * ABool;
82 
83  if(string(argv[1]) == string("input")) // input option
84  {
85  string filename(argv[2]);
86  ifstream inf;
87  inf.open(filename.c_str(), ios::in);
88  string header;
89  getline(inf,header);
90  bool isSymmetric = header.find("symmetric");
91  bool isUnweighted = header.find("pattern");
92  inf.close();
93 
94  ABool = new PSpMat_Bool();
95  ABool->ReadDistribute(filename, 0, isUnweighted); // unweighted
96  if(isSymmetric)
97  Symmetricize(*ABool);
98  SpParHelper::Print("Read input\n");
99  }
100  else if(string(argv[1]) == string("rmat"))
101  {
102  unsigned scale;
103  scale = static_cast<unsigned>(atoi(argv[2]));
104  double initiator[4] = {.57, .19, .19, .05};
106  DEL->GenGraph500Data(initiator, scale, EDGEFACTOR, true, false );
107  MPI_Barrier(MPI_COMM_WORLD);
108 
109  ABool = new PSpMat_Bool(*DEL, false);
110  Symmetricize(*ABool);
111  delete DEL;
112  }
113  else if(string(argv[1]) == string("er"))
114  {
115  unsigned scale;
116  scale = static_cast<unsigned>(atoi(argv[2]));
117  double initiator[4] = {.25, .25, .25, .25};
119  DEL->GenGraph500Data(initiator, scale, EDGEFACTOR, true, false );
120  MPI_Barrier(MPI_COMM_WORLD);
121 
122  ABool = new PSpMat_Bool(*DEL, false);
123  Symmetricize(*ABool);
124  delete DEL;
125  }
126  else
127  {
128  SpParHelper::Print("Unknown input option\n");
129  MPI_Finalize();
130  return -1;
131  }
132 
133  Symmetricize(*ABool);
134  PSpMat_Int64 A = *ABool;
135 
136 
137  MD(A);
138 
139  }
140  MPI_Finalize();
141  return 0;
142 }
143 
144 
145 
146 // Single source BFS
147 // Find all reachable vertices from surce via enodes
149 {
150 
154  x.SetElement(source, 1);
155  visited.SetElement(source, 1);
156  while(x.getnnz() > 0)
157  {
158  SpMV<SelectMinSR>(A, x, nx, false);
159  nx.Select(visited, [](int64_t visit){return visit==0;});
160  visited.Set(nx);
161  nx.Select(enodes, [](int64_t ev){return ev!=0;}); // newly visited enodes
162  x = nx;
163  }
164 
165  FullyDistSpVec<int64_t, int64_t> reach(visited, [](int64_t visit){return visit!=0;});
166  reach.Select(enodes, [](int64_t ev){return ev==0;});
167  reach.DelElement(source); // remove source
168  return reach;
169 }
170 
171 
172 
173 
174 
175 
176 /***************************************************************************
177  * Multisource (guided) BFS from all vertices in "sources" with SpGEMM
178  * Only Vertices in "enodes" are used in the traversal
179  * Inputs:
180  * sources: the sources of BFSs. Let nnz(sources) = k
181  * A: nxn adjacency matrix (n is the number of vertices)
182  * enodes: enodes[i]=1 if the ith vertex has already received an order.
183  * let the number of enodes is r
184  ****************************************************************************/
186 {
187 
188  // ---------------------------------------------------------------------------
189  // create an nxk initial frontier from sources
190  // ith column represent the current frontier of BFS tree rooted at the ith source vertex
191  //----------------------------------------------------------------------------
192  FullyDistVec<int64_t, int64_t> ri = sources.FindInds([](int64_t x){return true;});
194  ci.iota(sources.getnnz(), 0);
195 
196  PSpMat_Int64 fringe(A.getnrow(), sources.getnnz(), ri, ci, (int64_t) 1, false);
197  PSpMat_Int64 visited(A.getnrow(), sources.getnnz(), ri, ci, (int64_t) 1, false);
199 
200 
201  // ---------------------------------------------------------------------------
202  // create an nxn matrix E such that E(i,i) = 1 if i is an enode
203  // Note: we can avoid creating this matrix from scratch
204  // by incrementally adding nonzero entries in this matrix
205  //----------------------------------------------------------------------------
206  FullyDistVec<int64_t, int64_t> ri1 = enodes.FindInds([](int64_t x){return x>0;});
207  PSpMat_Int64 E(A.getnrow(), A.getnrow(), ri1, ri1, 1);
208 
209 
210 
211  while( fringe.getnnz() > 0 )
212  {
213  // get the next frontier
214  fringe = PSpGEMM<SelectMinSR>(A, fringe);
215 
216  // remove previously visited vertices
217  fringe = EWiseMult(fringe, visited, true);
218  visited += fringe;
219 
220  // keep enodes in the frontier
221  // this can be replaced by EWiseMult if we keep a matrix with repeated enodes (memory requirement can be high)
222  fringe = PSpGEMM<PTDD>(E, fringe);
223 
224  }
225 
226 
227  FullyDistVec<int64_t, int64_t> degrees(A.getcommgrid(), sources.getnnz(), 0);
228 
229  // create an nxn matrix NE such that E(i,i) = 1 if i is not an enode
230  // Note: NE and E together create a diagonal matrix
231  FullyDistVec<int64_t, int64_t> ri2 = enodes.FindInds([](int64_t x){return x==0;});
232  PSpMat_Int64 NE(A.getnrow(), A.getnrow(), ri2, ri2, 1);
233 
234  // keep only visited non enodes
235  visited = PSpGEMM<PTDD>(NE, visited);
236 
237  // count visited non enodes from each sources
238  visited.Reduce(degrees, Column, plus<int64_t>(), static_cast<int64_t>(0));
239  degrees.Apply([](int64_t val){return (val-1);}); // -1 to remove sources themselves
240 
241  // option 2, using maskedReduce
242  /*
243  FullyDistSpVec<int64_t, int64_t> nenodes(enodes,[](int64_t x){return x==0;});
244  visited.MaskedReduce(degrees, nenodes, Column, plus<int64_t>(), static_cast<int64_t>(0));
245 
246  // or use negative mask
247  //FullyDistSpVec<int64_t, int64_t> nenodes(enodes,[](int64_t x){return x>0;});
248  //visited.MaskedReduce(degrees1, nenodes, Column, plus<int64_t>(), static_cast<int64_t>(0), true);
249  degrees.Apply([](int64_t val){return (val-1);}); // -1 to remove sources themselves
250  */
251 
252  return FullyDistSpVec<int64_t, int64_t>(sources.TotalLength(), ri, degrees);
253 }
254 
255 
256 
257 
258 
259 // assume that source is an enode
261 {
262  int nprocs = sources.getcommgrid()->GetSize();
263  int myrank = sources.getcommgrid()->GetRank();
264 
265  FullyDistSpVec<int64_t, int64_t> degrees = sources;
266  vector<int64_t> locvals = sources.GetLocalInd();
267  int64_t j = 0;
268 
269  for(int i=0; i<nprocs; )
270  {
271  int64_t s = -1;
272  if(myrank==i && j<sources.getlocnnz())
273  {
274  s = locvals[j++] + sources.LengthUntil();
275  }
276  MPI_Bcast(&s, 1, MPIType<int64_t>(), i, sources.getcommgrid()->GetWorld());
277  if(s!=-1)
278  {
279  FullyDistSpVec<int64_t, int64_t> reach = getReach(s, A, enodes);
280  degrees.SetElement(s, reach.getnnz());
281  }
282  else i++;
283  }
284  return degrees;
285 }
286 
287 
288 
290 {
294  A.Reduce(degrees, Column, plus<int64_t>(), static_cast<int64_t>(0));
295  degrees.Apply([](int64_t x){return x-1;}); // magic
296 
297  FullyDistVec<int64_t, double> treach (A.getcommgrid(), A.getnrow(), (double) 0);
298  FullyDistVec<int64_t, double> treaches (A.getcommgrid(), A.getnrow(), (double) 0);
300 
301  int myrank, nprocs;
302  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
303  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
304  double time_beg = MPI_Wtime();
305 
306  //degrees.DebugPrint();
307  double time1=0, time2=0, time3=0;
308  for(int64_t i=0; i<A.getnrow(); i++)
309  {
310  //degrees.DebugPrint();
311  int64_t s = degrees.MinElement().first; // minimum degree vertex
312  enodes.SetElement(s, i+1);
313  mdOrder.SetElement(i, s+1);
314 
315  time1 = MPI_Wtime();
316  FullyDistSpVec<int64_t, int64_t> reach = getReach(s, A, enodes);
317  time2 += MPI_Wtime()-time1;
318 
319 
320  time1 = MPI_Wtime();
322  //updatedDeg = getReachesSPMV(reach, A, enodes);
323  updatedDeg = getReachesSPMM(reach, A, enodes);
324 
325 
326  time3 += MPI_Wtime()-time1;
327 
328  degrees.Set(updatedDeg);
329  degrees.SetElement(s, A.getnrow()); // set degree to infinite
330  //degrees.DebugPrint();
331 
332 
333  int nnz = reach.getnnz();
334  if(myrank==0)
335  {
336  if(i%20==0)
337  {
338  cout << i << " .................. " << nnz << " :: " << time2 << " + " << time3 << endl;
339  time2 = 0; time3 = 0;
340  }
341 
342  }
343 
344  }
345 
346 
347  double time_end = MPI_Wtime();
348 
349 
350  if(myrank==0)
351  cout << " Total time: " << time_end - time_beg << endl;
352 
353 
354 
355  //mdOrder.DebugPrint();
356 
357 
358 
359 }
360 
361 
int64_t T_promote
Definition: MD.cpp:34
void SetElement(IT indx, NT numx)
Indexing is performed 0-based.
FullyDistVec< IT, NT > Reduce(Dim dim, _BinaryOperation __binary_op, NT id, _UnaryOperation __unary_op) const
Definition: SpParMat.cpp:791
std::shared_ptr< CommGrid > getcommgrid() const
Definition: SpParMat.h:275
void SetElement(IT indx, NT numx)
static MPI_Op mpi_op()
Definition: MD.cpp:37
MPI_Datatype MPIType< int64_t >(void)
Definition: MPIType.cpp:64
void GenGraph500Data(double initiator[4], int log_numverts, int edgefactor, bool scramble=false, bool packed=false)
FullyDistSpVec< int64_t, int64_t > getReach(int64_t source, PSpMat_Int64 &A, FullyDistVec< int64_t, int64_t > &enodes)
Definition: MD.cpp:148
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
Return the indices where pred is true.
IT getnnz() const
Definition: SpParMat.cpp:676
Dcsc< IU, typename promote_trait< NU1, NU2 >::T_promote > EWiseMult(const Dcsc< IU, NU1 > &A, const Dcsc< IU, NU2 > *B, bool exclude)
Definition: Friends.h:694
#define EDGEFACTOR
Definition: MD.cpp:13
std::shared_ptr< CommGrid > getcommgrid() const
static T_promote add(const T_promote &arg1, const T_promote &arg2)
Definition: MD.cpp:39
FullyDistSpVec< int64_t, int64_t > getReachesSPMM(FullyDistSpVec< int64_t, int64_t > &sources, PSpMat_Int64 &A, FullyDistVec< int64_t, int64_t > &enodes)
Definition: MD.cpp:185
void Select(const FullyDistVec< IT, NT1 > &denseVec, _UnaryOperation unop)
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > PSpMat_Int64
Definition: MD.cpp:57
void ReadDistribute(const std::string &filename, int master, bool nonum, HANDLER handler, bool transpose=false, bool pario=false)
Definition: SpParMat.cpp:3560
static void axpy(bool a, const T_promote &x, T_promote &y)
Definition: MD.cpp:49
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
double A
static T_promote multiply(const bool &arg1, const T_promote &arg2)
Definition: MD.cpp:44
static void Print(const std::string &s)
void iota(IT globalsize, NT first)
long int64_t
Definition: compat.h:21
IT getncol() const
Definition: SpParMat.cpp:694
void MD(PSpMat_Int64 &A)
Definition: MD.cpp:289
int main(int argc, char *argv[])
Definition: MD.cpp:61
static T_promote id()
Definition: MD.cpp:35
Definition: CCGrid.h:4
void Symmetricize(PARMAT &A)
Definition: ReadMatDist.h:29
std::vector< IT > GetLocalInd()
FullyDistSpVec< int64_t, int64_t > getReachesSPMV(FullyDistSpVec< int64_t, int64_t > &sources, PSpMat_Int64 &A, FullyDistVec< int64_t, int64_t > &enodes)
Definition: MD.cpp:260
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > PSpMat_Bool
Definition: MD.cpp:56
static bool returnedSAID()
Definition: MD.cpp:36
IT getnrow() const
Definition: SpParMat.cpp:685