33#ifndef __STDC_CONSTANT_MACROS
34#define __STDC_CONSTANT_MACROS
36#ifndef __STDC_LIMIT_MACROS
37#define __STDC_LIMIT_MACROS
48#include "CombBLAS/CombBLAS.h"
118 param.ifilename =
"";
119 param.isInputMM =
false;
120 param.ofilename =
"";
126 param.remove_isolated =
false;
127 param.randpermute = 0;
130 param.inflation = 0.0;
133 param.prunelimit = 1.0/10000.0;
135 param.recover_num = 1400;
136 param.recover_pct = .9;
137 param.kselectVersion = 1;
141 param.perProcessMem = 0;
142 param.isDoublePrecision =
true;
143 param.is64bInt =
true;
152 runinfo <<
"\n======================================" <<
endl;
153 runinfo <<
"Running HipMCL with the parameters: " <<
endl;
154 runinfo <<
"======================================" <<
endl;
157 runinfo <<
" input file type: " ;
168 runinfo <<
" Remove isolated vertices? : ";
172 runinfo <<
" Randomly permute vertices? : ";
193 runinfo <<
" Memory avilable per process: ";
196 if(
param.isDoublePrecision)
runinfo <<
"Using double precision floating point" <<
endl;
197 else runinfo <<
"Using single precision floating point" <<
endl;
202 runinfo <<
" Show matrices after major steps? : ";
205 runinfo <<
"======================================" <<
endl;
211 for (
int i = 1; i <
argc; i++)
216 else if (
strcmp(
argv[i],
"--matrix-market")==0){
217 param.isInputMM =
true;
225 else if (
strcmp(
argv[i],
"--remove-isolated")==0){
226 param.remove_isolated =
true;
228 else if (
strcmp(
argv[i],
"--tournament-select")==0){
229 param.kselectVersion = 1;
231 else if (
strcmp(
argv[i],
"--quick-select")==0){
232 param.kselectVersion = 2;
250 if(
param.recover_pct>1)
param.recover_pct/=100.00;
260 else if (
strcmp(
argv[i],
"-per-process-mem")==0) {
263 else if (
strcmp(
argv[i],
"--single-precision")==0) {
264 param.isDoublePrecision =
false;
266 else if (
strcmp(
argv[i],
"--32bit-index")==0) {
267 param.is64bInt =
false;
271 if(
param.ofilename==
"")
273 param.ofilename =
param.ifilename +
".hipmcl";
283 runinfo <<
"Usage: ./hipmcl -M <input filename> -I <inlfation> (required)" <<
endl;
285 runinfo <<
"======================================" <<
endl;
287 runinfo <<
"======================================" <<
endl;
292 runinfo <<
" -M <input file name (labeled triples format)> (mandatory)" <<
endl;
293 runinfo <<
" --matrix-market : if provided, the input file is in the matrix market format (default: the file is in labeled triples format)" <<
endl;
294 runinfo <<
" -base <index of the first vertex in the matrix market file, 0|1> (default: 1) " <<
endl;
295 runinfo <<
" -o <output filename> (default: input_file_name.hipmcl )" <<
endl;
298 runinfo <<
"-I <inflation> (mandatory)\n";
301 runinfo <<
" -rand <randomly permute vertices> (default:0)\n";
302 runinfo <<
" --remove-isolated : if provided, remove isolated vertices (default: don't remove isolated vertices)\n";
306 runinfo <<
" -p <cutoff> (default: 1/10000)\n";
307 runinfo <<
" -R <recovery number> (default: 1400)\n";
308 runinfo <<
" -pct <recovery pct> (default: 90)\n";
309 runinfo <<
" -S <selection number> (default: 1100)\n";
313 runinfo <<
" -phases <number of phases> (default:1)\n";
314 runinfo <<
" -per-process-mem <memory (GB) available per process> (default:0, number of phases is not estimated)\n";
315 runinfo <<
" --single-precision (if not provided, use double precision floating point numbers)\n" <<
endl;
316 runinfo <<
" --32bit-index (if not provided, use 64 bit indexing for vertex ids)\n" <<
endl;
319 runinfo <<
" --show: show information about matrices after major steps (default: do not show matrices)" <<
endl;
323 runinfo <<
"======================================" <<
endl;
325 runinfo <<
"======================================" <<
endl;
326 runinfo <<
"Example with with a graph in labeled triples format on a laptop with 8GB memory and 8 cores:\nexport OMP_NUM_THREADS=8\nbin/hipmcl -M data/sevenvertexgraph.txt -I 2 -per-process-mem 8" <<
endl;
327 runinfo <<
"Same as above with 4 processes and 2 theaded per process cores:\nexport OMP_NUM_THREADS=2\nmpirun -np 4 bin/hipmcl -M data/sevenvertexgraph.txt -I 2 -per-process-mem 2" <<
endl;
328 runinfo <<
"Example with a graph in matrix market format:\nbin/hipmcl -M data/sevenvertex.mtx --matrix-market -base 1 -I 2 -per-process-mem 8" <<
endl;
330 runinfo <<
"Example on the NERSC/Edison system with 16 nodes and 24 threads per node: \nsrun -N 16 -n 16 -c 24 bin/hipmcl -M data/hep-th.mtx --matrix-market -base 1 -per-process-mem 64 -o hep-th.hipmcl" <<
endl;
337template <
typename IT,
typename NT,
typename DER>
352template <
typename IT,
typename NT,
typename DER>
360template <
typename IT,
typename NT,
typename DER>
376template <
typename IT,
typename NT,
typename DER>
385template <
typename IT,
typename NT,
typename DER>
398template <
typename IT,
typename NT,
typename DER>
418template <
typename IT,
typename NT,
typename DER>
422 if(
A.getnrow() ==
A.getncol())
425 p.iota(
A.getnrow(), 0);
436template <
typename IT,
typename NT,
typename DER>
439 if(
param.remove_isolated)
442 if(
param.randpermute)
470 A =
MemEfficientSpGEMM<PTFF, NT, DER>(
A,
A,
param.phases,
param.prunelimit, (
IT)
param.select, (
IT)
param.recover_num,
param.recover_pct,
param.kselectVersion,
param.perProcessMem);
520 cout <<
"================detailed timing==================" <<
endl;
533 cout <<
"=================================================" <<
endl;
543template <
typename IT,
typename NT,
typename DER>
555template <
typename GIT,
typename LIT,
typename NT>
584 GIT nnz =
A.getnnz();
586 outs <<
"Number of vertices: " <<
nv <<
" number of edges: "<< nnz <<
endl;
629 s2 <<
"=================================================\n" <<
endl ;
641 printf(
"ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
666 if(
param.ifilename==
"" ||
param.inflation == 0.0)
677 cout <<
"\nProcess Grid used (pr x pc x threads): " <<
sqrt(nprocs) <<
" x " <<
sqrt(nprocs) <<
" x " <<
nthreads <<
endl;
684 if(
param.perProcessMem==0)
688 cout <<
"******** Number of phases will not be estimated as -per-process-mem option is supplied. It is highly recommended that you provide -per-process-mem option for large-scale runs. *********** " <<
endl;
693 if(
param.isDoublePrecision)
698 else if(
param.is64bInt)
double mcl_localspgemmtime
double mcl_multiwaymergetime
double mcl_prunecolumntime
int main(int argc, char *argv[])
void ShowParam(HipMCLParam ¶m)
void RandPermute(SpParMat< IT, NT, DER > &A, HipMCLParam ¶m)
void MakeColStochastic(SpParMat< IT, NT, DER > &A)
void ProcessParam(int argc, char *argv[], HipMCLParam ¶m)
void InitParam(HipMCLParam ¶m)
FullyDistVec< IT, IT > HipMCL(SpParMat< IT, NT, DER > &A, HipMCLParam ¶m)
double cblas_allgathertime
void AdjustLoops(SpParMat< IT, NT, DER > &A)
void RemoveIsolated(SpParMat< IT, NT, DER > &A, HipMCLParam ¶m)
void Inflate(SpParMat< IT, NT, DER > &A, double power)
NT Chaos(SpParMat< IT, NT, DER > &A)
void MainBody(HipMCLParam ¶m)
FullyDistVec< IT, IT > Interpret(SpParMat< IT, NT, DER > &A)
FullyDistVec< int64_t, double > MPI_DenseVec
SpDCCols< int64_t, double > DCCols
SpParMat< int64_t, double, DCCols > MPI_DCCols
static void Print(const std::string &s)
void WriteMCLClusters(std::string ofName, FullyDistVec< IT, IT > clustIdForVtx, FullyDistVec< IT, std::array< char, MAXVERTNAME > > vtxLabels)
void Symmetricize(PARMAT &A)
FullyDistVec< IT, IT > CC(SpParMat< IT, NT, DER > &A, IT &nCC)