COMBINATORIAL_BLAS  1.6
Utility.h
Go to the documentation of this file.
1 #ifndef BP_UTILITY_H
2 #define BP_UTILITY_H
3 
4 #include "../CombBLAS.h"
5 
6 namespace combblas {
7 
8 /*
9  Remove isolated vertices and purmute
10  */
11 template <typename PARMAT>
12 void removeIsolated(PARMAT & A)
13 {
14 
15  int nprocs, myrank;
16  MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
17  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
18 
19 
20  FullyDistVec<int64_t, int64_t> * ColSums = new FullyDistVec<int64_t, int64_t>(A.getcommgrid());
21  FullyDistVec<int64_t, int64_t> * RowSums = new FullyDistVec<int64_t, int64_t>(A.getcommgrid());
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 
67  MPI_Barrier(MPI_COMM_WORLD);
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  */
79 template <class IT, class NT>
80 bool isMatching(FullyDistVec<IT,NT> & mateCol2Row, FullyDistVec<IT,NT> & mateRow2Col)
81 {
82 
83  int myrank;
84  MPI_Comm_rank(MPI_COMM_WORLD,&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 
112 template <class IT>
113 bool CheckMatching(FullyDistVec<IT,IT> & mateRow2Col, FullyDistVec<IT,IT> & mateCol2Row)
114 {
115 
116  int myrank;
117  MPI_Comm_rank(MPI_COMM_WORLD,&myrank);
118  int64_t nrow = mateRow2Col.TotalLength();
119  int64_t ncol = mateCol2Row.TotalLength();
120  FullyDistSpVec<IT,IT> mateRow2ColSparse (mateRow2Col, [](IT mate){return mate!=-1;});
121  FullyDistSpVec<IT,IT> mateCol2RowSparse (mateCol2Row, [](IT mate){return mate!=-1;});
122  FullyDistSpVec<IT,IT> mateRow2ColInverted = mateRow2ColSparse.Invert(ncol);
123  FullyDistSpVec<IT,IT> mateCol2RowInverted = mateCol2RowSparse.Invert(nrow);
124 
125 
126 
127  bool isMatching = false;
128  if((mateCol2RowSparse == mateRow2ColInverted) && (mateRow2ColSparse == mateCol2RowInverted))
129  isMatching = true;
130 
131  bool isPerfectMatching = false;
132  if((mateRow2ColSparse.getnnz()==nrow) && (mateCol2RowSparse.getnnz() == ncol))
133  isPerfectMatching = true;
134 
135 
136  if(myrank == 0)
137  {
138  std::cout << "-------------------------------" << std::endl;
139  if(isMatching)
140  {
141 
142  std::cout << "| This is a matching |" << std::endl;
143  if(isPerfectMatching)
144  std::cout << "| This is a perfect matching |" << std::endl;
145 
146 
147  }
148  else
149  std::cout << "| This is not a matching |" << std::endl;
150  std::cout << "-------------------------------" << std::endl;
151  }
152 
153  return isPerfectMatching;
154 
155 }
156 
157 
158 // Gievn a matrix and matching vectors, returns the weight of the matching
159 template <class IT, class NT, class DER>
161 {
162 
163  auto commGrid = A.getcommgrid();
164  int myrank=commGrid->GetRank();
165  MPI_Comm World = commGrid->GetWorld();
166  MPI_Comm ColWorld = commGrid->GetColWorld();
167  MPI_Comm RowWorld = commGrid->GetRowWorld();
168  int nprocs = commGrid->GetSize();
169  int pr = commGrid->GetGridRows();
170  int pc = commGrid->GetGridCols();
171  int rowrank = commGrid->GetRankInProcRow();
172  int colrank = commGrid->GetRankInProcCol();
173  int diagneigh = commGrid->GetComplementRank();
174 
175  //Information about the matrix distribution
176  //Assume that A is an nrow x ncol matrix
177  //The local submatrix is an lnrow x lncol matrix
178  IT nrows = A.getnrow();
179  IT ncols = A.getncol();
180  IT m_perproc = nrows / pr;
181  IT n_perproc = ncols / pc;
182  DER* spSeq = A.seqptr(); // local submatrix
183  Dcsc<IT, NT>* dcsc = spSeq->GetDCSC();
184  IT lnrow = spSeq->getnrow();
185  IT lncol = spSeq->getncol();
186  IT localRowStart = colrank * m_perproc; // first row in this process
187  IT localColStart = rowrank * n_perproc; // first col in this process
188 
189  // -----------------------------------------------------------
190  // replicate mate vectors for mateCol2Row
191  // -----------------------------------------------------------
192  int xsize = (int) mateCol2Row.LocArrSize();
193  int trxsize = 0;
194  MPI_Status status;
195  MPI_Sendrecv(&xsize, 1, MPI_INT, diagneigh, TRX, &trxsize, 1, MPI_INT, diagneigh, TRX, World, &status);
196  std::vector<IT> trxnums(trxsize);
197  MPI_Sendrecv(mateCol2Row.GetLocArr(), xsize, MPIType<IT>(), diagneigh, TRX, trxnums.data(), trxsize, MPIType<IT>(), diagneigh, TRX, World, &status);
198 
199 
200  std::vector<int> colsize(pc);
201  colsize[colrank] = trxsize;
202  MPI_Allgather(MPI_IN_PLACE, 1, MPI_INT, colsize.data(), 1, MPI_INT, ColWorld);
203  std::vector<int> dpls(pc,0); // displacements (zero initialized pid)
204  std::partial_sum(colsize.data(), colsize.data()+pc-1, dpls.data()+1);
205  int accsize = std::accumulate(colsize.data(), colsize.data()+pc, 0);
206  std::vector<IT> RepMateC2R(accsize);
207  MPI_Allgatherv(trxnums.data(), trxsize, MPIType<IT>(), RepMateC2R.data(), colsize.data(), dpls.data(), MPIType<IT>(), ColWorld);
208  // -----------------------------------------------------------
209 
210  /*
211  if(myrank==1)
212  {
213  for(int i=0; i<RepMateC2R.size(); i++)
214  cout << RepMateC2R[i] << ",";
215  }*/
216 
217  NT w = 0;
218  for(auto colit = spSeq->begcol(); colit != spSeq->endcol(); ++colit) // iterate over columns
219  {
220  IT lj = colit.colid(); // local numbering
221  IT mj = RepMateC2R[lj]; // mate of j
222  if(mj >= localRowStart && mj < (localRowStart+lnrow) )
223  {
224  for(auto nzit = spSeq->begnz(colit); nzit < spSeq->endnz(colit); ++nzit)
225  {
226  IT li = nzit.rowid();
227  IT i = li + localRowStart;
228  // TODO: use binary search to directly go to mj-th entry if more than 32 nonzero in this column
229  if( i == mj)
230  {
231  w += nzit.value();
232  //cout << myrank<< ":: row: " << i << " column: "<< lj+localColStart << " weight: " << nzit.value() << endl;
233  }
234  }
235  }
236 
237  }
238 
239  MPI_Barrier(World);
240  MPI_Allreduce(MPI_IN_PLACE, &w, 1, MPIType<NT>(), MPI_SUM, World);
241  //MPI_Allreduce(&w, &gw, 1, MPIType<NT>(), MPI_SUM, World);
242  //MPI_Reduce(&w, &gw, 1, MPIType<NT>(), MPI_SUM, 0, World);
243  //MPI_Allreduce(&w, &gw, 1, MPI_DOUBLE, MPI_SUM, World);
244  //cout << myrank << ": " << gw << endl;
245  return w;
246 }
247 
248 
249 
250 
251 
252 
253 template <typename PARMAT>
254 void Symmetricize(PARMAT & A)
255 {
256  // boolean addition is practically a "logical or"
257  // therefore this doesn't destruct any links
258  PARMAT AT = A;
259  AT.Transpose();
260  A += AT;
261 }
262 
263 }
264 
265 #endif
266 
#define TRX
Definition: SpDefs.h:102
std::shared_ptr< CommGrid > getcommgrid() const
Definition: SpParMat.h:275
FullyDistVec< IT, IT > FindInds(_Predicate pred) const
Return the indices where pred is true.
FullyDistSpVec< IT, NT > Invert(IT globallen)
double A
NT MatchingWeight(std::vector< NT > &RepMateWC2R, MPI_Comm RowWorld, NT &minw)
long int64_t
Definition: compat.h:21
IT getncol() const
Definition: SpParMat.cpp:694
bool isMatching(FullyDistVec< IT, NT > &mateCol2Row, FullyDistVec< IT, NT > &mateRow2Col)
Definition: Utility.h:80
Definition: CCGrid.h:4
bool CheckMatching(FullyDistVec< IT, IT > &mateRow2Col, FullyDistVec< IT, IT > &mateCol2Row)
Definition: Utility.h:113
NT Reduce(_BinaryOperation __binary_op, NT identity) const
void Symmetricize(PARMAT &A)
Definition: ReadMatDist.h:29
void removeIsolated(PARMAT &A)
Definition: Utility.h:12
const NT * GetLocArr() const
Definition: FullyDistVec.h:167
IT getnrow() const
Definition: SpParMat.cpp:685