COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
MatPermuteSave.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
22
23
24using namespace std;
25using namespace combblas;
26
27
28
29typedef SpParMat < int64_t, bool, SpDCCols<int64_t,bool> > Par_DCSC_Bool;
30typedef SpParMat < int64_t, int64_t, SpDCCols<int64_t, int64_t> > Par_DCSC_int64_t;
31typedef SpParMat < int64_t, double, SpDCCols<int64_t, double> > Par_DCSC_Double;
32
33
34
35int main(int argc, char* argv[])
36{
37 int provided;
38 MPI_Init_thread(&argc, &argv, MPI_THREAD_SERIALIZED, &provided);
39 if (provided < MPI_THREAD_SERIALIZED)
40 {
41 printf("ERROR: The MPI library does not have MPI_THREAD_SERIALIZED support\n");
42 MPI_Abort(MPI_COMM_WORLD, 1);
43 }
44
45 int nprocs, myrank;
46 MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
47 MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
48 if(argc < 2)
49 {
50 if(myrank == 0)
51 {
52 cout << "Usage: ./rcm <filename>" << endl;
53
54 }
55 MPI_Finalize();
56 return -1;
57 }
58 {
59 Par_DCSC_Bool * ABool;
60 ostringstream tinfo;
61
62 ABool = new Par_DCSC_Bool();
63 string filename(argv[1]);
64 tinfo.str("");
65 tinfo << "**** Reading input matrix: " << filename << " ******* " << endl;
66 SpParHelper::Print(tinfo.str());
67 double t01 = MPI_Wtime();
68 ABool->ParallelReadMM(filename, true, maximum<bool>());
69 double t02 = MPI_Wtime();
70 tinfo.str("");
71 tinfo << "Reader took " << t02-t01 << " seconds" << endl;
72 SpParHelper::Print(tinfo.str());
73
74
75 if(ABool->getnrow() == ABool->getncol())
76 {
78 p.iota(ABool->getnrow(), 0);
79 p.RandPerm();
80 (*ABool)(p,p,true);// in-place permute to save memory
81 SpParHelper::Print("Applied symmetric permutation.\n");
82 }
83 else
84 {
85 SpParHelper::Print("Rectangular matrix: Can not apply symmetric permutation.\n");
86 }
87
88 filename += "_permuted";
89 ABool->SaveGathered(filename);
90 tinfo.str("");
91 tinfo << "**** Saved to output matrix: " << filename << " ******* " << endl;
92 SpParHelper::Print(tinfo.str());
93 delete ABool;
94
95 }
96 MPI_Finalize();
97 return 0;
98}
99
int main(int argc, char *argv[])
SpParMat< int64_t, bool, SpDCCols< int64_t, bool > > Par_DCSC_Bool
SpParMat< int64_t, double, SpDCCols< int64_t, double > > Par_DCSC_Double
SpParMat< int64_t, int64_t, SpDCCols< int64_t, int64_t > > Par_DCSC_int64_t
void iota(IT globalsize, NT first)
static void Print(const std::string &s)
IT getncol() const
Definition SpParMat.cpp:694
void SaveGathered(std::string filename, HANDLER handler, bool transpose=false) const
void ParallelReadMM(const std::string &filename, bool onebased, _BinaryOperation BinOp)
std::shared_ptr< CommGrid > getcommgrid() const
Definition SpParMat.h:275
IT getnrow() const
Definition SpParMat.cpp:685
Compute the maximum of two values.
Definition Operations.h:155