COMBINATORIAL_BLAS  1.6
Roofline.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 #include "../Applications/TwitterEdge.h"
10 
11 using namespace std;
12 using namespace combblas;
13 
14 // One edge = 16 bytes (rounded up from 11)
15 // One parent = 8 bytes
16 #define INC 256
17 #define L1 4096 // maximum entries of edges+parents combined to fit 32 KB
18 #define REPEAT 1000
19 
20 
21 
22 
23 // all the local variables before the EWiseApply wouldn't be accessible,
24 // so we need to pass them as parameters to the functor constructor
25 class twitter_mult : public std::binary_function<ParentType, TwitterEdge, ParentType>
26 {
27 private:
28  time_t sincedate;
29 public:
30  twitter_mult(time_t since):sincedate(since) {};
31  ParentType operator()(const ParentType & arg1, const TwitterEdge & arg2) const
32  {
33  if(arg2.isFollower() && arg2.TweetSince(sincedate))
34  return arg1;
35  else
36  return ParentType(); // null-type parent id
37  }
38 };
39 
40 
41 int main(int argc, char* argv[])
42 {
43  int nprocs, myrank;
44  MPI_Init(&argc, &argv);
45  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
46  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
47 
48  {
49  int64_t len = INC;
50  time_t now;
51  time ( &now );
52  TwitterEdge twe(4, 1, now); // 4 retweets, latest now, following
53 
54  while(len < L1)
55  {
56  // FullyDistVec(IT globallen, NT initval)
57  FullyDistVec<int64_t,TwitterEdge> tvec(nprocs * len, twe);
59  pvec.iota(nprocs * len, ParentType());
60 
61  MPI_Barrier(MPI_COMM_WORLD);
62  double t1 = MPI_Wtime(); // initilize (wall-clock) timer
63 
64  time_t now = time(0);
65  struct tm * timeinfo = localtime( &now);
66  timeinfo->tm_mon = timeinfo->tm_mon-1;
67  time_t monthago = mktime(timeinfo);
68 
69  for(int i=0; i< REPEAT; ++i)
70  pvec.EWiseApply(tvec, twitter_mult(monthago));
71 
72  MPI_Barrier(MPI_COMM_WORLD);
73  double t2 = MPI_Wtime();
74 
75  if(myrank == 0)
76  {
77  cout<<"EWiseApply Iterations finished"<<endl;
78  double time = t2-t1;
79  double teps = (nprocs*len*REPEAT) / (time * 1000000);
80  printf("%.6lf seconds elapsed for %d iterations on vector of length %lld\n", time, REPEAT, nprocs*len);
81  printf("%.6lf million TEPS per second\n", teps);
82  }
83  len += INC;
84  }
85  }
86  MPI_Finalize();
87  return 0;
88 }
bool TweetSince(time_t begin) const
Definition: TwitterEdge.h:26
ParentType operator()(const ParentType &arg1, const TwitterEdge &arg2) const
Definition: Roofline.cpp:31
void EWiseApply(const FullyDistVec< IT, NT2 > &other, _BinaryOperation __binary_op, _BinaryPredicate _do_op, const bool useExtendedBinOp)
#define REPEAT
Definition: Roofline.cpp:18
#define L1
Definition: Roofline.cpp:17
void iota(IT globalsize, NT first)
long int64_t
Definition: compat.h:21
int main(int argc, char *argv[])
Definition: Roofline.cpp:41
twitter_mult(time_t since)
Definition: Roofline.cpp:30
Definition: CCGrid.h:4
#define INC
Definition: Roofline.cpp:16
bool isFollower() const
Definition: TwitterEdge.h:23