COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
auction.cpp
Go to the documentation of this file.
1#include "CombBLAS/CombBLAS.h"
2#include <mpi.h>
3#include <sys/time.h>
4#include <iostream>
5#include <functional>
6#include <algorithm>
7#include <vector>
8#include <string>
9#include <sstream>
10
11
12using namespace std;
13using namespace combblas;
14
16int init;
18bool fewexp;
19
20
21template <typename PARMAT>
23{
24 // boolean addition is practically a "logical or"
25 // therefore this doesn't destruct any links
26 PARMAT AT = A;
27 AT.Transpose();
28 A += AT;
29}
30
31
32
34{
35public:
36
37 // careful: secondMaxProfit should be -\max wij by default
38 // how can we pass that info to static VertexType id() so that SpMV can use it?
39 VertexType(int64_t oi=-1, double p = 0, double smp=-9999999)
40 {
41 objID = oi;
42 price = p;
44 };
45
46
47 friend bool operator<(const VertexType & vtx1, const VertexType & vtx2 ){return vtx1.price < vtx2.price;};
48 friend ostream& operator<<(ostream& os, const VertexType & vertex ){os << "(" << vertex.objID << ", " << vertex.price << ", " << vertex.secondMaxProfit << ")"; return os;};
49
50 // variables for the object
51 int64_t objID; // ultimately it will be the object with the max profit for a bidder
52 double price; //ultimately it will be the max profit for a bidder
53 double secondMaxProfit; //it will be the 2nd max profit for a bidder
54
55
56};
57
58
59// return the maximum and second maximum profits
60struct max2 : public std::binary_function<VertexType, VertexType, VertexType>
61{
62 const VertexType operator()(const VertexType& x, const VertexType& y) const
63 {
65 if(x.price > y.price)
66 {
67 ret = x;
68 if(x.secondMaxProfit < y.price) ret.secondMaxProfit = y.price;
69 }
70 else
71 {
72 ret = y;
73 if(y.secondMaxProfit < x.price) ret.secondMaxProfit = x.price;
74 }
75
76 return ret;
77 }
78};
79
80
81
83{
84 static VertexType id(){ return VertexType(); };
85 static bool returnedSAID() { return false; }
87
88
89
90 // addition compute the maximum and second maximum profit
91 static VertexType add(const VertexType & arg1, const VertexType & arg2)
92 {
93 return max2()(arg1, arg2);
94 }
95
96 // multiplication computes profit of an object to a bidder
97 // profit_j = cij - price_j
98 static VertexType multiply(const double & arg1, const VertexType & arg2)
99 {
100 // if the matrix is boolean, arg1=1
101 return VertexType(arg2.objID, arg1 - arg2.price, arg2.secondMaxProfit);
102 }
103
104 static void axpy(const double a, const VertexType & x, VertexType & y)
105 {
106 y = add(y, multiply(a, x));
107 }
108};
109
110
111
120template <class IT, class NT>
123
124
125
126
127/*
128 Remove isolated vertices and purmute
129 */
131{
132
133 int nprocs, myrank;
136
137
140 FullyDistVec<int64_t, int64_t> nonisoRowV; // id's of non-isolated (connected) Row vertices
141 FullyDistVec<int64_t, int64_t> nonisoColV; // id's of non-isolated (connected) Col vertices
142 FullyDistVec<int64_t, int64_t> nonisov; // id's of non-isolated (connected) vertices
143
144 A.Reduce(*ColSums, Column, plus<int64_t>(), static_cast<int64_t>(0));
145 A.Reduce(*RowSums, Row, plus<int64_t>(), static_cast<int64_t>(0));
146
147 // this steps for general graph
148 /*
149 ColSums->EWiseApply(*RowSums, plus<int64_t>()); not needed for bipartite graph
150 nonisov = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
151 nonisov.RandPerm(); // so that A(v,v) is load-balanced (both memory and time wise)
152 A.operator()(nonisov, nonisov, true); // in-place permute to save memory
153 */
154
155 // this steps for bipartite graph
156 nonisoColV = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
157 nonisoRowV = RowSums->FindInds(bind2nd(greater<int64_t>(), 0));
158 delete ColSums;
159 delete RowSums;
160
161
162 {
163 nonisoColV.RandPerm();
164 nonisoRowV.RandPerm();
165 }
166
167
168 int64_t nrows1=A.getnrow(), ncols1=A.getncol(), nnz1 = A.getnnz();
169 double avgDeg1 = (double) nnz1/(nrows1+ncols1);
170
171
172 A.operator()(nonisoRowV, nonisoColV, true);
173
174 int64_t nrows2=A.getnrow(), ncols2=A.getncol(), nnz2 = A.getnnz();
175 double avgDeg2 = (double) nnz2/(nrows2+ncols2);
176
177
178 if(myrank == 0)
179 {
180 cout << "ncol nrows nedges deg \n";
181 cout << nrows1 << " " << ncols1 << " " << nnz1 << " " << avgDeg1 << " \n";
182 cout << nrows2 << " " << ncols2 << " " << nnz2 << " " << avgDeg2 << " \n";
183 }
184
186
187
188}
189
190
191
192
193
194int main(int argc, char* argv[])
195{
196
197 // ------------ initialize MPI ---------------
198 int provided;
201 {
202 printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
204 }
205 int nprocs, myrank;
208 if(argc < 3)
209 {
210 //ShowUsage();
211 MPI_Finalize();
212 return -1;
213 }
214
215
216 // ------------ Process input arguments and build matrix ---------------
217 {
218
222 double t01, t02;
223 if(string(argv[1]) == string("input")) // input option
224 {
225 ABool = new PSpMat_Bool();
226
227 string filename(argv[2]);
228 tinfo.str("");
229 tinfo << "**** Reading input matrix: " << filename << " ******* " << endl;
231 t01 = MPI_Wtime();
232 ABool->ParallelReadMM(filename, false, maximum<bool>());
233 t02 = MPI_Wtime();
234 ABool->PrintInfo();
235 tinfo.str("");
236 tinfo << "Reader took " << t02-t01 << " seconds" << endl;
238 //GetOptions(argv+3, argc-3);
239
240 }
241 else if(argc < 4)
242 {
243 //ShowUsage();
244 MPI_Finalize();
245 return -1;
246 }
247 else
248 {
249 unsigned scale = (unsigned) atoi(argv[2]);
250 unsigned EDGEFACTOR = (unsigned) atoi(argv[3]);
251 double initiator[4];
252 if(string(argv[1]) == string("er"))
253 {
254 initiator[0] = .25;
255 initiator[1] = .25;
256 initiator[2] = .25;
257 initiator[3] = .25;
258 cout << "ER ******** \n";
259 }
260 else if(string(argv[1]) == string("g500"))
261 {
262 initiator[0] = .57;
263 initiator[1] = .19;
264 initiator[2] = .19;
265 initiator[3] = .05;
266 cout << "g500 ******** \n";
267 }
268 else if(string(argv[1]) == string("ssca"))
269 {
270 initiator[0] = .6;
271 initiator[1] = .4/3;
272 initiator[2] = .4/3;
273 initiator[3] = .4/3;
274 cout << "ER ******** \n";
275 }
276 else
277 {
278 if(myrank == 0)
279 printf("The input type - %s - is not recognized.\n", argv[2]);
281 }
282
283 SpParHelper::Print("Generating input matrix....\n");
284 t01 = MPI_Wtime();
287 ABool = new PSpMat_Bool(*DEL, false);
288 delete DEL;
289 t02 = MPI_Wtime();
290 ABool->PrintInfo();
291 tinfo.str("");
292 tinfo << "Generator took " << t02-t01 << " seconds" << endl;
294
296 //removeIsolated(*ABool);
297 SpParHelper::Print("Generated matrix symmetricized....\n");
298 ABool->PrintInfo();
299
300 //GetOptions(argv+4, argc-4);
301
302 }
303
304
305 // randomly permute for load balance
306 SpParHelper::Print("Performing random permutation of matrix.\n");
309 prow.iota(ABool->getnrow(), 0);
310 pcol.iota(ABool->getncol(), 0);
311 prow.RandPerm();
312 pcol.RandPerm();
313 (*ABool)(prow, pcol, true);
314 SpParHelper::Print("Performed random permutation of matrix.\n");
315
316
318
319
320 FullyDistVec<int64_t, int64_t> mateRow2Col ( A.getcommgrid(), A.getnrow(), (int64_t) -1);
321 FullyDistVec<int64_t, int64_t> mateCol2Row ( A.getcommgrid(), A.getncol(), (int64_t) -1);
322
324
325
326
327 }
328 MPI_Finalize();
329 return 0;
330}
331
332
333
334
335
336
337
340{
341
342 int nprocs, myrank;
345
346 int64_t nrow = A.getnrow();
347 int64_t ncol = A.getncol();
350 FullyDistVec<int64_t, int64_t> leaves ( A.getcommgrid(), ncol, (int64_t) -1);
351
353 vector<int> layers;
356
357 bool matched = true;
358 int phase = 0;
359 int totalLayer = 0;
361
362
363
365 // set index to value
366 objects.ApplyInd([](VertexType vtx, int64_t idx){return VertexType(idx, vtx.price, vtx.secondMaxProfit);});
368
369 // future: I need an SPMV that could return a dense vector of type other than the input
371
372
373 // Remove bidders without any profitable objects
375
376 // remove matched bidders
377 // we have done unncessary computation for already matched bidders
378 // option1: bottom up style
379 // option2: masked SpMV
380
382 [](VertexType vtx, int64_t mate){return vtx;},
383 [](VertexType vtx, int64_t mate){return mate==-1;},
384 false, VertexType());
385
386
387 // compute bid
388 activeBidders.Apply([](VertexType bidder){return VertexType(bidder.objID, bidder.price - bidder.secondMaxProfit);}); // I don't care secondMaxProfit anymore
389
390
391 // place bid
392 // objects need to select the best bidder
393 activeBidders.DebugPrint();
395 activeBidders.Invert(ncol,
396 [](VertexType bidder, int64_t idx){return bidder.objID;},
397 [](VertexType bidder, int64_t idx){return VertexType(idx, bidder.price);},
398 [](VertexType bid1, VertexType bid2){return bid1.price>bid2.price? bid1: bid2;}); // I don't care secondMaxProfit anymore
399
400 bidObject.DebugPrint();
401
402 // just creating a simplified object with the highest bidder
403 // mateCol2Row is used just as a mask
405 [](VertexType vtx, int64_t mate){return vtx.objID;},
406 [](VertexType vtx, int64_t mate){return true;},
407 false, VertexType());
408
409
410 //mateCol2Row.DebugPrint();
411 // bidders previously matched to current successfull bids will become unmatched
413 [](int64_t newbidder, int64_t mate){return mate;},
414 [](int64_t newbidder, int64_t mate){return mate!=-1;},
415 false, (int64_t)-1);
416
417
419
420 //cout << " djkfhksjdfh \n";
421 //successfullBids.DebugPrint();
422 //revokedBids.DebugPrint();
423
424 // previously unmatched bidders that will be matched
426
427 // previously matched bidders that will be unmatched
429 // they are mutually exclusive
430 successfullBidders.DebugPrint();
431 revokedBidders.DebugPrint();
432
434 revokedBidders.Apply([](int64_t prevmate){return (int64_t)-1;});
435
436 //mateRow2Col.Set(revokedBidders);
437
438
439
440
441
442 //objects.DebugPrint();
443
444}
445
#define EDGEFACTOR
Definition DirOptBFS.cpp:81
int main()
Definition Driver.cpp:12
bool moreSplit
Definition auction.cpp:15
bool mvInvertMate
Definition auction.cpp:15
void auction(PSpMat_s32p64 &A, FullyDistVec< int64_t, int64_t > &mateRow2Col, FullyDistVec< int64_t, int64_t > &mateCol2Row)
Definition auction.cpp:338
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > PSpMat_Bool
Definition auction.cpp:112
bool randMaximal
Definition auction.cpp:17
bool fewexp
Definition auction.cpp:18
bool prune
Definition auction.cpp:15
int init
Definition auction.cpp:16
SpParMat< int64_t, bool, SpDCCols< int32_t, bool > > PSpMat_s32p64
Definition auction.cpp:113
bool randMM
Definition auction.cpp:15
bool isMaximalmatching(PSpMat_Int64 &A, FullyDistVec< IT, NT > &mateRow2Col, FullyDistVec< IT, NT > &mateCol2Row, FullyDistSpVec< int64_t, int64_t > unmatchedRow, FullyDistSpVec< int64_t, int64_t > unmatchedCol)
SpParMat< int64_t, float, SpDCCols< int64_t, float > > PSpMat_float
Definition auction.cpp:115
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > PSpMat_Int64
Definition auction.cpp:114
void maximumMatching(PSpMat_s32p64 &Aeff, FullyDistVec< int64_t, int64_t > &mateRow2Col, FullyDistVec< int64_t, int64_t > &mateCol2Row)
void removeIsolated(PSpMat_Bool &A)
Definition auction.cpp:130
void GenGraph500Data(double initiator[4], int log_numverts, int edgefactor, bool scramble=false, bool packed=false)
static void Print(const std::string &s)
int nprocs
Definition comms.cpp:55
long int64_t
Definition compat.h:21
@ Column
Definition SpDefs.h:115
void Symmetricize(PARMAT &A)
Definition ReadMatDist.h:29
double A
static VertexType multiply(const double &arg1, const VertexType &arg2)
Definition auction.cpp:98
static bool returnedSAID()
Definition auction.cpp:85
static VertexType add(const VertexType &arg1, const VertexType &arg2)
Definition auction.cpp:91
static void axpy(const double a, const VertexType &x, VertexType &y)
Definition auction.cpp:104
static MPI_Op mpi_op()
Definition auction.cpp:86
static VertexType id()
Definition auction.cpp:84
friend bool operator<(const VertexType &vtx1, const VertexType &vtx2)
Definition auction.cpp:47
int64_t objID
Definition auction.cpp:51
double price
Definition auction.cpp:52
VertexType(int64_t oi=-1, double p=0, double smp=-9999999)
Definition auction.cpp:39
double secondMaxProfit
Definition auction.cpp:53
friend ostream & operator<<(ostream &os, const VertexType &vertex)
Definition auction.cpp:48
static MPI_Op op()
Definition MPIOp.h:80
const VertexType operator()(const VertexType &x, const VertexType &y) const
Definition auction.cpp:62