COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
FastSV.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 <chrono>
49#include "CombBLAS/CombBLAS.h"
50#include "FastSV.h"
51
52using namespace std;
53using namespace combblas;
54
59using Int = int64_t;
60
61class Dist
62{
63public:
66};
67
68
69
70int main(int argc, char* argv[])
71{
72 int provided;
75 {
76 printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
78 }
79
80 int nthreads = 1;
81#ifdef THREADED
82#pragma omp parallel
83 {
84 nthreads = omp_get_num_threads();
85 }
86#endif
87
88 int nprocs, myrank;
91 if(myrank == 0)
92 {
93 cout << "Process Grid (p x p x t): " << sqrt(nprocs) << " x " << sqrt(nprocs) << " x " << nthreads << endl;
94 }
95
96 if(argc < 3)
97 {
98 if(myrank == 0)
99 {
100 cout << "Usage: ./fastsv -I <mm|triples> -M <FILENAME> (required)\n";
101 cout << "-I <INPUT FILE TYPE> (mm: matrix market, triples: (vtx1, vtx2, edge_weight) triples. default:mm)\n";
102 cout << "-base <BASE OF MATRIX MARKET> (default:1)\n";
103 cout << "-rand <RANDOMLY PERMUTE VERTICES> (default:0)\n";
104 cout << "Example (0-indexed mtx with random permutation): ./fastsv -M input.mtx -base 0 -rand 1" << endl;
105 cout << "Example (triples format): ./fastsv -I triples -M input.txt" << endl;
106 }
107 MPI_Finalize();
108 return -1;
109 }
110 {
111 string ifilename = "";
112 int base = 1;
113 int randpermute = 0;
114 int step = 1;
115 bool isMatrixMarket = true;
116
117 for (int i = 1; i < argc; i++)
118 {
119 if (strcmp(argv[i],"-I")==0)
120 {
121 string ifiletype = string(argv[i+1]);
122 if(ifiletype == "triples") isMatrixMarket = false;
123 }
124 if (strcmp(argv[i],"-M")==0)
125 {
126 ifilename = string(argv[i+1]);
127 if(myrank == 0) printf("filename: %s\n",ifilename.c_str());
128 }
129 else if (strcmp(argv[i],"-base")==0)
130 {
131 base = atoi(argv[i + 1]);
132 if(myrank == 0) printf("\nBase of MM (1 or 0):%d\n",base);
133 }
134 else if (strcmp(argv[i],"-rand")==0)
135 {
136 randpermute = atoi(argv[i + 1]);
137 if(myrank == 0) printf("\nRandomly permute the matrix? (1 or 0):%d\n",randpermute);
138 }
139 else if (strcmp(argv[i],"-step")==0)
140 {
141 step = atoi(argv[i + 1]);
142 if(myrank == 0) printf("\nGraph size:%.2f%%\n", 100.0 / step);
143 }
144 }
145
146 double tIO = MPI_Wtime();
147 Dist::MPI_DCCols A; // construct object
148
150 A.ParallelReadMM(ifilename, base, maximum<bool>());
151 else
152 A.ReadGeneralizedTuples(ifilename, maximum<bool>());
153 A.PrintInfo();
154
156 AT.Transpose();
157 if(!(AT == A))
158 {
159 SpParHelper::Print("Symmatricizing an unsymmetric input matrix.\n");
160 A += AT;
161 }
162 A.PrintInfo();
163
165 outs << "File Read time: " << MPI_Wtime() - tIO << endl;
167
168 if(randpermute && isMatrixMarket) // no need for GeneralizedTuples
169 {
170 double tPerm = MPI_Wtime();
171 // randomly permute for load balance
172 if(A.getnrow() == A.getncol())
173 {
174 FullyDistVec<Int, Int> p( A.getcommgrid());
175 p.iota(A.getnrow(), 0);
176 p.RandPerm();
177 (A)(p,p,true);// in-place permute to save memory
178 SpParHelper::Print("Applied symmetric permutation.\n");
179 }
180 else
181 {
182 SpParHelper::Print("Rectangular matrix: Can not apply symmetric permutation.\n");
183 }
184 double dur = MPI_Wtime() - tPerm;
185 outs.str("");
186 outs.clear();
187 outs << "Random permutation time: " << dur << endl;
189 }
190
193 outs.str("");
194 outs.clear();
195 outs << "isolated vertice: " << isov.TotalLength() << endl;
197
198 float balance = A.LoadImbalance();
199 Int nnz = A.getnnz();
200 outs.str("");
201 outs.clear();
202 outs << "Load balance: " << balance << endl;
203 outs << "Nonzeros: " << nnz << endl;
205 double t1 = MPI_Wtime();
206
207 Int nCC = 0;
209
210 double t2 = MPI_Wtime();
211 //outs.str("");
212 //outs.clear();
213 //outs << "Total time: " << t2 - t1 << endl;
214 //SpParHelper::Print(outs.str());
215 //HistCC(cclabels, nCC);
216
217 Int nclusters = cclabels.Reduce(maximum<Int>(), (Int) 0 ) ;
218 nclusters ++; // because of zero based indexing for clusters
219
220 double tend = MPI_Wtime();
222 s2 << "Number of components: " << nclusters << endl;
223 s2 << "Total time: " << (t2-t1) << endl;
224 s2 << "=================================================\n" << endl ;
225 SpParHelper::Print(s2.str());
226
227 }
228
229 MPI_Finalize();
230 return 0;
231}
int main()
Definition Driver.cpp:12
double tIO
Definition MCL.cpp:86
SpParMat< Int, Int, DCCols > MPI_DCCols
Definition FastSV.cpp:65
SpDCCols< Int, Int > DCCols
Definition FastSV.cpp:64
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
FullyDistVec< IT, IT > SV(SpParMat< IT, NT, DER > &A, IT &nCC)
Definition FastSV.h:336
double A