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