COMBINATORIAL_BLAS  1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SpAsgnTiming.cpp
Go to the documentation of this file.
1 #include <mpi.h>
2 #include <sys/time.h>
3 #include <iostream>
4 #include <functional>
5 #include <algorithm>
6 #include <vector>
7 #include <sstream>
8 #include "../CombBLAS.h"
9 
10 using namespace std;
11 
12 #define ITERATIONS 10
13 #define EDGEFACTOR 8
14 
15 int main(int argc, char* argv[])
16 {
17  int nprocs, myrank;
18  MPI_Init(&argc, &argv);
19  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
20  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
21 
22  if(argc < 2)
23  {
24  if(myrank == 0)
25  {
26  cout << "Usage: ./IndexingTiming <Scale>" << endl;
27  }
28  MPI_Finalize();
29  return -1;
30  }
31  {
33  PARDBMAT *A, *B; // declare objects
34  double initiator[4] = {.6, .4/3, .4/3, .4/3};
36 
37  int scale = static_cast<unsigned>(atoi(argv[1]));
38  ostringstream outs, outs2, outs3;
39  outs << "Forcing scale to : " << scale << endl;
40  SpParHelper::Print(outs.str());
41  DEL->GenGraph500Data(initiator, scale, EDGEFACTOR, true, true ); // generate packed edges
42  SpParHelper::Print("Generated renamed edge lists\n");
43 
44  // conversion from distributed edge list, keeps self-loops, sums duplicates
45  A = new PARDBMAT(*DEL, false); // already creates renumbered vertices (hence balanced)
46  delete DEL; // free memory before symmetricizing
47  SpParHelper::Print("Created double Sparse Matrix\n");
48 
49  float balance = A->LoadImbalance();
50  outs2 << "Load balance: " << balance << endl;
51  SpParHelper::Print(outs2.str());
52  A->PrintInfo();
53 
54  for(unsigned i=1; i<4; i++)
55  {
56  DEL = new DistEdgeList<int64_t>();
57  DEL->GenGraph500Data(initiator, scale-i, ((double) EDGEFACTOR) / pow(2.0,i) , true, true ); // "i" scale smaller
58  B = new PARDBMAT(*DEL, false);
59  delete DEL;
60  SpParHelper::Print("Created RHS Matrix\n");
61  B->PrintInfo();
62  FullyDistVec<int,int> perm; // get a different permutation
63  perm.iota(A->getnrow(), 0);
64  perm.RandPerm();
65 
66  //void FullyDistVec::iota(IT globalsize, NT first)
68  sel.iota(B->getnrow(), 0);
69  perm = perm(sel); // just get the first B->getnrow() entries of the permutation
70  perm.PrintInfo("Index vector");
71 
72  A->SpAsgn(perm,perm,*B); // overriding A with a structurally similar piece.
73  A->PrintInfo();
74 
75  double t1 = MPI_Wtime();
76  for(int j=0; j< ITERATIONS; ++j)
77  {
78  A->SpAsgn(perm,perm,*B);
79  }
80  double t2 = MPI_Wtime();
81 
82  if(myrank == 0)
83  {
84  cout<< "Scale " << scale-i << " assignment iterations finished"<<endl;
85  printf("%.6lf seconds elapsed per iteration\n", (t2-t1)/(double)ITERATIONS);
86 
87  }
88  delete B;
89  }
90  }
91  MPI_Finalize();
92  return 0;
93 }