COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
gathertest.cpp
Go to the documentation of this file.
1//#define DETERMINISTIC 1
2
3#ifdef THREADED
4#ifndef _OPENMP
5#define _OPENMP // should be defined before any COMBBLAS header is included
6#endif
7#include <omp.h>
8#endif
9
10#include "CombBLAS/CombBLAS.h"
11#include <mpi.h>
12#include <sys/time.h>
13#include <iostream>
14#include <functional>
15#include <algorithm>
16#include <vector>
17#include <string>
18#include <sstream>
19
20
21#define EDGEFACTOR 16
22#define RAND_PERMUTE 1
23
24#ifdef DETERMINISTIC
26#else
27MTRand GlobalMT; // generate random numbers with Mersenne Twister
28#endif
29
35
36
37
38using namespace std;
39using namespace combblas;
40
41
42
43template <typename PARMAT>
45{
46 // boolean addition is practically a "logical or"
47 // therefore this doesn't destruct any links
48 PARMAT AT = A;
49 AT.Transpose();
50 AT.RemoveLoops(); // not needed for boolean matrices, but no harm in keeping it
51 A += AT;
52}
53
54
55
56
61
62
63
64
65int main(int argc, char* argv[])
66{
67 int nprocs, myrank;
68 MPI_Init(&argc, &argv);
71 if(argc < 3)
72 {
73 if(myrank == 0)
74 {
75 cout << "Usage: ./rcm <rmat|er|input> <scale|filename> <splitPerThread>" << endl;
76 cout << "Example: mpirun -np 4 ./rcm rmat 20" << endl;
77 cout << "Example: mpirun -np 4 ./rcm er 20" << endl;
78 cout << "Example: mpirun -np 4 ./rcm input a.mtx" << endl;
79
80 }
82 return -1;
83 }
84 {
87 bool unsym=false;
89
90 if(string(argv[1]) == string("input")) // input option
91 {
92 ABool = new Par_DCSC_Bool();
93 string filename(argv[2]);
94 tinfo.str("");
95 tinfo << "**** Reading input matrix: " << filename << " ******* " << endl;
96
97
99 double t01 = MPI_Wtime();
100 ABool->ParallelReadMM(filename, false, maximum<bool>());
101 double t02 = MPI_Wtime();
102 int64_t bw = ABool->Bandwidth();
103 int64_t pf = ABool->Profile();
104 tinfo.str("");
105 tinfo << "Reader took " << t02-t01 << " seconds" << endl;
106 tinfo << "Bandwidth before random permutation " << bw << endl;
107 tinfo << "Profile before random permutation " << pf << endl;
109
110
111
112
113 }
114 else if(string(argv[1]) == string("rmat"))
115 {
116 unsigned scale;
117 scale = static_cast<unsigned>(atoi(argv[2]));
118 double initiator[4] = {.57, .19, .19, .05};
120 DEL->GenGraph500Data(initiator, scale, EDGEFACTOR, true, false );
122
123 ABool = new Par_DCSC_Bool(*DEL, false);
125 delete DEL;
126 }
127 else if(string(argv[1]) == string("er"))
128 {
129 unsigned scale;
130 scale = static_cast<unsigned>(atoi(argv[2]));
131 double initiator[4] = {.25, .25, .25, .25};
133 DEL->GenGraph500Data(initiator, scale, EDGEFACTOR, true, false );
135
136 ABool = new Par_DCSC_Bool(*DEL, false);
138 delete DEL;
139 }
140 else
141 {
142 SpParHelper::Print("Unknown input option\n");
143 MPI_Finalize();
144 return -1;
145 }
146
147
148
149
150 MPI_Comm com = ABool->getcommgrid()->GetWorld();
151 double gtime = MPI_Wtime();
152 SpParHelper::GatherMatrix(com, ABool->seq(), (int)0);
153 if(myrank==0)
154 {
155 cout << "gathertime " << MPI_Wtime() - gtime << endl;
156 }
157
158
159
160
161
162
163 delete ABool;
164
165
166 }
167 MPI_Finalize();
168 return 0;
169}
170
int main()
Definition Driver.cpp:12
MTRand GlobalMT
double cblas_transvectime
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > Par_DCSC_Bool
double cblas_localspmvtime
double cblas_allgathertime
double cblas_alltoalltime
SpParMat< int64_t, bool, SpCCols< int64_t, bool > > Par_CSC_Bool
SpParMat< int64_t, double, SpDCCols< int64_t, double > > Par_DCSC_Double
#define EDGEFACTOR
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > Par_DCSC_int64_t
double cblas_mergeconttime
void GenGraph500Data(double initiator[4], int log_numverts, int edgefactor, bool scramble=false, bool packed=false)
static void GatherMatrix(MPI_Comm &comm1d, SpMat< IT, NT, DER > &Matrix, int root)
static void Print(const std::string &s)
int nprocs
Definition comms.cpp:55
void Symmetricize(PARMAT &A)
Definition ReadMatDist.h:29
double A