COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
RandomParentBFS.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#ifdef THREADED
12 #ifndef _OPENMP
13 #define _OPENMP
14 #endif
15#endif
16
17#ifdef _OPENMP
18 #include <omp.h>
19#endif
20
21
27#ifdef _OPENMP
29#else
31#endif
32
33#define ITERS 16
34#define EDGEFACTOR 16
35using namespace std;
36using namespace combblas;
37
38
39MTRand GlobalMT(123); // for reproducable result
40
41
42template <typename PARMAT>
44{
45 // boolean addition is practically a "logical or"
46 // therefore this doesn't destruct any links
47 PARMAT AT = A;
48 AT.Transpose();
49 A += AT;
50}
51
52
53
55{
56public:
57 ParentType(){parent=-1; p = 0;};
59 friend ostream& operator<<(ostream& os, const ParentType & vertex ){os << "Parent=" << vertex.parent << " p=" << vertex.p; return os;};
60 //private:
62 float p;
63
64};
65
66
67
68
69
70
71struct Edge_randomizer : public std::unary_function<std::pair<bool, float>, std::pair<bool, float>>
72{
73 const std::pair<bool, float> operator()(const std::pair<bool, float> & x) const
74 {
75 float edgeRand = static_cast<float>(rand()); // random range(0,1)
76 return std::pair<bool, float>(x.first, edgeRand);
77 }
78};
79
80
81
82static void MPI_randuniq(void * invec, void * inoutvec, int * len, MPI_Datatype *datatype)
83{
87 for (int i=0; i<*len; i++ )
89}
90
91
93{
96 static ParentType id(){ return ParentType(); };
97 static bool returnedSAID() { return false; }
98 //static MPI_Op mpi_op() { return MPI_MAX; }; // do we need this?
99
100 static ParentType add(const ParentType & arg1, const ParentType & arg2)
101 {
102 //cout << arg1 << " ;;; " << arg2 << endl;
103 if(arg1.p < arg2.p) return arg1;
104 else return arg2;
105 }
106
108 {
110 temp.parent = arg2.parent;
111 temp.p = GlobalMT.rand();
112 return temp;
113 }
114
115 static void axpy(T_promote a, const ParentType & x, ParentType & y)
116 {
117 y = add(y, multiply(a, x));
118 }
119};
120
121
122
123// This one is used for BFS iteration
125{
127 static T_promote id(){ return -1; };
128 static bool returnedSAID() { return false; }
129 //static MPI_Op mpi_op() { return MPI_MAX; };
130
131 static T_promote add(const T_promote & arg1, const T_promote & arg2)
132 {
133 cout << arg1 << " a " << arg2 << endl;
134 return std::max(arg1, arg2);
135 }
136
137 static T_promote multiply(const bool & arg1, const T_promote & arg2)
138 {
139 cout << arg1 << " m " << arg2 << endl;
140 return arg2;
141 }
142 /*
143 static void axpy(bool a, const T_promote & x, T_promote & y)
144 {
145 y = std::max(y, x);
146 }*/
147};
148
152
153
155{
156
157 int nprocs, myrank;
160
161
162 FullyDistSpVec<int64_t, ParentType> fringe(Aeff.getcommgrid(), Aeff.getncol());
163
164 fringe.SetElement(0, ParentType(0));
165 fringe.SetElement(1, ParentType(1));
166 fringe.SetElement(5, ParentType(5));
167 fringe.SetElement(6, ParentType(6));
168 fringe.SetElement(7, ParentType(7));
169
171 //A.PrintInfo();
173 fringe.DebugPrint();
174}
175
176
177
178int main(int argc, char* argv[])
179{
180 int nprocs, myrank;
181 MPI_Init(&argc, &argv);
184 if(argc < 2)
185 {
186 if(myrank == 0)
187 {
188 cout << "Usage: ./rpbfs <Scale>" << endl;
189 cout << "Example: mpirun -np 4 ./spbfs 20" << endl;
190 }
191 MPI_Finalize();
192 return -1;
193 }
194 {
195
196 // Declare objects
198 FullyDistVec<int64_t, int64_t> nonisov; // id's of non-isolated (connected) vertices
199 unsigned scale;
200
201 scale = static_cast<unsigned>(atoi(argv[1]));
202 double initiator[4] = {.57, .19, .19, .05};
203
204
208
209
210 PSpMat_Bool * ABool = new PSpMat_Bool(*DEL, false);
211 delete DEL;
212 int64_t removed = ABool->RemoveLoops();
213 ABool->PrintInfo();
214
217
218 ABool->Reduce(*ColSums, Column, plus<int64_t>(), static_cast<int64_t>(0));
219 ABool->Reduce(*RowSums, Row, plus<int64_t>(), static_cast<int64_t>(0));
220 ColSums->EWiseApply(*RowSums, plus<int64_t>());
221 delete RowSums;
222 nonisov = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
223 delete ColSums;
224 nonisov.RandPerm(); // so that A(v,v) is load-balanced (both memory and time wise)
225#ifndef NOPERMUTE
226 ABool->operator()(nonisov, nonisov, true); // in-place permute to save memory
227#endif
228
229
231
232
233 }
234 MPI_Finalize();
235 return 0;
236}
237
238
239
240
241
242
243
244
245
int main()
Definition Driver.cpp:12
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > PSpMat_Bool
double cblas_transvectime
double cblas_localspmvtime
int cblas_splits
double cblas_allgathertime
double cblas_alltoalltime
MTRand GlobalMT(123)
SpParMat< int64_t, bool, SpDCCols< int32_t, bool > > PSpMat_s32p64
void RandomParentBFS(PSpMat_Bool &Aeff)
#define EDGEFACTOR
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > PSpMat_Int64
double cblas_mergeconttime
double rand()
void GenGraph500Data(double initiator[4], int log_numverts, int edgefactor, bool scramble=false, bool packed=false)
int nprocs
Definition comms.cpp:55
@ Column
Definition SpDefs.h:115
void Symmetricize(PARMAT &A)
Definition ReadMatDist.h:29
double A
const std::pair< bool, float > operator()(const std::pair< bool, float > &x) const
friend ostream & operator<<(ostream &os, const ParentType &vertex)
ParentType(int64_t x)
static bool returnedSAID()
static T_promote id()
static T_promote multiply(const bool &arg1, const T_promote &arg2)
static T_promote add(const T_promote &arg1, const T_promote &arg2)
static ParentType multiply(const T_promote &arg1, const ParentType &arg2)
static ParentType id()
static bool returnedSAID()
static ParentType add(const ParentType &arg1, const ParentType &arg2)
static void axpy(T_promote a, const ParentType &x, ParentType &y)
static MPI_Op MPI_BFSRAND