COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
Utility.h
Go to the documentation of this file.
1#ifndef BP_UTILITY_H
2#define BP_UTILITY_H
3
4#include "CombBLAS/CombBLAS.h"
5
6namespace combblas {
7
8/*
9 Remove isolated vertices and purmute
10 */
11template <typename PARMAT>
13{
14
15 int nprocs, myrank;
18
19
22 FullyDistVec<int64_t, int64_t> nonisoRowV; // id's of non-isolated (connected) Row vertices
23 FullyDistVec<int64_t, int64_t> nonisoColV; // id's of non-isolated (connected) Col vertices
24 FullyDistVec<int64_t, int64_t> nonisov; // id's of non-isolated (connected) vertices
25
26 A.Reduce(*ColSums, Column, std::plus<int64_t>(), static_cast<int64_t>(0));
27 A.Reduce(*RowSums, Row, std::plus<int64_t>(), static_cast<int64_t>(0));
28
29 // this steps for general graph
30 /*
31 ColSums->EWiseApply(*RowSums, plus<int64_t>()); not needed for bipartite graph
32 nonisov = ColSums->FindInds(bind2nd(greater<int64_t>(), 0));
33 nonisov.RandPerm(); // so that A(v,v) is load-balanced (both memory and time wise)
34 A.operator()(nonisov, nonisov, true); // in-place permute to save memory
35 */
36
37 // this steps for bipartite graph
38 nonisoColV = ColSums->FindInds(bind2nd(std::greater<int64_t>(), 0));
39 nonisoRowV = RowSums->FindInds(bind2nd(std::greater<int64_t>(), 0));
40 delete ColSums;
41 delete RowSums;
42
43
44 {
45 nonisoColV.RandPerm();
46 nonisoRowV.RandPerm();
47 }
48
49
50 int64_t nrows1=A.getnrow(), ncols1=A.getncol(), nnz1 = A.getnnz();
51 double avgDeg1 = (double) nnz1/(nrows1+ncols1);
52
53
54 A.operator()(nonisoRowV, nonisoColV, true);
55
56 int64_t nrows2=A.getnrow(), ncols2=A.getncol(), nnz2 = A.getnnz();
57 double avgDeg2 = (double) nnz2/(nrows2+ncols2);
58
59
60 if(myrank == 0)
61 {
62 std::cout << "ncol nrows nedges deg \n";
63 std::cout << nrows1 << " " << ncols1 << " " << nnz1 << " " << avgDeg1 << " \n";
64 std::cout << nrows2 << " " << ncols2 << " " << nnz2 << " " << avgDeg2 << " \n";
65 }
66
68
69
70}
71
72
73
74
75/*
76 * Serial: Check the validity of the matching solution;
77 we need a better solution using invert
78 */
79template <class IT, class NT>
81{
82
83 int myrank;
85 for(int i=0; i< mateRow2Col.glen ; i++)
86 {
87 int t = mateRow2Col[i];
88
89 if(t!=-1 && mateCol2Row[t]!=i)
90 {
91 if(myrank == 0)
92 std::cout << "Does not satisfy the matching constraints\n";
93 return false;
94 }
95 }
96
97 for(int i=0; i< mateCol2Row.glen ; i++)
98 {
99 int t = mateCol2Row[i];
100 if(t!=-1 && mateRow2Col[t]!=i)
101 {
102 if(myrank == 0)
103 std::cout << "Does not satisfy the matching constraints\n";
104 return false;
105 }
106 }
107 return true;
108}
109
110
111
112template <class IT>
114{
115
116 int myrank;
118 int64_t nrow = mateRow2Col.TotalLength();
119 int64_t ncol = mateCol2Row.TotalLength();
124
125
126
127 bool isMatching = false;
129 isMatching = true;
130 else
131 {
132 SpParHelper::Print("Warning: This is not a matching! Need to check the correctness of the matching (HWPM) code\n");
133 }
134
135 bool isPerfectMatching = false;
136 if((mateRow2ColSparse.getnnz()==nrow) && (mateCol2RowSparse.getnnz() == ncol))
137 isPerfectMatching = true;
138 return isPerfectMatching;
139}
140
141
142// Gievn a matrix and matching vectors, returns the weight of the matching
143template <class IT, class NT, class DER>
145{
146
147 auto commGrid = A.getcommgrid();
148 int myrank=commGrid->GetRank();
149 MPI_Comm World = commGrid->GetWorld();
150 MPI_Comm ColWorld = commGrid->GetColWorld();
151 MPI_Comm RowWorld = commGrid->GetRowWorld();
152 int nprocs = commGrid->GetSize();
153 int pr = commGrid->GetGridRows();
154 int pc = commGrid->GetGridCols();
155 int rowrank = commGrid->GetRankInProcRow();
156 int colrank = commGrid->GetRankInProcCol();
157 int diagneigh = commGrid->GetComplementRank();
158
159 //Information about the matrix distribution
160 //Assume that A is an nrow x ncol matrix
161 //The local submatrix is an lnrow x lncol matrix
162 IT nrows = A.getnrow();
163 IT ncols = A.getncol();
164 IT m_perproc = nrows / pr;
165 IT n_perproc = ncols / pc;
166 DER* spSeq = A.seqptr(); // local submatrix
167 Dcsc<IT, NT>* dcsc = spSeq->GetDCSC();
168 IT lnrow = spSeq->getnrow();
169 IT lncol = spSeq->getncol();
170 IT localRowStart = colrank * m_perproc; // first row in this process
171 IT localColStart = rowrank * n_perproc; // first col in this process
172
173 // -----------------------------------------------------------
174 // replicate mate vectors for mateCol2Row
175 // -----------------------------------------------------------
176 int xsize = (int) mateCol2Row.LocArrSize();
177 int trxsize = 0;
179 MPI_Sendrecv(&xsize, 1, MPI_INT, diagneigh, TRX, &trxsize, 1, MPI_INT, diagneigh, TRX, World, &status);
180 std::vector<IT> trxnums(trxsize);
181 MPI_Sendrecv(mateCol2Row.GetLocArr(), xsize, MPIType<IT>(), diagneigh, TRX, trxnums.data(), trxsize, MPIType<IT>(), diagneigh, TRX, World, &status);
182
183
184 std::vector<int> colsize(pc);
185 colsize[colrank] = trxsize;
187 std::vector<int> dpls(pc,0); // displacements (zero initialized pid)
188 std::partial_sum(colsize.data(), colsize.data()+pc-1, dpls.data()+1);
189 int accsize = std::accumulate(colsize.data(), colsize.data()+pc, 0);
190 std::vector<IT> RepMateC2R(accsize);
191 MPI_Allgatherv(trxnums.data(), trxsize, MPIType<IT>(), RepMateC2R.data(), colsize.data(), dpls.data(), MPIType<IT>(), ColWorld);
192 // -----------------------------------------------------------
193
194 /*
195 if(myrank==1)
196 {
197 for(int i=0; i<RepMateC2R.size(); i++)
198 cout << RepMateC2R[i] << ",";
199 }*/
200
201 NT w = 0;
202 for(auto colit = spSeq->begcol(); colit != spSeq->endcol(); ++colit) // iterate over columns
203 {
204 IT lj = colit.colid(); // local numbering
205 IT mj = RepMateC2R[lj]; // mate of j
206 if(mj >= localRowStart && mj < (localRowStart+lnrow) )
207 {
208 for(auto nzit = spSeq->begnz(colit); nzit < spSeq->endnz(colit); ++nzit)
209 {
210 IT li = nzit.rowid();
211 IT i = li + localRowStart;
212 // TODO: use binary search to directly go to mj-th entry if more than 32 nonzero in this column
213 if( i == mj)
214 {
215 w += nzit.value();
216 //cout << myrank<< ":: row: " << i << " column: "<< lj+localColStart << " weight: " << nzit.value() << endl;
217 }
218 }
219 }
220
221 }
222
225 //MPI_Allreduce(&w, &gw, 1, MPIType<NT>(), MPI_SUM, World);
226 //MPI_Reduce(&w, &gw, 1, MPIType<NT>(), MPI_SUM, 0, World);
227 //MPI_Allreduce(&w, &gw, 1, MPI_DOUBLE, MPI_SUM, World);
228 //cout << myrank << ": " << gw << endl;
229 return w;
230}
231
232
233
234
235
236
237template <typename PARMAT>
238void Symmetricize(PARMAT & A)
239{
240 // boolean addition is practically a "logical or"
241 // therefore this doesn't destruct any links
242 PARMAT AT = A;
243 AT.Transpose();
244 A += AT;
245}
246
247}
248
249#endif
250
int64_t IT
double NT
static void Print(const std::string &s)
int nprocs
Definition comms.cpp:55
long int64_t
Definition compat.h:21
#define TRX
Definition SpDefs.h:103
@ Column
Definition SpDefs.h:115
bool CheckMatching(FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row)
Definition Utility.h:113
void removeIsolated(PARMAT &A)
Definition Utility.h:12
NT MatchingWeight(std::vector< NT > &RepMateWC2R, MPI_Comm RowWorld, NT &minw)
void Symmetricize(PARMAT &A)
Definition ReadMatDist.h:29
bool isMatching(FullyDistVec< IT, NT > &mateCol2Row, FullyDistVec< IT, NT > &mateRow2Col)
Definition Utility.h:80
double A