COMBINATORIAL_BLAS  1.6
IndexingTest.cpp
Go to the documentation of this file.
1 /****************************************************************/
2 /* Parallel Combinatorial BLAS Library (for Graph Computations) */
3 /* version 1.5 -------------------------------------------------*/
4 /* date: 10/09/2015 ---------------------------------------------*/
5 /* authors: Ariful Azad, Aydin Buluc, Adam Lugowski ------------*/
6 /****************************************************************/
7 /*
8  Copyright (c) 2010-2015, 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 #include <mpi.h>
30 #include <sys/time.h>
31 #include <iostream>
32 #include <functional>
33 #include <algorithm>
34 #include <vector>
35 #include <sstream>
36 #include "CombBLAS/CombBLAS.h"
37 
38 using namespace std;
39 using namespace combblas;
40 
41 template <typename IT, typename NT>
42 pair< FullyDistVec<IT,IT>, FullyDistVec<IT,NT> > TopK(FullyDistSpVec<IT,NT> & v, IT k)
43 {
44  // FullyDistVec::FullyDistVec(shared_ptr<CommGrid> commgrid, IT glen, NT initval)
45  FullyDistVec<IT,IT> sel(v.getcommgrid(), k, 0);
46 
47  //void FullyDistVec::iota(IT globalsize, NT first)
48  sel.iota(k, v.TotalLength() - k);
49 
50  FullyDistSpVec<IT,NT> sorted(v);
51  FullyDistSpVec<IT,IT> perm = sorted.sort();
52 
53  // FullyDistVec FullyDistSpVec::operator(FullyDistVec & v)
54  FullyDistVec<IT,IT> topkind = perm(sel);
55  FullyDistVec<IT,NT> topkele = v(topkind);
56  return make_pair(topkind, topkele);
57 }
58 
59 
60 int main(int argc, char* argv[])
61 {
62  int nprocs, myrank;
63  MPI_Init(&argc, &argv);
64  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
65  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
66 
67  if(argc < 6)
68  {
69  if(myrank == 0)
70  {
71  cout << "Usage: ./IndexingTest <BASEADDRESS> <Matrix> <IndexedMatrix> <VectorOne> <VectorTwo>" << endl;
72  cout << "Example: ./IndexingTest ../mfiles B_100x100.txt B_10x30_Indexed.txt rand10outta100.txt rand30outta100.txt" << endl;
73  cout << "Input files should be under <BASEADDRESS> in tuples format" << endl;
74  }
75  MPI_Finalize();
76  return -1;
77  }
78  {
79  string directory(argv[1]);
80  string normalname(argv[2]);
81  string indexdname(argv[3]);
82  string vec1name(argv[4]);
83  string vec2name(argv[5]);
84  normalname = directory+"/"+normalname;
85  indexdname = directory+"/"+indexdname;
86  vec1name = directory+"/"+vec1name;
87  vec2name = directory+"/"+vec2name;
88 
89  ifstream inputvec1(vec1name.c_str());
90  ifstream inputvec2(vec2name.c_str());
91 
92  if(myrank == 0)
93  {
94  if(inputvec1.fail() || inputvec2.fail())
95  {
96  cout << "One of the input vector files do not exist, aborting" << endl;
97  MPI_Abort(MPI_COMM_WORLD, NOFILE);
98  return -1;
99  }
100  }
101  MPI_Barrier(MPI_COMM_WORLD);
103  shared_ptr<CommGrid> fullWorld;
104  fullWorld.reset( new CommGrid(MPI_COMM_WORLD, 0, 0) );
105 
106  PARDBMAT A(fullWorld);
107  PARDBMAT AID(fullWorld);
108  PARDBMAT ACID(fullWorld);
109  FullyDistVec<int,int> vec1(fullWorld);
110  FullyDistVec<int,int> vec2(fullWorld);
111 
112  A.ReadDistribute(normalname, 0);
113  AID.ReadDistribute(indexdname, 0);
114  vec1.ReadDistribute(inputvec1, 0);
115  vec2.ReadDistribute(inputvec2, 0);
116 
117  vec1.Apply(bind2nd(minus<int>(), 1)); // For 0-based indexing
118  vec2.Apply(bind2nd(minus<int>(), 1));
119  ACID = A(vec1, vec2);
120 
121  if (ACID == AID)
122  {
123  SpParHelper::Print("Indexing working correctly\n");
124  }
125  else
126  {
127  SpParHelper::Print("ERROR in indexing, go fix it!\n");
128  }
129 
130  FullyDistVec<int,int> crow(fullWorld);
131  FullyDistVec<int,int> ccol(fullWorld);
132  FullyDistVec<int,double> cval(fullWorld);
133  A.Find(crow, ccol, cval);
134  FullyDistSpVec<int, double> sval = cval;
135  sval.DebugPrint();
136 
137  pair< FullyDistVec<int,int> , FullyDistVec<int,double> > ptopk;
138  ptopk = TopK(sval, 3);
139  //ptopk.first.DebugPrint();
140  //ptopk.second.DebugPrint();
141 
142  // generate random permutations
143  FullyDistVec<int,int> p(fullWorld);
144  FullyDistVec<int,int> q(fullWorld);
145  p.iota(A.getnrow(), 0);
146  q.iota(A.getncol(), 0);
147  p.RandPerm();
148  q.RandPerm();
149 
150  PARDBMAT B = A(p,q);
151  A.PrintInfo();
152  B.PrintInfo();
153 
154  float oldbalance = A.LoadImbalance();
155  float newbalance = B.LoadImbalance();
156  ostringstream outs;
157  outs << "Old balance: " << oldbalance << endl;
158  outs << "New balance: " << newbalance << endl;
159  SpParHelper::Print(outs.str());
160  SpParHelper::Print(outs.str());
161 
162  // B = P A Q
163  // get the inverse permutations
164  FullyDistVec<int, int> pinv = p.sort();
165  FullyDistVec<int, int> qinv = q.sort();
166  SpParHelper::Print("Sorts are done\n");
167  PARDBMAT C = B(pinv,qinv);
168  if(C == A)
169  {
170  SpParHelper::Print("Double permutation successfully restored the original\n");
171  }
172  else
173  {
174  SpParHelper::Print("Error in permutation\n");
175  }
176 
177  inputvec1.clear();
178  inputvec1.close();
179  inputvec2.clear();
180  inputvec2.close();
181  }
182  MPI_Finalize();
183  return 0;
184 }
double B
int main(int argc, char *argv[])
#define NOFILE
Definition: SpDefs.h:75
pair< FullyDistVec< IT, IT >, FullyDistVec< IT, NT > > TopK(FullyDistSpVec< IT, NT > &v, IT k)
std::shared_ptr< CommGrid > getcommgrid() const
double A
void iota(IT globalsize, NT first)
FullyDistVec< IT, IT > sort()
double C
void Apply(_UnaryOperation __unary_op)
Definition: FullyDistVec.h:182
Definition: CCGrid.h:4
std::ifstream & ReadDistribute(std::ifstream &infile, int master, HANDLER handler)
FullyDistSpVec< IT, IT > sort()
sort the vector itself, return the permutation vector (0-based)