COMBINATORIAL_BLAS  1.6
FindSparse.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 PARMAT>
42 void Symmetricize(PARMAT & A)
43 {
44  // boolean addition is practically a "logical or"
45  // therefore this doesn't destruct any links
46  PARMAT AT = A;
47  AT.Transpose();
48  A += AT;
49 }
50 
51 int main(int argc, char* argv[])
52 {
53  int nprocs, myrank;
54  MPI_Init(&argc, &argv);
55  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
56  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
57 
58  if(argc < 3)
59  {
60  if(myrank == 0)
61  {
62  cout << "Usage: ./FindSparse <BASEADDRESS> <Matrix>" << endl;
63  cout << "Input files should be under <BASEADDRESS> in appropriate format" << endl;
64  }
65  MPI_Finalize();
66  return -1;
67  }
68  {
69  string directory(argv[1]);
70  string matrixname(argv[2]);
71  matrixname = directory+"/"+matrixname;
72 
74  shared_ptr<CommGrid> fullWorld;
75  fullWorld.reset( new CommGrid(MPI_COMM_WORLD, 0, 0) );
76 
77  PARDBMAT A(fullWorld); // declare objects
78  FullyDistVec<int,int> crow(fullWorld);
79  FullyDistVec<int,int> ccol(fullWorld);
80  FullyDistVec<int,double> cval(fullWorld);
81 
82  A.ReadDistribute(matrixname, 0);
83 
84  A.Find(crow, ccol, cval);
85  PARDBMAT B(A.getnrow(), A.getncol(), crow, ccol, cval); // Sparse()
86 
87  if (A == B)
88  {
89  SpParHelper::Print("Find and Sparse working correctly\n");
90  }
91  else
92  {
93  SpParHelper::Print("ERROR in Find(), go fix it!\n");
94 
95  SpParHelper::Print("Rows array: \n");
96  crow.DebugPrint();
97 
98  SpParHelper::Print("Columns array: \n");
99  ccol.DebugPrint();
100 
101  SpParHelper::Print("Values array: \n");
102  cval.DebugPrint();
103  }
104 
105  // Begin testing a Matlab like use case
106  // Example provided by Arne De Coninck <arne.deconinck@ugent.be>
107  // X = ones(size(A,2),1)
108  // M = [X' * X, X' * A; A'* X, A' * A]
109 
110  A.PrintInfo();
111  Symmetricize(A);
112  FullyDistVec<int,int> rowsym(fullWorld);
113  FullyDistVec<int,int> colsym(fullWorld);
114  FullyDistVec<int,double> valsym(fullWorld);
115  A.Find(rowsym, colsym, valsym);
116 
117  FullyDistVec<int,double> colsums(fullWorld);
118  A.Reduce(colsums, Column, plus<double>(), 0.0);
119 
120  FullyDistVec<int,double> numcols(A.getcommgrid(), 1, A.getncol());
121 
122  #if defined(COMBBLAS_TR1) || defined(COMBBLAS_BOOST) || defined(NOTGNU)
123  vector< FullyDistVec<int,double> > vals2concat;
124  vals2concat.push_back(numcols);
125  vals2concat.push_back(colsums);
126  vals2concat.push_back(colsums);
127  vals2concat.push_back(valsym);
128  #else
129  vector< FullyDistVec<int,double> > vals2concat{numcols, colsums, colsums, valsym};
130  #endif
131  FullyDistVec<int,double> nval = Concatenate(vals2concat);
132  nval.PrintInfo("Values:");
133 
134  // sparse() expects a zero based index
135  FullyDistVec<int,int> firstrowrids(A.getcommgrid(), A.getncol()+1, 0); // M(1,:)
136  FullyDistVec<int,int> firstcolrids(A.getcommgrid(), A.getncol(), 0); // M(2:end,1)
137  firstcolrids.iota(A.getncol(),1); // fill M(2:end,1)'s row ids
138  rowsym.Apply(bind2nd(plus<int>(), 1));
139 
140  #if defined(COMBBLAS_TR1) || defined(COMBBLAS_BOOST) || defined(NOTGNU)
141  vector< FullyDistVec<int,int> > rows2concat;
142  rows2concat.push_back(firstrowrids);
143  rows2concat.push_back(firstcolrids);
144  rows2concat.push_back(rowsym);
145  #else
146  vector< FullyDistVec<int,int> > rows2concat{firstrowrids, firstcolrids, rowsym}; // C++11 style
147  #endif
148  FullyDistVec<int,int> nrow = Concatenate(rows2concat);
149  nrow.PrintInfo("Row ids:");
150 
151  FullyDistVec<int,int> firstrowcids(A.getcommgrid(), A.getncol()+1, 0); // M(1,:)
152  firstrowcids.iota(A.getncol()+1,0); // fill M(1,:)'s column ids
153  FullyDistVec<int,int> firstcolcids(A.getcommgrid(), A.getncol(), 0); // M(2:end,1)
154  colsym.Apply(bind2nd(plus<int>(), 1));
155 
156  #if defined(COMBBLAS_TR1) || defined(COMBBLAS_BOOST) || defined(NOTGNU)
157  vector< FullyDistVec<int,int> > cols2concat;
158  cols2concat.push_back(firstrowcids);
159  cols2concat.push_back(firstcolcids);
160  cols2concat.push_back(colsym);
161  #else
162  vector< FullyDistVec<int,int> > cols2concat{firstrowcids, firstcolcids, colsym}; // C++11 style
163  #endif
164  FullyDistVec<int,int> ncol = Concatenate(cols2concat);
165  ncol.PrintInfo("Column ids:");
166 
167  PARDBMAT M(A.getnrow()+1, A.getncol()+1, nrow, ncol, nval); // Sparse()
168  M.PrintInfo();
169 
170  // End Arne's test case
171  }
172  MPI_Finalize();
173  return 0;
174 }
double B
FullyDistVec< IT, NT > Concatenate(std::vector< FullyDistVec< IT, NT > > &vecs)
Definition: ParFriends.h:59
double A
void iota(IT globalsize, NT first)
void PrintInfo(std::string vectorname) const
void Apply(_UnaryOperation __unary_op)
Definition: FullyDistVec.h:182
Definition: CCGrid.h:4
void Symmetricize(PARMAT &A)
Definition: FindSparse.cpp:42
int main(int argc, char *argv[])
Definition: FindSparse.cpp:51