COMBINATORIAL_BLAS  1.6
CC.cpp
Go to the documentation of this file.
1 /****************************************************************/
2 /* Parallel Combinatorial BLAS Library (for Graph Computations) */
3 /* version 1.6 -------------------------------------------------*/
4 /* date: 6/15/2017 ---------------------------------------------*/
5 /* authors: Ariful Azad, Aydin Buluc --------------------------*/
6 /****************************************************************/
7 /*
8  Copyright (c) 2010-2017, The Regents of the University of California
9 
10  Permission is hereby granted, free of charge, to any person obtaining a copy
11  of this software and associated documentation files (the "Software"), to deal
12  in the Software without restriction, including without limitation the rights
13  to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
14  copies of the Software, and to permit persons to whom the Software is
15  furnished to do so, subject to the following conditions:
16 
17  The above copyright notice and this permission notice shall be included in
18  all copies or substantial portions of the Software.
19 
20  THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
21  IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
22  FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
23  AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
24  LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
25  OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
26  THE SOFTWARE.
27  */
28 
29 
30 #include <mpi.h>
31 
32 // These macros should be defined before stdint.h is included
33 #ifndef __STDC_CONSTANT_MACROS
34 #define __STDC_CONSTANT_MACROS
35 #endif
36 #ifndef __STDC_LIMIT_MACROS
37 #define __STDC_LIMIT_MACROS
38 #endif
39 #include <stdint.h>
40 
41 #include <sys/time.h>
42 #include <iostream>
43 #include <fstream>
44 #include <string>
45 #include <sstream>
46 #include <ctime>
47 #include <cmath>
48 #include "CombBLAS/CombBLAS.h"
49 #include "CC.h"
50 
51 using namespace std;
52 using namespace combblas;
53 
59 class Dist
60 {
61 public:
64 };
65 
66 
67 
68 int main(int argc, char* argv[])
69 {
70  int provided;
71  MPI_Init_thread(&argc, &argv, MPI_THREAD_SERIALIZED, &provided);
72  if (provided < MPI_THREAD_SERIALIZED)
73  {
74  printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
75  MPI_Abort(MPI_COMM_WORLD, 1);
76  }
77 
78  int nthreads = 1;
79 #ifdef THREADED
80 #pragma omp parallel
81  {
82  nthreads = omp_get_num_threads();
83  }
84 #endif
85 
86  int nprocs, myrank;
87  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
88  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
89  if(myrank == 0)
90  {
91  cout << "Process Grid (p x p x t): " << sqrt(nprocs) << " x " << sqrt(nprocs) << " x " << nthreads << endl;
92  }
93 
94  if(argc < 3)
95  {
96  if(myrank == 0)
97  {
98  cout << "Usage: ./cc -I <mm|triples> -M <FILENAME_MATRIX_MARKET> (required)\n";
99  cout << "-I <INPUT FILE TYPE> (mm: matrix market, triples: (vtx1, vtx2, edge_weight) triples. default:mm)\n";
100  cout << "-base <BASE OF MATRIX MARKET> (default:1)\n";
101  cout << "-rand <RANDOMLY PERMUTE VERTICES> (default:0)\n";
102  cout << "Example (0-indexed mtx with random permutation): ./cc -M input.mtx -base 0 -rand 1" << endl;
103  cout << "Example (triples format): ./cc -I triples -M input.txt" << endl;
104  }
105  MPI_Finalize();
106  return -1;
107  }
108  {
109  string ifilename = "";
110  int base = 1;
111  int randpermute = 0;
112  bool isMatrixMarket = true;
113 
114  for (int i = 1; i < argc; i++)
115  {
116  if (strcmp(argv[i],"-I")==0)
117  {
118  string ifiletype = string(argv[i+1]);
119  if(ifiletype == "triples") isMatrixMarket = false;
120  }
121  if (strcmp(argv[i],"-M")==0)
122  {
123  ifilename = string(argv[i+1]);
124  if(myrank == 0) printf("filename: %s",ifilename.c_str());
125  }
126  else if (strcmp(argv[i],"-base")==0)
127  {
128  base = atoi(argv[i + 1]);
129  if(myrank == 0) printf("\nBase of MM (1 or 0):%d",base);
130  }
131  else if (strcmp(argv[i],"-rand")==0)
132  {
133  randpermute = atoi(argv[i + 1]);
134  if(myrank == 0) printf("\nRandomly permute the matrix? (1 or 0):%d",randpermute);
135  }
136  }
137 
138  double tIO = MPI_Wtime();
139  Dist::MPI_DCCols A; // construct object
140 
141  if(isMatrixMarket)
142  A.ParallelReadMM(ifilename, base, maximum<bool>());
143  else
144  A.ReadGeneralizedTuples(ifilename, maximum<bool>());
145  A.PrintInfo();
146 
147  Dist::MPI_DCCols AT = A;
148  AT.Transpose();
149  if(!(AT == A))
150  {
151  SpParHelper::Print("Symmatricizing an unsymmetric input matrix.\n");
152  A += AT;
153  }
154  A.PrintInfo();
155 
156 
157  ostringstream outs;
158  outs << "File Read time: " << MPI_Wtime() - tIO << endl;
159  SpParHelper::Print(outs.str());
160 
161  if(randpermute && isMatrixMarket) // no need for GeneralizedTuples
162  {
163  // randomly permute for load balance
164  if(A.getnrow() == A.getncol())
165  {
167  p.iota(A.getnrow(), 0);
168  p.RandPerm();
169  (A)(p,p,true);// in-place permute to save memory
170  SpParHelper::Print("Applied symmetric permutation.\n");
171  }
172  else
173  {
174  SpParHelper::Print("Rectangular matrix: Can not apply symmetric permutation.\n");
175  }
176  }
177 
178  FullyDistVec<int64_t,double> ColSums = A.Reduce(Column, plus<double>(), 0.0);
179  FullyDistVec<int64_t, int64_t> isov = ColSums.FindInds(bind2nd(equal_to<double>(), 0));
180  outs.str("");
181  outs.clear();
182  outs << "isolated vertice: " << isov.TotalLength() << endl;
183  SpParHelper::Print(outs.str());
184 
185  float balance = A.LoadImbalance();
186  int64_t nnz = A.getnnz();
187  outs.str("");
188  outs.clear();
189  outs << "Load balance: " << balance << endl;
190  outs << "Nonzeros: " << nnz << endl;
191  SpParHelper::Print(outs.str());
192  double t1 = MPI_Wtime();
193 
194  int64_t nCC = 0;
195  FullyDistVec<int64_t, int64_t> cclabels = CC(A, nCC);
196 
197  double t2 = MPI_Wtime();
198  //outs.str("");
199  //outs.clear();
200  //outs << "Total time: " << t2 - t1 << endl;
201  //SpParHelper::Print(outs.str());
202  //HistCC(cclabels, nCC);
203 
204  int64_t nclusters = cclabels.Reduce(maximum<int64_t>(), (int64_t) 0 ) ;
205  nclusters ++; // because of zero based indexing for clusters
206 
207  double tend = MPI_Wtime();
208  stringstream s2;
209  s2 << "Number of clusters: " << nclusters << endl;
210  s2 << "Total time: " << (t2-t1) << endl;
211  s2 << "=================================================\n" << endl ;
212  SpParHelper::Print(s2.str());
213 
214  }
215 
216  MPI_Finalize();
217  return 0;
218 }
FullyDistVec< IT, NT > Reduce(Dim dim, _BinaryOperation __binary_op, NT id, _UnaryOperation __unary_op) const
Definition: SpParMat.cpp:791
std::shared_ptr< CommGrid > getcommgrid() const
Definition: SpParMat.h:275
FullyDistVec< IT, std::array< char, MAXVERTNAME > > ReadGeneralizedTuples(const std::string &, _BinaryOperation)
Definition: SpParMat.cpp:3319
Compute the maximum of two values.
Definition: Operations.h:154
SpDCCols< int64_t, double > DCCols
Definition: CC.cpp:62
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
Return the indices where pred is true.
IT getnnz() const
Definition: SpParMat.cpp:676
float LoadImbalance() const
Definition: SpParMat.cpp:665
FullyDistVec< IT, IT > CC(SpParMat< IT, NT, DER > &A, IT &nCC)
Definition: CC.h:276
double A
void iota(IT globalsize, NT first)
long int64_t
Definition: compat.h:21
IT getncol() const
Definition: SpParMat.cpp:694
void PrintInfo() const
Definition: SpParMat.cpp:2387
double tIO
Definition: MCL.cpp:65
SpParMat< int64_t, double, DCCols > MPI_DCCols
Definition: CC.cpp:63
Definition: CCGrid.h:4
int main(int argc, char *argv[])
Definition: CC.cpp:68
NT Reduce(_BinaryOperation __binary_op, NT identity) const
void ParallelReadMM(const std::string &filename, bool onebased, _BinaryOperation BinOp)
Definition: SpParMat.cpp:3417
IT getnrow() const
Definition: SpParMat.cpp:685