ColPack
|
00001 /************************************************************************************ 00002 Copyright (C) 2005-2008 Assefaw H. Gebremedhin, Arijit Tarafdar, Duc Nguyen, 00003 Alex Pothen 00004 00005 This file is part of ColPack. 00006 00007 ColPack is free software: you can redistribute it and/or modify 00008 it under the terms of the GNU Lesser General Public License as published 00009 by the Free Software Foundation, either version 3 of the License, or 00010 (at your option) any later version. 00011 00012 ColPack is distributed in the hope that it will be useful, 00013 but WITHOUT ANY WARRANTY; without even the implied warranty of 00014 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00015 GNU Lesser General Public License for more details. 00016 00017 You should have received a copy of the GNU Lesser General Public License 00018 along with ColPack. If not, see <http://www.gnu.org/licenses/>. 00019 ************************************************************************************/ 00020 00021 #include "ColPackHeaders.h" 00022 00023 using namespace std; 00024 00025 namespace ColPack 00026 { 00027 //Private Function 2301 00028 int BipartiteGraphPartialOrdering::CheckVertexOrdering(string s_VertexOrderingVariant) 00029 { 00030 if(m_s_VertexOrderingVariant.compare(s_VertexOrderingVariant) == 0) 00031 { 00032 return(_TRUE); 00033 } 00034 00035 if(m_s_VertexOrderingVariant.compare("ALL") != 0) 00036 { 00037 m_s_VertexOrderingVariant = s_VertexOrderingVariant; 00038 } 00039 00040 return(_FALSE); 00041 } 00042 00043 00044 //Public Constructor 2351 00045 BipartiteGraphPartialOrdering::BipartiteGraphPartialOrdering() 00046 { 00047 Clear(); 00048 } 00049 00050 00051 //Public Destructor 2352 00052 BipartiteGraphPartialOrdering::~BipartiteGraphPartialOrdering() 00053 { 00054 Clear(); 00055 } 00056 00057 00058 //Virtual Function 2353 00059 void BipartiteGraphPartialOrdering::Clear() 00060 { 00061 BipartiteGraphInputOutput::Clear(); 00062 00063 m_d_OrderingTime = _UNKNOWN; 00064 00065 m_s_VertexOrderingVariant.clear(); 00066 00067 m_vi_OrderedVertices.clear(); 00068 00069 return; 00070 } 00071 00072 00073 //Virtual Function 2354 00074 void BipartiteGraphPartialOrdering::Reset() 00075 { 00076 m_d_OrderingTime = _UNKNOWN; 00077 00078 m_s_VertexOrderingVariant.clear(); 00079 00080 m_vi_OrderedVertices.clear(); 00081 00082 return; 00083 } 00084 00085 00086 //Public Function 2355 00087 int BipartiteGraphPartialOrdering::RowNaturalOrdering() 00088 { 00089 if(CheckVertexOrdering("ROW_NATURAL")) 00090 { 00091 return(_TRUE); 00092 } 00093 00094 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00095 00096 m_vi_OrderedVertices.clear(); 00097 m_vi_OrderedVertices.reserve((unsigned) i_LeftVertexCount); 00098 00099 for(int i = 0; i<i_LeftVertexCount; i++) 00100 { 00101 m_vi_OrderedVertices.push_back(i); 00102 } 00103 00104 return(_TRUE); 00105 } 00106 00107 00108 //Public Function 2356 00109 int BipartiteGraphPartialOrdering::ColumnNaturalOrdering() 00110 { 00111 if(CheckVertexOrdering("COLUMN_NATURAL")) 00112 { 00113 return(_TRUE); 00114 } 00115 00116 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00117 int i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size()); 00118 00119 m_vi_OrderedVertices.clear(); 00120 m_vi_OrderedVertices.reserve((unsigned) i_RightVertexCount); 00121 00122 for(int i = 0; i < i_RightVertexCount; i++) 00123 { 00124 m_vi_OrderedVertices.push_back(i + i_LeftVertexCount); 00125 } 00126 00127 return(_TRUE); 00128 } 00129 00130 int BipartiteGraphPartialOrdering::RowRandomOrdering() { 00131 if(CheckVertexOrdering("ROW_RANDOM")) 00132 { 00133 return(_TRUE); 00134 } 00135 00136 m_s_VertexOrderingVariant = "ROW_RANDOM"; 00137 00138 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00139 00140 m_vi_OrderedVertices.clear(); 00141 m_vi_OrderedVertices.resize((unsigned) i_LeftVertexCount); 00142 00143 for(unsigned int i = 0; i<i_LeftVertexCount; i++) { 00144 m_vi_OrderedVertices[i] = i; 00145 } 00146 00147 randomOrdering(m_vi_OrderedVertices); 00148 00149 return(_TRUE); 00150 } 00151 00152 int BipartiteGraphPartialOrdering::ColumnRandomOrdering() { 00153 if(CheckVertexOrdering("COLUMN_RANDOM")) 00154 { 00155 return(_TRUE); 00156 } 00157 00158 m_s_VertexOrderingVariant = "COLUMN_RANDOM"; 00159 00160 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00161 int i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size()); 00162 00163 m_vi_OrderedVertices.clear(); 00164 m_vi_OrderedVertices.resize((unsigned) i_RightVertexCount); 00165 00166 for(unsigned int i = 0; i<i_RightVertexCount; i++) { 00167 m_vi_OrderedVertices[i] = i + i_LeftVertexCount; 00168 } 00169 00170 randomOrdering(m_vi_OrderedVertices); 00171 00172 return(_TRUE); 00173 } 00174 00175 //Public Function 2357 00176 int BipartiteGraphPartialOrdering::RowLargestFirstOrdering() 00177 { 00178 if(CheckVertexOrdering("ROW_LARGEST_FIRST")) 00179 { 00180 return(_TRUE); 00181 } 00182 00183 int i, j, k; 00184 int i_DegreeCount = 0; 00185 int i_VertexCount, i_Current; 00186 vector<int> vi_Visited; 00187 vector< vector<int> > vvi_GroupedVertexDegree; 00188 00189 i_VertexCount = (int)m_vi_LeftVertices.size () - 1; 00190 m_i_MaximumVertexDegree = 0; 00191 m_i_MinimumVertexDegree = i_VertexCount; 00192 vvi_GroupedVertexDegree.resize ( i_VertexCount ); 00193 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 00194 00195 // enter code here 00196 for ( i=0; i<i_VertexCount; ++i ) 00197 { 00198 // clear the visted nodes 00199 //vi_VistedNodes.clear (); 00200 // reset the degree count 00201 i_DegreeCount = 0; 00202 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 00203 for ( j=m_vi_LeftVertices [i]; j<m_vi_LeftVertices [i+1]; ++j ) 00204 { 00205 i_Current = m_vi_Edges [j]; 00206 00207 for ( k=m_vi_RightVertices [i_Current]; k<m_vi_RightVertices [i_Current+1]; ++k ) 00208 { 00209 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 00210 00211 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 00212 { 00213 ++i_DegreeCount; 00214 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 00215 vi_Visited [m_vi_Edges [k]] = i; 00216 } 00217 } 00218 } 00219 00220 vvi_GroupedVertexDegree [i_DegreeCount].push_back ( i ); 00221 00222 if ( m_i_MaximumVertexDegree < i_DegreeCount ) 00223 { 00224 m_i_MaximumVertexDegree = i_DegreeCount; 00225 } 00226 else if (m_i_MinimumVertexDegree > i_DegreeCount) 00227 { 00228 m_i_MinimumVertexDegree = i_DegreeCount; 00229 } 00230 } 00231 00232 if(i_VertexCount <2) m_i_MinimumVertexDegree = i_DegreeCount; 00233 00234 // clear the vertexorder 00235 m_vi_OrderedVertices.clear (); 00236 // take the bucket and place it in the vertexorder 00237 for ( i=m_i_MaximumVertexDegree; i>=m_i_MinimumVertexDegree; --i ) 00238 { 00239 // j = size of the bucket 00240 // get the size of the i-th bucket 00241 j = (int)vvi_GroupedVertexDegree [i].size (); 00242 // place it into vertex ordering 00243 for ( k=0; k<j; ++k ) 00244 { 00245 m_vi_OrderedVertices.push_back ( vvi_GroupedVertexDegree [i][k] ); 00246 } 00247 } 00248 00249 return(_TRUE); 00250 } 00251 00252 00253 00254 //Public Function 2358 00255 int BipartiteGraphPartialOrdering::ColumnLargestFirstOrdering() 00256 { 00257 if(CheckVertexOrdering("COLUMN_LARGEST_FIRST")) 00258 { 00259 return(_TRUE); 00260 } 00261 00262 int i, j, k; 00263 int i_DegreeCount = 0; 00264 int i_VertexCount, i_Current; 00265 vector<int> vi_Visited; 00266 vector< vector<int> > vvi_GroupedVertexDegree; 00267 00268 i_VertexCount = (int)m_vi_RightVertices.size () - 1; 00269 m_i_MaximumVertexDegree = 0; 00270 m_i_MinimumVertexDegree = i_VertexCount; 00271 vvi_GroupedVertexDegree.resize ( i_VertexCount ); 00272 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 00273 00274 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00275 00276 // enter code here 00277 for ( i=0; i<i_VertexCount; ++i ) 00278 { 00279 // clear the visted nodes 00280 //vi_VistedNodes.clear (); 00281 // reset the degree count 00282 i_DegreeCount = 0; 00283 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 00284 for ( j=m_vi_RightVertices [i]; j<m_vi_RightVertices [i+1]; ++j ) 00285 { 00286 i_Current = m_vi_Edges [j]; 00287 00288 for ( k=m_vi_LeftVertices [i_Current]; k<m_vi_LeftVertices [i_Current+1]; ++k ) 00289 { 00290 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 00291 00292 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 00293 { 00294 ++i_DegreeCount; 00295 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 00296 vi_Visited [m_vi_Edges [k]] = i; 00297 } 00298 } 00299 } 00300 vvi_GroupedVertexDegree [i_DegreeCount].push_back ( i ); 00301 00302 if ( m_i_MaximumVertexDegree < i_DegreeCount ) 00303 { 00304 m_i_MaximumVertexDegree = i_DegreeCount; 00305 } 00306 else if (m_i_MinimumVertexDegree > i_DegreeCount) 00307 { 00308 m_i_MinimumVertexDegree = i_DegreeCount; 00309 } 00310 } 00311 00312 if(i_VertexCount <2) m_i_MinimumVertexDegree = i_DegreeCount; 00313 00314 // clear the vertexorder 00315 m_vi_OrderedVertices.clear (); 00316 // take the bucket and place it in the vertexorder 00317 00318 for ( i=m_i_MaximumVertexDegree; i>=m_i_MinimumVertexDegree; --i ) 00319 { 00320 // j = size of the bucket 00321 // get the size of the i-th bucket 00322 j = (int)vvi_GroupedVertexDegree [i].size (); 00323 // place it into vertex ordering 00324 for ( k=0; k<j; ++k ) 00325 { 00326 m_vi_OrderedVertices.push_back ( vvi_GroupedVertexDegree [i][k] + i_LeftVertexCount); 00327 } 00328 } 00329 00330 return(_TRUE); 00331 } 00332 00333 00334 int BipartiteGraphPartialOrdering::RowSmallestLastOrdering() { 00335 return BipartiteGraphPartialOrdering::RowSmallestLastOrdering_serial(); 00336 /*#ifdef _OPENMP 00337 return BipartiteGraphPartialOrdering::RowSmallestLastOrdering_OMP(); 00338 #else 00339 return BipartiteGraphPartialOrdering::RowSmallestLastOrdering_serial(); 00340 #endif*/ 00341 } 00342 00343 //Line 1: procedure SMALLESTLASTORDERING-EASY(G = (V;E)) 00344 int BipartiteGraphPartialOrdering::RowSmallestLastOrdering_OMP() 00345 { 00346 // cout<<"IN ROW_SMALLEST_LAST_OMP()"<<endl<<flush; 00347 if(CheckVertexOrdering("ROW_SMALLEST_LAST_OMP")) 00348 { 00349 return(_TRUE); 00350 } 00351 00352 // PrintBipartiteGraph(); 00353 00354 int j, k, l, u; 00355 int i_LeftVertexCount = (signed) m_vi_LeftVertices.size() - 1; 00356 int i_RightVertexCount = (signed) m_vi_RightVertices.size() - 1; 00357 vector<int> vi_Visited; 00358 vi_Visited.clear(); 00359 vi_Visited.resize ( i_LeftVertexCount, _UNKNOWN ); 00360 m_vi_OrderedVertices.clear(); 00361 vector<int> d; // current distance-2 degree of each Left vertex 00362 d.resize(i_LeftVertexCount, _UNKNOWN); 00363 vector<int> VertexThreadGroup; 00364 VertexThreadGroup.resize(i_LeftVertexCount, _UNKNOWN); 00365 int i_MaxNumThreads; 00366 #ifdef _OPENMP 00367 i_MaxNumThreads = omp_get_max_threads(); 00368 #else 00369 i_MaxNumThreads = 1; 00370 #endif 00371 int i_MaxDegree = 0; 00372 int* i_MaxDegree_Private = new int[i_MaxNumThreads]; 00373 int* i_MinDegree_Private = new int[i_MaxNumThreads]; 00374 // ??? is this really neccessary ? => #pragma omp parallel for default(none) schedule(static) shared() 00375 for(int i=0; i < i_MaxNumThreads; i++) { 00376 i_MaxDegree_Private[i] = 0; 00377 i_MinDegree_Private[i] = i_LeftVertexCount; 00378 } 00379 int* delta = new int[i_MaxNumThreads]; 00380 00381 vector<int>** B; //private buckets. Each thread i will have their own buckets B[i][] 00382 B = new vector<int>*[i_MaxNumThreads]; 00383 #ifdef _OPENMP 00384 #pragma omp parallel for default(none) schedule(static) shared(B, i_MaxDegree, i_MaxNumThreads) 00385 #endif 00386 for(int i=0; i < i_MaxNumThreads; i++) { 00387 B[i] = new vector<int>[i_MaxDegree]; 00388 } 00389 00390 //DONE Line 2: for each vertex v in V in parallel do 00391 #ifdef _OPENMP 00392 #pragma omp parallel for default(none) schedule(static) shared(B, d, i_LeftVertexCount, i_MaxDegree_Private, i_MinDegree_Private) firstprivate(vi_Visited) 00393 #endif 00394 for(int v=0; v < i_LeftVertexCount; v++) { 00395 //DONE Line 3: d(v) <- d2(v,G) . Also find i_MaxDegree_Private 00396 d[v] = 0; 00397 for(int i=m_vi_LeftVertices[v]; i<m_vi_LeftVertices[v+1];i++) { 00398 int i_Current = m_vi_Edges[i]; 00399 for(int j=m_vi_RightVertices [i_Current]; j<m_vi_RightVertices [i_Current+1]; j++) { 00400 if ( m_vi_Edges [j] != v && vi_Visited [m_vi_Edges [j]] != v ) { 00401 d[v]++; 00402 vi_Visited [m_vi_Edges [j]] = v; 00403 } 00404 } 00405 } 00406 00407 int i_thread_num; 00408 #ifdef _OPENMP 00409 i_thread_num = omp_get_thread_num(); 00410 #else 00411 i_thread_num = 0; 00412 #endif 00413 if(i_MaxDegree_Private[i_thread_num]<d[v]) i_MaxDegree_Private[i_thread_num]=d[v]; 00414 if(i_MinDegree_Private[i_thread_num]>d[v]) { 00415 i_MinDegree_Private[i_thread_num]=d[v]; 00416 } 00417 } 00418 // find i_MaxDegree; populate delta 00419 for(int i=0; i < i_MaxNumThreads; i++) { 00420 if(i_MaxDegree<i_MaxDegree_Private[i] ) i_MaxDegree = i_MaxDegree_Private[i]; 00421 delta[i] = i_MinDegree_Private[i]; 00422 } 00423 00424 #ifdef _OPENMP 00425 #pragma omp parallel for default(none) schedule(static) shared(B, i_MaxDegree, i_MaxNumThreads) 00426 #endif 00427 for(int i=0; i < i_MaxNumThreads; i++) { 00428 int i_thread_num; 00429 #ifdef _OPENMP 00430 i_thread_num = omp_get_thread_num(); 00431 #else 00432 i_thread_num = 0; 00433 #endif 00434 B[i_thread_num] = new vector<int>[i_MaxDegree+1]; 00435 } 00436 00437 //DONE Line 2: for each vertex v in V in parallel do 00438 #ifdef _OPENMP 00439 #pragma omp parallel for default(none) schedule(static) shared(B, d, i_LeftVertexCount, VertexThreadGroup) 00440 #endif 00441 for(int v=0; v < i_LeftVertexCount; v++) { 00442 int i_thread_num; 00443 #ifdef _OPENMP 00444 i_thread_num = omp_get_thread_num(); 00445 #else 00446 i_thread_num = 0; 00447 #endif 00448 //DONE Line 4: add v to B_t(v) [d (v)] 00449 B[ i_thread_num ][ d[v] ].push_back(v); 00450 00451 //record that v is in B_t(v) 00452 VertexThreadGroup[v] = i_thread_num; 00453 00454 } 00455 00456 //DONE Line 5: i_NumOfVerticesToBeColored <- |V| 00457 int i_NumOfVerticesToBeColored = i_LeftVertexCount; 00458 00459 //Line 6: for k = 1 to p in parallel do 00460 #ifdef _OPENMP 00461 #pragma omp parallel for default(none) schedule(static) shared(i_MaxNumThreads, i_NumOfVerticesToBeColored, B, delta, VertexThreadGroup, d, cout, i_MaxDegree, i_MaxDegree_Private ) firstprivate(vi_Visited) 00462 #endif 00463 for(int k=0; k< i_MaxNumThreads; k++) { 00464 //reset vi_Visited 00465 for(int i=0; i< vi_Visited.size();i++) vi_Visited[i] = _UNKNOWN; 00466 00467 //Line 7: while i_NumOfVerticesToBeColored >= 0 do // !!! ??? why not i_NumOfVerticesToBeColored > 0 00468 while(i_NumOfVerticesToBeColored > 0) { 00469 int i_thread_num; 00470 #ifdef _OPENMP 00471 i_thread_num = omp_get_thread_num(); 00472 #else 00473 i_thread_num = 0; 00474 #endif 00475 //Line 8: Let delta be the smallest index j such that B_k [j] is non-empty 00476 //update delta 00477 // cout<<"delta[i_thread_num] 1="<< delta[i_thread_num] <<endl; 00478 // cout<<"B[i_thread_num]:"<<endl<<'\t'; 00479 // for(int i=0; i<=i_MaxDegree; i++) {cout<<B[i_thread_num][i].size()<<' ';} 00480 // cout<<'*'<<endl; 00481 if(delta[i_thread_num]!=0 && B[ i_thread_num ][ delta[i_thread_num] - 1 ].size() != 0) delta[i_thread_num] --; 00482 // cout<<"delta[i_thread_num] 2="<< delta[i_thread_num] <<endl; 00483 00484 //Line 9: Let v be a vertex drawn from B_k [delta] 00485 int v; 00486 00487 for(int i=delta[i_thread_num] ; i<i_MaxDegree_Private[i_thread_num]; i++) { 00488 if(B[ i_thread_num ][ i ].size()!=0) { 00489 v = B[ i_thread_num ][ delta[i_thread_num] ][ B[ i_thread_num ][ delta[i_thread_num] ].size() - 1 ]; 00490 d[v]= -1; // mark v as selected 00491 00492 //Line 10: remove v from B_k [delta] 00493 B[ i_thread_num ][ delta[i_thread_num] ].pop_back(); 00494 00495 break; 00496 } 00497 else delta[i_thread_num]++; 00498 } 00499 // cout<<"Select vertex v="<<v<<" ; d[v]="<< d[v]<<endl; 00500 // cout<<"delta[i_thread_num] 3="<< delta[i_thread_num] <<endl; 00501 00502 //Line 11: for each vertex w which is distance-2-neighbor of (v) do 00503 for(int l = m_vi_LeftVertices[v]; l < m_vi_LeftVertices[v+1]; l++) { 00504 int i_D1Neighbor = m_vi_Edges[l]; 00505 for(int m = m_vi_RightVertices[i_D1Neighbor]; m < m_vi_RightVertices[i_D1Neighbor+1]; m++ ) { 00506 int w = m_vi_Edges[m]; 00507 00508 //Line 12: if w in B_k then 00509 if( VertexThreadGroup[w] != i_thread_num || vi_Visited [w] == v || d[w] < 1 || w == v ) continue; 00510 00511 //Line 13: remove w from B_k [d (w)] 00512 // find location of w in B_k [d (w)] and pop it . !!! naive, improvable by keeping extra data. See if the extra data affacts concurrency 00513 int i_w_location = B[ i_thread_num ][ d[w] ].size() - 1; 00514 // cout<<"d[w]="<<d[w]<<endl; 00515 // cout<<"i_w_location before="<<i_w_location<<endl; 00516 // for(int ii = 0; ii <= i_w_location; ii++) {cout<<' '<< B[ i_thread_num ][ d[w] ][ii] ;} 00517 // cout<<"find w="<<w<<endl; 00518 while(i_w_location>=0 && B[ i_thread_num ][ d[w] ][i_w_location] != w) i_w_location--; 00519 // if(i_w_location<0) { 00520 // cout<<"*** i_w_location<0"<<endl<<flush; 00521 // } 00522 // cout<<"i_w_location after="<<i_w_location<<endl; 00523 if(i_w_location != (B[ i_thread_num ][ d[w] ].size() - 1) ) B[ i_thread_num ] [ d[w] ][i_w_location] = B[ i_thread_num ] [ d[w] ][ B[ i_thread_num ][ d[w] ].size() - 1 ]; 00524 B[ i_thread_num ] [ d[w] ].pop_back(); 00525 00526 //Line 14: d (w) <- d (w) - 1 00527 d[w]--; 00528 00529 //Line 15: add w to B_k [d (w)] 00530 B[ i_thread_num ] [ d[w] ].push_back(w); 00531 00532 } 00533 } 00534 00535 //DONE Line 16: W [i_NumOfVerticesToBeColored] <- v; i_NumOfVerticesToBeColored <- i_NumOfVerticesToBeColored - 1 . critical statements 00536 #ifdef _OPENMP 00537 #pragma omp critical 00538 #endif 00539 { 00541 m_vi_OrderedVertices.push_back(v); 00542 i_NumOfVerticesToBeColored--; 00543 // cout<<"i_NumOfVerticesToBeColored="<<i_NumOfVerticesToBeColored<<endl; 00544 } 00545 } 00546 } 00547 // cout<<"OUT ROW_SMALLEST_LAST_OMP()"<<endl<<flush; 00548 00549 return(_TRUE); 00550 } 00551 00552 //Public Function 2359 00553 int BipartiteGraphPartialOrdering::RowSmallestLastOrdering_serial() 00554 { 00555 if(CheckVertexOrdering("ROW_SMALLEST_LAST")) 00556 { 00557 return(_TRUE); 00558 } 00559 00560 int i, j, k, u, l; 00561 int i_Current; 00562 int i_SelectedVertex, i_SelectedVertexCount; 00563 int i_VertexCount; 00564 int i_VertexCountMinus1; // = i_VertexCount - 1, used when inserting selected vertices into m_vi_OrderedVertices 00565 int i_HighestInducedVertexDegree, i_InducedVertexDegree; 00566 vector<int> vi_InducedVertexDegree; 00567 vector<int> vi_Visited; 00568 vector< vector<int> > vvi_GroupedInducedVertexDegree; 00569 vector< int > vi_VertexLocation; 00570 00571 // initialize 00572 00573 i_SelectedVertex = _UNKNOWN; 00574 i_VertexCount = (int)m_vi_LeftVertices.size () - 1; 00575 i_VertexCountMinus1 = i_VertexCount - 1; 00576 i_HighestInducedVertexDegree = 0; 00577 vi_Visited.clear(); 00578 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 00579 m_vi_OrderedVertices.clear(); 00580 m_vi_OrderedVertices.resize(i_VertexCount, _UNKNOWN); 00581 00582 vi_InducedVertexDegree.clear(); 00583 vi_InducedVertexDegree.reserve((unsigned) i_VertexCount); 00584 vvi_GroupedInducedVertexDegree.clear(); 00585 vvi_GroupedInducedVertexDegree.resize((unsigned) i_VertexCount); 00586 vi_VertexLocation.clear(); 00587 vi_VertexLocation.reserve((unsigned) i_VertexCount); 00588 00589 for ( i=0; i<i_VertexCount; ++i ) 00590 { 00591 // clear the visted nodes 00592 //vi_VistedNodes.clear (); 00593 // reset the degree count 00594 i_InducedVertexDegree = 0; 00595 // let's loop from mvi_LeftVertices[i] to mvi_LeftVertices[i+1] for the i'th column 00596 for ( j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[i+1]; ++j ) 00597 { 00598 i_Current = m_vi_Edges [j]; 00599 00600 for (k=m_vi_RightVertices[i_Current]; k<m_vi_RightVertices[i_Current+1]; ++k) 00601 { 00602 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 00603 00604 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 00605 { 00606 ++i_InducedVertexDegree; 00607 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 00608 vi_Visited [m_vi_Edges [k]] = i; 00609 } 00610 } 00611 } 00612 00613 //vi_InducedVertexDegree[i] = vertex degree of vertex i 00614 vi_InducedVertexDegree.push_back ( i_InducedVertexDegree ); 00615 // vector vvi_GroupedInducedVertexDegree[i] = all the vertices with degree i 00616 // for every new vertex with degree i, it will be pushed to the back of vector vvi_GroupedInducedVertexDegree[i] 00617 vvi_GroupedInducedVertexDegree [i_InducedVertexDegree].push_back ( i ); 00618 //vi_VertexLocation[i] = location of vertex i in vvi_GroupedInducedVertexDegree[i_InducedVertexDegree] 00619 vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1); 00620 00621 //get max degree (i_HighestInducedVertexDegree) 00622 if ( i_HighestInducedVertexDegree < i_InducedVertexDegree ) 00623 { 00624 i_HighestInducedVertexDegree = i_InducedVertexDegree; 00625 } 00626 } 00627 00628 i_SelectedVertexCount = 0; 00629 // first clear the visited nodes 00630 vi_Visited.clear (); 00631 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 00632 // end clear nodes 00633 00634 int iMin = 1; 00635 00636 // just counting the number of vertices that we have worked with, 00637 // stop when i_SelectedVertexCount == i_VertexCount, i.e. we have looked through all the vertices 00638 while ( i_SelectedVertexCount < i_VertexCount ) 00639 { 00640 if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin - 1].size() != _FALSE) 00641 iMin--; 00642 00643 // selecte first item from the bucket 00644 for ( i=iMin; i<(i_HighestInducedVertexDegree+1); ++i ) 00645 { 00646 00647 if ( vvi_GroupedInducedVertexDegree[i].size () != 0 ) 00648 { 00649 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back (); 00650 //remove the i_SelectedVertex from vvi_GroupedInducedVertexDegree 00651 vvi_GroupedInducedVertexDegree[i].pop_back(); 00652 break; 00653 } 00654 else 00655 iMin++; 00656 } 00657 // end select first nonzero item from the bucket 00658 00659 // go to the neighbors of i_SelectedVertex and decrease the degrees by 1 00660 for ( i=m_vi_LeftVertices [i_SelectedVertex]; i<m_vi_LeftVertices [i_SelectedVertex+1]; ++i ) 00661 { 00662 // which Column element is Row (i_SelectedVertex) pointing to? 00663 i_Current = m_vi_Edges [i]; 00664 // go to each neighbor of Col (i_SelectedVertex), decrease their degree by 1 00665 // and then update their position in vvi_GroupedInducedVertexDegree and vi_VertexLocation 00666 for ( j=m_vi_RightVertices [i_Current]; j<m_vi_RightVertices [i_Current+1]; ++j ) 00667 { 00668 // note: m_vi_Edges [j] is the neighbor of Col (i_SelectedVertex) 00669 // make sure it's not pointing back to itself, or if we already visited this node 00670 if ( m_vi_Edges [j] == i_SelectedVertex || vi_Visited [m_vi_Edges [j]] == i_SelectedVertex ) 00671 { 00672 // it is pointed to itself or already visited 00673 continue; 00674 } 00675 00676 u = m_vi_Edges[j]; 00677 00678 // now check to make sure that the neighbor's degree isn't UNKNOWN 00679 if ( vi_InducedVertexDegree [u] == _UNKNOWN ) 00680 { 00681 // our neighbor doesn't have any neighbors, so skip it. 00682 continue; 00683 } 00684 00685 // update that we visited this 00686 vi_Visited [u] = i_SelectedVertex; 00687 // end up update that we visited this 00688 00689 // move the last element in this bucket to u's position to get rid of expensive erase operation 00690 if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1) 00691 { 00692 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back(); 00693 00694 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l; 00695 00696 00697 vi_VertexLocation[l] = vi_VertexLocation[u]; 00698 } 00699 // remove last element from this bucket 00700 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back(); 00701 00702 // reduce degree of u by 1 00703 vi_InducedVertexDegree[u]--; 00704 00705 // move u to appropriate bucket 00706 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u); 00707 00708 // update vi_VertexLocation[u] since it has now been changed 00709 vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1; 00710 00711 00712 00713 } 00714 } 00715 // end of go to the neighbors of i_SelectedVertex and decrease the degrees by 1 00716 00717 //Mark the i_SelectedVertex as _UNKNOWN, so that we don't look at it again 00718 vi_InducedVertexDegree [i_SelectedVertex] = _UNKNOWN; 00719 00720 m_vi_OrderedVertices[i_VertexCountMinus1 - i_SelectedVertexCount] = i_SelectedVertex; 00721 00722 // go to next vertex 00723 ++i_SelectedVertexCount; 00724 } 00725 00726 // clear the buffer 00727 vi_InducedVertexDegree.clear(); 00728 vi_VertexLocation.clear(); 00729 vvi_GroupedInducedVertexDegree.clear(); 00730 00731 return(_TRUE); 00732 } 00733 00734 int BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering() { 00735 return BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering_serial(); 00736 /*#ifdef _OPENMP 00737 return BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering_OMP(); 00738 #else 00739 return BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering_serial(); 00740 #endif*/ 00741 } 00742 00743 //Line 1: procedure SMALLESTLASTORDERING-EASY(G = (V;E)) 00744 int BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering_OMP() 00745 { 00746 // cout<<"IN COLUMN_SMALLEST_LAST_OMP()"<<endl<<flush; 00747 if(CheckVertexOrdering("COLUMN_SMALLEST_LAST_OMP")) 00748 { 00749 return(_TRUE); 00750 } 00751 00752 // PrintBipartiteGraph(); 00753 00754 int j, k, l, u; 00755 int i_LeftVertexCount = (signed) m_vi_LeftVertices.size() - 1; 00756 int i_RightVertexCount = (signed) m_vi_RightVertices.size() - 1; 00757 vector<int> vi_Visited; 00758 vi_Visited.clear(); 00759 vi_Visited.resize ( i_RightVertexCount, _UNKNOWN ); 00760 m_vi_OrderedVertices.clear(); 00761 vector<int> d; // current distance-2 degree of each right vertex 00762 d.resize(i_RightVertexCount, _UNKNOWN); 00763 vector<int> VertexThreadGroup; 00764 VertexThreadGroup.resize(i_RightVertexCount, _UNKNOWN); 00765 int i_MaxNumThreads; 00766 #ifdef _OPENMP 00767 i_MaxNumThreads = omp_get_max_threads(); 00768 #else 00769 i_MaxNumThreads = 1; 00770 #endif 00771 int i_MaxDegree = 0; 00772 int* i_MaxDegree_Private = new int[i_MaxNumThreads]; 00773 int* i_MinDegree_Private = new int[i_MaxNumThreads]; 00774 // ??? is this really neccessary ? => #pragma omp parallel for default(none) schedule(static) shared() 00775 for(int i=0; i < i_MaxNumThreads; i++) { 00776 i_MaxDegree_Private[i] = 0; 00777 i_MinDegree_Private[i] = i_RightVertexCount; 00778 } 00779 int* delta = new int[i_MaxNumThreads]; 00780 00781 vector<int>** B; //private buckets. Each thread i will have their own buckets B[i][] 00782 B = new vector<int>*[i_MaxNumThreads]; 00783 #ifdef _OPENMP 00784 #pragma omp parallel for default(none) schedule(static) shared(B, i_MaxDegree, i_MaxNumThreads) 00785 #endif 00786 for(int i=0; i < i_MaxNumThreads; i++) { 00787 B[i] = new vector<int>[i_MaxDegree]; 00788 } 00789 00790 //DONE Line 2: for each vertex v in V in parallel do 00791 #ifdef _OPENMP 00792 #pragma omp parallel for default(none) schedule(static) shared(B, d, i_RightVertexCount, i_MaxDegree_Private, i_MinDegree_Private) firstprivate(vi_Visited) 00793 #endif 00794 for(int v=0; v < i_RightVertexCount; v++) { 00795 //DONE Line 3: d(v) <- d2(v,G) . Also find i_MaxDegree_Private 00796 d[v] = 0; 00797 for(int i=m_vi_RightVertices[v]; i<m_vi_RightVertices[v+1];i++) { 00798 int i_Current = m_vi_Edges[i]; 00799 for(int j=m_vi_LeftVertices [i_Current]; j<m_vi_LeftVertices [i_Current+1]; j++) { 00800 if ( m_vi_Edges [j] != v && vi_Visited [m_vi_Edges [j]] != v ) { 00801 d[v]++; 00802 vi_Visited [m_vi_Edges [j]] = v; 00803 } 00804 } 00805 } 00806 00807 int i_thread_num; 00808 #ifdef _OPENMP 00809 i_thread_num = omp_get_thread_num(); 00810 #else 00811 i_thread_num = 0; 00812 #endif 00813 if(i_MaxDegree_Private[i_thread_num]<d[v]) i_MaxDegree_Private[i_thread_num]=d[v]; 00814 if(i_MinDegree_Private[i_thread_num]>d[v]) { 00815 i_MinDegree_Private[i_thread_num]=d[v]; 00816 } 00817 } 00818 // find i_MaxDegree; populate delta 00819 for(int i=0; i < i_MaxNumThreads; i++) { 00820 if(i_MaxDegree<i_MaxDegree_Private[i] ) i_MaxDegree = i_MaxDegree_Private[i]; 00821 delta[i] = i_MinDegree_Private[i]; 00822 } 00823 00824 #ifdef _OPENMP 00825 #pragma omp parallel for default(none) schedule(static) shared(B, i_MaxDegree, i_MaxNumThreads) 00826 #endif 00827 for(int i=0; i < i_MaxNumThreads; i++) { 00828 int i_thread_num; 00829 #ifdef _OPENMP 00830 i_thread_num = omp_get_thread_num(); 00831 #else 00832 i_thread_num = 0; 00833 #endif 00834 B[i_thread_num] = new vector<int>[i_MaxDegree+1]; 00835 } 00836 00837 //DONE Line 2: for each vertex v in V in parallel do 00838 #ifdef _OPENMP 00839 #pragma omp parallel for default(none) schedule(static) shared(B, d, i_RightVertexCount, VertexThreadGroup) 00840 #endif 00841 for(int v=0; v < i_RightVertexCount; v++) { 00842 int i_thread_num; 00843 #ifdef _OPENMP 00844 i_thread_num = omp_get_thread_num(); 00845 #else 00846 i_thread_num = 0; 00847 #endif 00848 //DONE Line 4: add v to B_t(v) [d (v)] 00849 B[ i_thread_num ][ d[v] ].push_back(v); 00850 00851 //record that v is in B_t(v) 00852 VertexThreadGroup[v] = i_thread_num; 00853 00854 } 00855 00856 //DONE Line 5: i_NumOfVerticesToBeColored <- |V| 00857 int i_NumOfVerticesToBeColored = i_RightVertexCount; 00858 00859 //Line 6: for k = 1 to p in parallel do 00860 #ifdef _OPENMP 00861 #pragma omp parallel for default(none) schedule(static) shared(i_LeftVertexCount, i_MaxNumThreads, i_NumOfVerticesToBeColored, B, delta, VertexThreadGroup, d, cout, i_MaxDegree, i_MaxDegree_Private ) firstprivate(vi_Visited) 00862 #endif 00863 for(int k=0; k< i_MaxNumThreads; k++) { 00864 //reset vi_Visited 00865 for(int i=0; i< vi_Visited.size();i++) vi_Visited[i] = _UNKNOWN; 00866 00867 //Line 7: while i_NumOfVerticesToBeColored >= 0 do // !!! ??? why not i_NumOfVerticesToBeColored > 0 00868 while(i_NumOfVerticesToBeColored > 0) { 00869 int i_thread_num; 00870 #ifdef _OPENMP 00871 i_thread_num = omp_get_thread_num(); 00872 #else 00873 i_thread_num = 0; 00874 #endif 00875 //Line 8: Let delta be the smallest index j such that B_k [j] is non-empty 00876 //update delta 00877 // cout<<"delta[i_thread_num] 1="<< delta[i_thread_num] <<endl; 00878 // cout<<"B[i_thread_num]:"<<endl<<'\t'; 00879 // for(int i=0; i<=i_MaxDegree; i++) {cout<<B[i_thread_num][i].size()<<' ';} 00880 // cout<<'*'<<endl; 00881 if(delta[i_thread_num]!=0 && B[ i_thread_num ][ delta[i_thread_num] - 1 ].size() != 0) delta[i_thread_num] --; 00882 // cout<<"delta[i_thread_num] 2="<< delta[i_thread_num] <<endl; 00883 00884 //Line 9: Let v be a vertex drawn from B_k [delta] 00885 int v; 00886 00887 for(int i=delta[i_thread_num] ; i<i_MaxDegree_Private[i_thread_num]; i++) { 00888 if(B[ i_thread_num ][ i ].size()!=0) { 00889 v = B[ i_thread_num ][ delta[i_thread_num] ][ B[ i_thread_num ][ delta[i_thread_num] ].size() - 1 ]; 00890 d[v]= -1; // mark v as selected 00891 00892 //Line 10: remove v from B_k [delta] 00893 B[ i_thread_num ][ delta[i_thread_num] ].pop_back(); 00894 00895 break; 00896 } 00897 else delta[i_thread_num]++; 00898 } 00899 // cout<<"Select vertex v="<<v<<" ; d[v]="<< d[v]<<endl; 00900 // cout<<"delta[i_thread_num] 3="<< delta[i_thread_num] <<endl; 00901 00902 //Line 11: for each vertex w which is distance-2-neighbor of (v) do 00903 for(int l = m_vi_RightVertices[v]; l < m_vi_RightVertices[v+1]; l++) { 00904 int i_D1Neighbor = m_vi_Edges[l]; 00905 for(int m = m_vi_LeftVertices[i_D1Neighbor]; m < m_vi_LeftVertices[i_D1Neighbor+1]; m++ ) { 00906 int w = m_vi_Edges[m]; 00907 00908 //Line 12: if w in B_k then 00909 if( VertexThreadGroup[w] != i_thread_num || vi_Visited [w] == v || d[w] < 1 || w == v ) continue; 00910 00911 //Line 13: remove w from B_k [d (w)] 00912 // find location of w in B_k [d (w)] and pop it . !!! naive, improvable by keeping extra data. See if the extra data affacts concurrency 00913 int i_w_location = B[ i_thread_num ][ d[w] ].size() - 1; 00914 // cout<<"d[w]="<<d[w]<<endl; 00915 // cout<<"i_w_location before="<<i_w_location<<endl; 00916 // for(int ii = 0; ii <= i_w_location; ii++) {cout<<' '<< B[ i_thread_num ][ d[w] ][ii] ;} 00917 // cout<<"find w="<<w<<endl; 00918 while(i_w_location>=0 && B[ i_thread_num ][ d[w] ][i_w_location] != w) i_w_location--; 00919 // if(i_w_location<0) { 00920 // cout<<"*** i_w_location<0"<<endl<<flush; 00921 // } 00922 // cout<<"i_w_location after="<<i_w_location<<endl; 00923 if(i_w_location != (B[ i_thread_num ][ d[w] ].size() - 1) ) B[ i_thread_num ] [ d[w] ][i_w_location] = B[ i_thread_num ] [ d[w] ][ B[ i_thread_num ][ d[w] ].size() - 1 ]; 00924 B[ i_thread_num ] [ d[w] ].pop_back(); 00925 00926 //Line 14: d (w) <- d (w) - 1 00927 d[w]--; 00928 00929 //Line 15: add w to B_k [d (w)] 00930 B[ i_thread_num ] [ d[w] ].push_back(w); 00931 00932 } 00933 } 00934 00935 //DONE Line 16: W [i_NumOfVerticesToBeColored] <- v; i_NumOfVerticesToBeColored <- i_NumOfVerticesToBeColored - 1 . critical statements 00936 #ifdef _OPENMP 00937 #pragma omp critical 00938 #endif 00939 { 00941 m_vi_OrderedVertices.push_back(v + i_LeftVertexCount); 00942 i_NumOfVerticesToBeColored--; 00943 // cout<<"i_NumOfVerticesToBeColored="<<i_NumOfVerticesToBeColored<<endl; 00944 } 00945 } 00946 } 00947 // cout<<"OUT COLUMN_SMALLEST_LAST_OMP()"<<endl<<flush; 00948 00949 return(_TRUE); 00950 } 00951 00952 00953 //Public Function 2360 00954 int BipartiteGraphPartialOrdering::ColumnSmallestLastOrdering_serial() 00955 { 00956 if(CheckVertexOrdering("COLUMN_SMALLEST_LAST")) 00957 { 00958 return(_TRUE); 00959 } 00960 00961 int i, j, k, u, l; 00962 int i_Current; 00963 int i_SelectedVertex, i_SelectedVertexCount; 00964 int i_VertexCount; 00965 int i_VertexCountMinus1; // = i_VertexCount - 1, used when inserting selected vertices into m_vi_OrderedVertices 00966 int i_HighestInducedVertexDegree, i_InducedVertexDegree; 00967 vector<int> vi_InducedVertexDegree; 00968 vector<int> vi_Visited; 00969 vector< vector<int> > vvi_GroupedInducedVertexDegree; 00970 vector< int > vi_VertexLocation; 00971 00972 // initialize 00973 i_SelectedVertex = _UNKNOWN; 00974 i_VertexCount = (int)m_vi_RightVertices.size () - 1; 00975 i_VertexCountMinus1 = i_VertexCount - 1; 00976 i_HighestInducedVertexDegree = 0; 00977 vi_Visited.clear(); 00978 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 00979 m_vi_OrderedVertices.clear(); 00980 m_vi_OrderedVertices.resize(i_VertexCount, _UNKNOWN); 00981 00982 vi_InducedVertexDegree.clear(); 00983 vi_InducedVertexDegree.reserve((unsigned) i_VertexCount); 00984 vvi_GroupedInducedVertexDegree.clear(); 00985 vvi_GroupedInducedVertexDegree.resize((unsigned) i_VertexCount); 00986 vi_VertexLocation.clear(); 00987 vi_VertexLocation.reserve((unsigned) i_VertexCount); 00988 00989 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 00990 00991 for ( i=0; i<i_VertexCount; ++i ) 00992 { 00993 // clear the visted nodes 00994 //vi_VistedNodes.clear (); 00995 // reset the degree count 00996 i_InducedVertexDegree = 0; 00997 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 00998 for ( j=m_vi_RightVertices [i]; j<m_vi_RightVertices [i+1]; ++j ) 00999 { 01000 i_Current = m_vi_Edges [j]; 01001 01002 for ( k=m_vi_LeftVertices [i_Current]; k<m_vi_LeftVertices [i_Current+1]; ++k ) 01003 { 01004 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 01005 01006 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 01007 { 01008 ++i_InducedVertexDegree; 01009 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 01010 vi_Visited [m_vi_Edges [k]] = i; 01011 } 01012 } 01013 } 01014 01015 //vi_InducedVertexDegree[i] = vertex degree of vertex i 01016 vi_InducedVertexDegree.push_back ( i_InducedVertexDegree ); 01017 // vector vvi_GroupedInducedVertexDegree[i] = all the vertices with degree i 01018 // for every new vertex with degree i, it will be pushed to the back of vector vvi_GroupedInducedVertexDegree[i] 01019 vvi_GroupedInducedVertexDegree [i_InducedVertexDegree].push_back ( i ); 01020 //vi_VertexLocation[i] = location of vertex i in vvi_GroupedInducedVertexDegree[i_InducedVertexDegree] 01021 vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1); 01022 01023 //get max degree (i_HighestInducedVertexDegree) 01024 if ( i_HighestInducedVertexDegree < i_InducedVertexDegree ) 01025 { 01026 i_HighestInducedVertexDegree = i_InducedVertexDegree; 01027 } 01028 } 01029 01030 i_SelectedVertexCount = 0; 01031 // first clear the visited nodes 01032 vi_Visited.clear (); 01033 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01034 // end clear nodes 01035 01036 int iMin = 1; 01037 01038 // just counting the number of vertices that we have worked with, 01039 // stop when i_SelectedVertexCount == i_VertexCount, i.e. we have looked through all the vertices 01040 while ( i_SelectedVertexCount < i_VertexCount ) 01041 { 01042 if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin - 1].size() != _FALSE) 01043 iMin--; 01044 01045 // selecte first item from the bucket 01046 for ( i=iMin; i<(i_HighestInducedVertexDegree+1); ++i ) 01047 { 01048 01049 if ( vvi_GroupedInducedVertexDegree[i].size () != 0 ) 01050 { 01051 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back (); 01052 //remove the i_SelectedVertex from vvi_GroupedInducedVertexDegree 01053 vvi_GroupedInducedVertexDegree[i].pop_back(); 01054 break; 01055 } 01056 else 01057 iMin++; 01058 } 01059 // end select first nonzero item from the bucket 01060 01061 // go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01062 for ( i=m_vi_RightVertices [i_SelectedVertex]; i<m_vi_RightVertices [i_SelectedVertex+1]; ++i ) 01063 { 01064 // which Row element is Col (i_SelectedVertex) pointing to? 01065 i_Current = m_vi_Edges [i]; 01066 // go to each neighbor of Col (i_SelectedVertex), decrease their degree by 1 01067 // and then update their position in vvi_GroupedInducedVertexDegree and vi_VertexLocation 01068 for ( j=m_vi_LeftVertices [i_Current]; j<m_vi_LeftVertices [i_Current+1]; ++j ) 01069 { 01070 // note: m_vi_Edges [j] is the neighbor of Col (i_SelectedVertex) 01071 // make sure it's not pointing back to itself, or if we already visited this node 01072 if ( m_vi_Edges [j] == i_SelectedVertex || vi_Visited [m_vi_Edges [j]] == i_SelectedVertex ) 01073 { 01074 // it is pointed to itself or already visited 01075 continue; 01076 } 01077 01078 u = m_vi_Edges[j]; 01079 01080 // now check to make sure that the neighbor's degree isn't UNKNOWN 01081 if ( vi_InducedVertexDegree [u] == _UNKNOWN ) 01082 { 01083 // our neighbor doesn't have any neighbors, so skip it. 01084 continue; 01085 } 01086 01087 // update that we visited this 01088 vi_Visited [u] = i_SelectedVertex; 01089 // end up update that we visited this 01090 01091 // move the last element in this bucket to u's position to get rid of expensive erase operation 01092 if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1) 01093 { 01094 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back(); 01095 01096 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l; 01097 01098 01099 vi_VertexLocation[l] = vi_VertexLocation[u]; 01100 } 01101 // remove last element from this bucket 01102 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back(); 01103 01104 // reduce degree of u by 1 01105 vi_InducedVertexDegree[u]--; 01106 01107 // move u to appropriate bucket 01108 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u); 01109 01110 // update vi_VertexLocation[u] since it has now been changed 01111 vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1; 01112 } 01113 } 01114 // end of go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01115 01116 //Mark the i_SelectedVertex as _UNKNOWN, so that we don't look at it again 01117 vi_InducedVertexDegree [i_SelectedVertex] = _UNKNOWN; 01118 01119 m_vi_OrderedVertices[i_VertexCountMinus1 - i_SelectedVertexCount] = i_SelectedVertex + i_LeftVertexCount; 01120 01121 // go to next vertex 01122 ++i_SelectedVertexCount; 01123 } 01124 01125 // clear the buffer 01126 vi_InducedVertexDegree.clear(); 01127 vi_VertexLocation.clear(); 01128 vvi_GroupedInducedVertexDegree.clear(); 01129 01130 return(_TRUE); 01131 } 01132 01133 int BipartiteGraphPartialOrdering::ColumnDynamicLargestFirstOrdering() 01134 { 01135 if(CheckVertexOrdering("COLUMN_DYNAMIC_LARGEST_FIRST")) 01136 { 01137 return(_TRUE); 01138 } 01139 01140 int i, j, k, u, l; 01141 int i_Current; 01142 int i_SelectedVertex, i_SelectedVertexCount; 01143 int i_VertexCount; 01144 01145 int i_HighestInducedVertexDegree, i_InducedVertexDegree; 01146 01147 vector < int > vi_InducedVertexDegree; 01148 vector < int > vi_Visited; 01149 vector < vector < int > > vvi_GroupedInducedVertexDegree; 01150 vector < int > vi_VertexLocation; 01151 01152 // initialize 01153 i_SelectedVertex = _UNKNOWN; 01154 i_VertexCount = (int)m_vi_RightVertices.size () - 1; 01155 01156 i_HighestInducedVertexDegree = 0; 01157 01158 vi_Visited.clear(); 01159 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01160 01161 m_vi_OrderedVertices.clear(); 01162 m_vi_OrderedVertices.reserve(i_VertexCount); 01163 01164 vi_InducedVertexDegree.clear(); 01165 vi_InducedVertexDegree.reserve((unsigned) i_VertexCount); 01166 01167 vvi_GroupedInducedVertexDegree.clear(); 01168 vvi_GroupedInducedVertexDegree.resize((unsigned) i_VertexCount); 01169 01170 vi_VertexLocation.clear(); 01171 vi_VertexLocation.reserve((unsigned) i_VertexCount); 01172 01173 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 01174 01175 for ( i=0; i<i_VertexCount; ++i) 01176 { 01177 // clear the visted nodes 01178 //vi_VistedNodes.clear (); 01179 // reset the degree count 01180 i_InducedVertexDegree = 0; 01181 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 01182 for ( j=m_vi_RightVertices [i]; j<m_vi_RightVertices [i+1]; ++j ) 01183 { 01184 i_Current = m_vi_Edges [j]; 01185 01186 for ( k=m_vi_LeftVertices [i_Current]; k<m_vi_LeftVertices [i_Current+1]; ++k ) 01187 { 01188 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 01189 01190 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 01191 { 01192 ++i_InducedVertexDegree; 01193 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 01194 vi_Visited [m_vi_Edges [k]] = i; 01195 } 01196 } 01197 } 01198 01199 //vi_InducedVertexDegree[i] = vertex degree of vertex i 01200 vi_InducedVertexDegree.push_back ( i_InducedVertexDegree ); 01201 // vector vvi_GroupedInducedVertexDegree[i] = all the vertices with degree i 01202 // for every new vertex with degree i, it will be pushed to the back of vector vvi_GroupedInducedVertexDegree[i] 01203 vvi_GroupedInducedVertexDegree [i_InducedVertexDegree].push_back ( i ); 01204 //vi_VertexLocation[i] = location of vertex i in vvi_GroupedInducedVertexDegree[i_InducedVertexDegree] 01205 vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1); 01206 01207 //get max degree (i_HighestInducedVertexDegree) 01208 if ( i_HighestInducedVertexDegree < i_InducedVertexDegree ) 01209 { 01210 i_HighestInducedVertexDegree = i_InducedVertexDegree; 01211 } 01212 } 01213 01214 i_SelectedVertexCount = 0; 01215 // first clear the visited nodes 01216 vi_Visited.clear (); 01217 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01218 // end clear nodes 01219 01220 //int iMin = 1; 01221 01222 // just counting the number of vertices that we have worked with, 01223 // stop when i_SelectedVertexCount == i_VertexCount, i.e. we have looked through all the vertices 01224 while ( i_SelectedVertexCount < i_VertexCount ) 01225 { 01226 //if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin - 1].size() != _FALSE) 01227 // iMin--; 01228 01229 // selecte first item from the bucket 01230 for ( i= i_HighestInducedVertexDegree; i>= 0; i-- ) 01231 { 01232 if ( vvi_GroupedInducedVertexDegree[i].size () != 0 ) 01233 { 01234 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back (); 01235 //remove the i_SelectedVertex from vvi_GroupedInducedVertexDegree 01236 vvi_GroupedInducedVertexDegree[i].pop_back(); 01237 break; 01238 } 01239 else 01240 i_HighestInducedVertexDegree--; 01241 } 01242 // end select first nonzero item from the bucket 01243 01244 // go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01245 for ( i=m_vi_RightVertices [i_SelectedVertex]; i<m_vi_RightVertices [i_SelectedVertex+1]; ++i ) 01246 { 01247 // which Row element is Col (i_SelectedVertex) pointing to? 01248 i_Current = m_vi_Edges [i]; 01249 // go to each neighbor of Col (i_SelectedVertex), decrease their degree by 1 01250 // and then update their position in vvi_GroupedInducedVertexDegree and vi_VertexLocation 01251 for ( j=m_vi_LeftVertices [i_Current]; j<m_vi_LeftVertices [i_Current+1]; ++j ) 01252 { 01253 // note: m_vi_Edges [j] is the neighbor of Col (i_SelectedVertex) 01254 // make sure it's not pointing back to itself, or if we already visited this node 01255 if ( m_vi_Edges [j] == i_SelectedVertex || vi_Visited [m_vi_Edges [j]] == i_SelectedVertex ) 01256 { 01257 // it is pointed to itself or already visited 01258 continue; 01259 } 01260 01261 u = m_vi_Edges[j]; 01262 01263 // now check to make sure that the neighbor's degree isn't UNKNOWN 01264 if ( vi_InducedVertexDegree [u] == _UNKNOWN ) 01265 { 01266 // our neighbor doesn't have any neighbors, so skip it. 01267 continue; 01268 } 01269 01270 // update that we visited this 01271 vi_Visited [u] = i_SelectedVertex; 01272 // end up update that we visited this 01273 01274 // move the last element in this bucket to u's position to get rid of expensive erase operation 01275 if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1) 01276 { 01277 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back(); 01278 01279 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l; 01280 01281 01282 vi_VertexLocation[l] = vi_VertexLocation[u]; 01283 } 01284 // remove last element from this bucket 01285 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back(); 01286 01287 // reduce degree of u by 1 01288 vi_InducedVertexDegree[u]--; 01289 01290 // move u to appropriate bucket 01291 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u); 01292 01293 // update vi_VertexLocation[u] since it has now been changed 01294 vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1; 01295 } 01296 } 01297 // end of go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01298 01299 //Mark the i_SelectedVertex as _UNKNOWN, so that we don't look at it again 01300 vi_InducedVertexDegree [i_SelectedVertex] = _UNKNOWN; 01301 01302 m_vi_OrderedVertices.push_back(i_SelectedVertex + i_LeftVertexCount); 01303 01304 // go to next vertex 01305 ++i_SelectedVertexCount; 01306 } 01307 01308 // clear the buffer 01309 vi_InducedVertexDegree.clear(); 01310 vi_VertexLocation.clear(); 01311 vvi_GroupedInducedVertexDegree.clear(); 01312 vi_Visited.clear (); 01313 return(_TRUE); 01314 } 01315 01316 int BipartiteGraphPartialOrdering::RowDynamicLargestFirstOrdering() 01317 { 01318 if(CheckVertexOrdering("ROW_DYNAMIC_LARGEST_FIRST")) 01319 { 01320 return(_TRUE); 01321 } 01322 01323 int i, j, k, u, l; 01324 int i_Current; 01325 int i_SelectedVertex, i_SelectedVertexCount; 01326 int i_VertexCount; 01327 01328 int i_HighestInducedVertexDegree, i_InducedVertexDegree; 01329 vector<int> vi_InducedVertexDegree; 01330 vector<int> vi_Visited; 01331 vector< vector<int> > vvi_GroupedInducedVertexDegree; 01332 vector< int > vi_VertexLocation; 01333 01334 // initialize 01335 01336 i_SelectedVertex = _UNKNOWN; 01337 i_VertexCount = (int)m_vi_LeftVertices.size () - 1; 01338 i_HighestInducedVertexDegree = 0; 01339 01340 vi_Visited.clear(); 01341 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01342 01343 m_vi_OrderedVertices.clear(); 01344 m_vi_OrderedVertices.reserve(i_VertexCount); 01345 01346 vi_InducedVertexDegree.clear(); 01347 vi_InducedVertexDegree.reserve((unsigned) i_VertexCount); 01348 vvi_GroupedInducedVertexDegree.clear(); 01349 vvi_GroupedInducedVertexDegree.resize((unsigned) i_VertexCount); 01350 vi_VertexLocation.clear(); 01351 vi_VertexLocation.reserve((unsigned) i_VertexCount); 01352 01353 for ( i=0; i<i_VertexCount; ++i ) 01354 { 01355 // clear the visted nodes 01356 //vi_VistedNodes.clear (); 01357 // reset the degree count 01358 i_InducedVertexDegree = 0; 01359 // let's loop from mvi_LeftVertices[i] to mvi_LeftVertices[i+1] for the i'th column 01360 for ( j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[i+1]; ++j ) 01361 { 01362 i_Current = m_vi_Edges [j]; 01363 01364 for (k=m_vi_RightVertices[i_Current]; k<m_vi_RightVertices[i_Current+1]; ++k) 01365 { 01366 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 01367 01368 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 01369 { 01370 ++i_InducedVertexDegree; 01371 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 01372 vi_Visited [m_vi_Edges [k]] = i; 01373 } 01374 } 01375 } 01376 01377 //vi_InducedVertexDegree[i] = vertex degree of vertex i 01378 vi_InducedVertexDegree.push_back ( i_InducedVertexDegree ); 01379 // vector vvi_GroupedInducedVertexDegree[i] = all the vertices with degree i 01380 // for every new vertex with degree i, it will be pushed to the back of vector vvi_GroupedInducedVertexDegree[i] 01381 vvi_GroupedInducedVertexDegree [i_InducedVertexDegree].push_back ( i ); 01382 //vi_VertexLocation[i] = location of vertex i in vvi_GroupedInducedVertexDegree[i_InducedVertexDegree] 01383 vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1); 01384 01385 //get max degree (i_HighestInducedVertexDegree) 01386 if ( i_HighestInducedVertexDegree < i_InducedVertexDegree ) 01387 { 01388 i_HighestInducedVertexDegree = i_InducedVertexDegree; 01389 } 01390 } 01391 01392 i_SelectedVertexCount = 0; 01393 // first clear the visited nodes 01394 vi_Visited.clear (); 01395 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01396 // end clear nodes 01397 01398 //int iMin = 1; 01399 01400 // just counting the number of vertices that we have worked with, 01401 // stop when i_SelectedVertexCount == i_VertexCount, i.e. we have looked through all the vertices 01402 while ( i_SelectedVertexCount < i_VertexCount ) 01403 { 01404 //if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin - 1].size() != _FALSE) 01405 // iMin--; 01406 01407 // selecte first item from the bucket 01408 for ( i= i_HighestInducedVertexDegree; i>= 0; i-- ) 01409 { 01410 01411 if ( vvi_GroupedInducedVertexDegree[i].size () != 0 ) 01412 { 01413 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back (); 01414 //remove the i_SelectedVertex from vvi_GroupedInducedVertexDegree 01415 vvi_GroupedInducedVertexDegree[i].pop_back(); 01416 break; 01417 } 01418 else 01419 i_HighestInducedVertexDegree--; 01420 } 01421 // end select first nonzero item from the bucket 01422 01423 // go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01424 for ( i=m_vi_LeftVertices [i_SelectedVertex]; i<m_vi_LeftVertices [i_SelectedVertex+1]; ++i ) 01425 { 01426 // which Column element is Row (i_SelectedVertex) pointing to? 01427 i_Current = m_vi_Edges [i]; 01428 // go to each neighbor of Col (i_SelectedVertex), decrease their degree by 1 01429 // and then update their position in vvi_GroupedInducedVertexDegree and vi_VertexLocation 01430 for ( j=m_vi_RightVertices [i_Current]; j<m_vi_RightVertices [i_Current+1]; ++j ) 01431 { 01432 // note: m_vi_Edges [j] is the neighbor of Col (i_SelectedVertex) 01433 // make sure it's not pointing back to itself, or if we already visited this node 01434 if ( m_vi_Edges [j] == i_SelectedVertex || vi_Visited [m_vi_Edges [j]] == i_SelectedVertex ) 01435 { 01436 // it is pointed to itself or already visited 01437 continue; 01438 } 01439 01440 u = m_vi_Edges[j]; 01441 01442 // now check to make sure that the neighbor's degree isn't UNKNOWN 01443 if ( vi_InducedVertexDegree [u] == _UNKNOWN ) 01444 { 01445 // our neighbor doesn't have any neighbors, so skip it. 01446 continue; 01447 } 01448 01449 // update that we visited this 01450 vi_Visited [u] = i_SelectedVertex; 01451 // end up update that we visited this 01452 01453 // move the last element in this bucket to u's position to get rid of expensive erase operation 01454 if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1) 01455 { 01456 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back(); 01457 01458 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l; 01459 01460 01461 vi_VertexLocation[l] = vi_VertexLocation[u]; 01462 } 01463 // remove last element from this bucket 01464 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back(); 01465 01466 // reduce degree of u by 1 01467 vi_InducedVertexDegree[u]--; 01468 01469 // move u to appropriate bucket 01470 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u); 01471 01472 // update vi_VertexLocation[u] since it has now been changed 01473 vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1; 01474 01475 } 01476 } 01477 // end of go to the neighbors of i_SelectedVertex and decrease the degrees by 1 01478 01479 //Mark the i_SelectedVertex as _UNKNOWN, so that we don't look at it again 01480 vi_InducedVertexDegree [i_SelectedVertex] = _UNKNOWN; 01481 01482 m_vi_OrderedVertices.push_back(i_SelectedVertex); 01483 01484 // go to next vertex 01485 ++i_SelectedVertexCount; 01486 } 01487 01488 // clear the buffer 01489 vi_InducedVertexDegree.clear(); 01490 vi_VertexLocation.clear(); 01491 vvi_GroupedInducedVertexDegree.clear(); 01492 vi_Visited.clear(); 01493 01494 return(_TRUE); 01495 } 01496 01497 //Public Function 2361 01498 int BipartiteGraphPartialOrdering::RowIncidenceDegreeOrdering() 01499 { 01500 if(CheckVertexOrdering("ROW_INCIDENCE_DEGREE")) 01501 { 01502 return(_TRUE); 01503 } 01504 01505 int i, j, k, l, u; 01506 int i_Current; 01507 int i_HighestDegreeVertex, m_i_MaximumVertexDegree; 01508 int i_VertexCount, i_VertexDegree, i_IncidenceVertexDegree; 01509 int i_SelectedVertex, i_SelectedVertexCount; 01510 vector<int> vi_IncidenceVertexDegree; 01511 vector<int> vi_Visited; 01512 vector< vector <int> > vvi_GroupedIncidenceVertexDegree; 01513 vector< int > vi_VertexLocation; 01514 01515 // initialize 01516 i_SelectedVertex = _UNKNOWN; 01517 i_VertexCount = (int)m_vi_LeftVertices.size () - 1; 01518 vvi_GroupedIncidenceVertexDegree.clear(); 01519 vvi_GroupedIncidenceVertexDegree.resize ( i_VertexCount ); 01520 i_HighestDegreeVertex = _UNKNOWN; 01521 m_i_MaximumVertexDegree = _UNKNOWN; 01522 i_IncidenceVertexDegree = 0; 01523 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01524 m_vi_OrderedVertices.clear(); 01525 m_vi_OrderedVertices.reserve(i_VertexCount); 01526 vi_VertexLocation.clear(); 01527 vi_VertexLocation.reserve(i_VertexCount); 01528 vi_IncidenceVertexDegree.clear(); 01529 vi_IncidenceVertexDegree.reserve(i_VertexCount); 01530 01531 for ( i=0; i<i_VertexCount; ++i ) 01532 { 01533 // clear the visted nodes 01534 //vi_VistedNodes.clear (); 01535 // reset the degree count 01536 i_VertexDegree = 0; 01537 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 01538 for ( j=m_vi_LeftVertices [i]; j<m_vi_LeftVertices [i+1]; ++j ) 01539 { 01540 i_Current = m_vi_Edges [j]; 01541 01542 for ( k=m_vi_RightVertices [i_Current]; k<m_vi_RightVertices [i_Current+1]; ++k ) 01543 { 01544 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 01545 01546 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 01547 { 01548 ++i_VertexDegree; 01549 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 01550 vi_Visited [m_vi_Edges [k]] = i; 01551 } 01552 } 01553 } 01554 vi_IncidenceVertexDegree.push_back ( i_IncidenceVertexDegree ); 01555 vvi_GroupedIncidenceVertexDegree [i_IncidenceVertexDegree].push_back ( i ); 01556 vi_VertexLocation.push_back ( vvi_GroupedIncidenceVertexDegree [i_IncidenceVertexDegree].size () -1 ); 01557 01558 if ( m_i_MaximumVertexDegree < i_VertexDegree ) 01559 { 01560 m_i_MaximumVertexDegree = i_VertexDegree; 01561 i_HighestDegreeVertex = i; 01562 } 01563 } 01564 01565 // initialize more things for the bucket "moving" 01566 i_SelectedVertexCount = 0; 01567 vi_Visited.clear (); 01568 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01569 01570 int iMax = 0; 01571 01572 while ( i_SelectedVertexCount < i_VertexCount ) 01573 { 01574 if(iMax != m_i_MaximumVertexDegree && vvi_GroupedIncidenceVertexDegree[iMax + 1].size() != _FALSE) 01575 iMax++; 01576 01577 for ( i=iMax; i>=0; i-- ) 01578 { 01579 if ( (int)vvi_GroupedIncidenceVertexDegree [i].size () != 0 ) 01580 { 01581 i_SelectedVertex = vvi_GroupedIncidenceVertexDegree [i].back (); 01582 // now we remove i_SelectedVertex from the bucket 01583 vvi_GroupedIncidenceVertexDegree [i].pop_back(); 01584 break; 01585 } 01586 } 01587 // end select the vertex with the highest Incidence degree 01588 01589 // tell the neighbors of i_SelectedVertex that we are removing it 01590 for ( i=m_vi_LeftVertices [i_SelectedVertex]; i<m_vi_LeftVertices [i_SelectedVertex+1]; ++i ) 01591 { 01592 i_Current = m_vi_Edges [i]; 01593 01594 for ( j=m_vi_RightVertices [i_Current]; j<m_vi_RightVertices [i_Current+1]; ++j ) 01595 { 01596 u = m_vi_Edges[j]; 01597 01598 // now check if the degree of i_SelectedVertex's neighbor isn't unknow to us 01599 if ( vi_IncidenceVertexDegree [u] == _UNKNOWN ) 01600 { 01601 // i_SelectedVertex's neighbor's degree is unknown. skip! 01602 continue; 01603 } 01604 01605 // note: u = neighbor of i_SelectedVertex 01606 // first check if we are pointing to itself or if we already visited this neighbor 01607 if ( u == i_SelectedVertex || vi_Visited [u] == i_SelectedVertex ) 01608 { 01609 // we already visited or pointing to itself. skip. 01610 continue; 01611 } 01612 01613 // move the last element in this bucket to u's position to get rid of expensive erase operation 01614 if(vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].size() > 1) 01615 { 01616 l = vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].back(); 01617 01618 vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]][vi_VertexLocation[u]] = l; 01619 01620 vi_VertexLocation[l] = vi_VertexLocation[u]; 01621 } 01622 01623 // remove the last element from vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]] 01624 vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].pop_back(); 01625 01626 // now update that we are visiting u 01627 vi_Visited [u] = i_SelectedVertex; 01628 // end now update that we are visiting u 01629 01630 // now we tell that neighbor to increase its degree by 1 01631 ++vi_IncidenceVertexDegree [u]; 01632 // place the neighbor in 1 degree up in the bucket 01633 vvi_GroupedIncidenceVertexDegree [vi_IncidenceVertexDegree [u]].push_back ( u ); 01634 // update the location of the now moved neighbor 01635 vi_VertexLocation [u] = vvi_GroupedIncidenceVertexDegree [vi_IncidenceVertexDegree [u]].size () -1; 01636 } 01637 } 01638 // end tell the neighbors of i_SelectedVertex that we are removing it 01639 01640 // now we say that i_SelectedVertex's degree is unknown to us 01641 vi_IncidenceVertexDegree [i_SelectedVertex] = _UNKNOWN; 01642 // now we push i_SelectedVertex to the back or VertexOrder 01643 m_vi_OrderedVertices.push_back ( i_SelectedVertex ); 01644 // loop it again 01645 ++i_SelectedVertexCount; 01646 } 01647 01648 // clear the buffer 01649 vi_IncidenceVertexDegree.clear(); 01650 vi_VertexLocation.clear(); 01651 vvi_GroupedIncidenceVertexDegree.clear(); 01652 01653 return(_TRUE); 01654 01655 } 01656 01657 01658 01659 //Public Function 2362 01660 int BipartiteGraphPartialOrdering::ColumnIncidenceDegreeOrdering() 01661 { 01662 if(CheckVertexOrdering("COLUMN_INCIDENCE_DEGREE")) 01663 { 01664 return(_TRUE); 01665 } 01666 01667 int i, j, k, l, u; 01668 int i_Current; 01669 int i_HighestDegreeVertex, m_i_MaximumVertexDegree; 01670 int i_VertexCount, i_VertexDegree, i_IncidenceVertexDegree; 01671 int i_SelectedVertex, i_SelectedVertexCount; 01672 vector<int> vi_IncidenceVertexDegree; 01673 vector<int> vi_Visited; 01674 vector< vector<int> > vvi_GroupedIncidenceVertexDegree; 01675 vector< int > vi_VertexLocation; 01676 01677 // initialize 01678 i_SelectedVertex = _UNKNOWN; 01679 i_VertexCount = (int)m_vi_RightVertices.size () - 1; 01680 vvi_GroupedIncidenceVertexDegree.resize ( i_VertexCount ); 01681 i_HighestDegreeVertex = _UNKNOWN; 01682 m_i_MaximumVertexDegree = _UNKNOWN; 01683 i_IncidenceVertexDegree = 0; 01684 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01685 m_vi_OrderedVertices.clear(); 01686 01687 int i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size()); 01688 01689 // enter code here 01690 for ( i=0; i<i_VertexCount; ++i ) 01691 { 01692 // clear the visted nodes 01693 //vi_VistedNodes.clear (); 01694 // reset the degree count 01695 i_VertexDegree = 0; 01696 // let's loop from mvi_RightVertices[i] to mvi_RightVertices[i+1] for the i'th column 01697 for ( j=m_vi_RightVertices [i]; j<m_vi_RightVertices [i+1]; ++j ) 01698 { 01699 i_Current = m_vi_Edges [j]; 01700 01701 for ( k=m_vi_LeftVertices [i_Current]; k<m_vi_LeftVertices [i_Current+1]; ++k ) 01702 { 01703 // b_visited = visitedAlready ( vi_VistedNodes, m_vi_Edges [k] ); 01704 01705 if ( m_vi_Edges [k] != i && vi_Visited [m_vi_Edges [k]] != i ) 01706 { 01707 ++i_VertexDegree; 01708 // vi_VistedNodes.push_back ( m_vi_Edges [k] ); 01709 vi_Visited [m_vi_Edges [k]] = i; 01710 } 01711 } 01712 } 01713 vi_IncidenceVertexDegree.push_back ( i_IncidenceVertexDegree ); 01714 vvi_GroupedIncidenceVertexDegree [i_IncidenceVertexDegree].push_back ( i ); 01715 vi_VertexLocation.push_back ( vvi_GroupedIncidenceVertexDegree [i_IncidenceVertexDegree].size () - 1); 01716 01717 if ( m_i_MaximumVertexDegree < i_VertexDegree ) 01718 { 01719 m_i_MaximumVertexDegree = i_VertexDegree; 01720 i_HighestDegreeVertex = i; 01721 } 01722 } 01723 01724 // initialize more things for the bucket "moving" 01725 i_SelectedVertexCount = 0; 01726 vi_Visited.clear (); 01727 vi_Visited.resize ( i_VertexCount, _UNKNOWN ); 01728 01729 int iMax = 0; 01730 01731 while ( i_SelectedVertexCount < i_VertexCount ) 01732 { 01733 if(iMax != m_i_MaximumVertexDegree && vvi_GroupedIncidenceVertexDegree[iMax + 1].size() != _FALSE) 01734 iMax++; 01735 01736 // select the vertex with the highest Incidence degree 01737 for ( i=m_i_MaximumVertexDegree; i>=0; i-- ) 01738 { 01739 if ( (int)vvi_GroupedIncidenceVertexDegree [i].size () != 0 ) 01740 { 01741 i_SelectedVertex = vvi_GroupedIncidenceVertexDegree [i].back (); 01742 // now we remove i_SelectedVertex from the bucket 01743 vvi_GroupedIncidenceVertexDegree [i].pop_back(); 01744 break; 01745 } 01746 } 01747 // end select the vertex with the highest Incidence degree 01748 01749 // tell the neighbors of i_SelectedVertex that we are removing it 01750 for ( i=m_vi_RightVertices [i_SelectedVertex]; i<m_vi_RightVertices [i_SelectedVertex+1]; ++i ) 01751 { 01752 i_Current = m_vi_Edges [i]; 01753 01754 for ( j=m_vi_LeftVertices [i_Current]; j<m_vi_LeftVertices [i_Current+1]; ++j ) 01755 { 01756 u = m_vi_Edges[j]; 01757 01758 // now check if the degree of i_SelectedVertex's neighbor isn't unknow to us 01759 if ( vi_IncidenceVertexDegree [u] == _UNKNOWN ) 01760 { 01761 // i_SelectedVertex's neighbor's degree is unknown. skip! 01762 continue; 01763 } 01764 01765 // note: u = neighbor of i_SelectedVertex 01766 // first check if we are pointing to itself or if we already visited a neighbor 01767 if ( u == i_SelectedVertex || vi_Visited [u] == i_SelectedVertex ) 01768 { 01769 // we already visited or pointing to itself. skip. 01770 continue; 01771 } 01772 01773 // move the last element in this bucket to u's position to get rid of expensive erase operation 01774 if(vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].size() > 1) 01775 { 01776 l = vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].back(); 01777 01778 vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]][vi_VertexLocation[u]] = l; 01779 01780 vi_VertexLocation[l] = vi_VertexLocation[u]; 01781 } 01782 01783 // remove the last element from vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]] 01784 vvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].pop_back(); 01785 01786 // now update that we are visiting u 01787 vi_Visited [u] = i_SelectedVertex; 01788 // end now update that we are visiting u 01789 01790 // now we tell that neighbor to increase its degree by 1 01791 ++vi_IncidenceVertexDegree [u]; 01792 // place the neighbor in 1 degree up in the bucket 01793 vvi_GroupedIncidenceVertexDegree [vi_IncidenceVertexDegree [u]].push_back ( u ); 01794 // update the location of the now moved neighbor 01795 vi_VertexLocation [u] = vvi_GroupedIncidenceVertexDegree [vi_IncidenceVertexDegree [u]].size () -1; 01796 } 01797 } 01798 // end tell the neighbors of i_SelectedVertex that we are removing it 01799 01800 // now we say that i_SelectedVertex's degree is unknown to us 01801 vi_IncidenceVertexDegree [i_SelectedVertex] = _UNKNOWN; 01802 // now we push i_SelectedVertex to the back or VertexOrder 01803 m_vi_OrderedVertices.push_back ( i_SelectedVertex + i_LeftVertexCount ); 01804 // loop it again 01805 ++i_SelectedVertexCount; 01806 } 01807 01808 // clear the buffer 01809 vi_IncidenceVertexDegree.clear(); 01810 vi_VertexLocation.clear(); 01811 vvi_GroupedIncidenceVertexDegree.clear(); 01812 01813 return(_TRUE); 01814 } 01815 01816 01817 //Public Function 2363 01818 string BipartiteGraphPartialOrdering::GetVertexOrderingVariant() 01819 { 01820 if(m_s_VertexOrderingVariant.compare("ROW_NATURAL") == 0) 01821 { 01822 return("Row Natural"); 01823 } 01824 else 01825 if(m_s_VertexOrderingVariant.compare("COLUMN_NATURAL") == 0) 01826 { 01827 return("Column Natural"); 01828 } 01829 else 01830 if(m_s_VertexOrderingVariant.compare("ROW_LARGEST_FIRST") == 0) 01831 { 01832 return("Row Largest First"); 01833 } 01834 else 01835 if(m_s_VertexOrderingVariant.compare("COLUMN_LARGEST_FIRST") == 0) 01836 { 01837 return("Column Largest First"); 01838 } 01839 else 01840 if(m_s_VertexOrderingVariant.compare("ROW_SMALLEST_LAST") == 0) 01841 { 01842 return("Row Smallest Last"); 01843 } 01844 else 01845 if(m_s_VertexOrderingVariant.compare("COLUMN_SMALLEST_LAST") == 0) 01846 { 01847 return("Column Smallest Last"); 01848 } 01849 else 01850 if(m_s_VertexOrderingVariant.compare("ROW_INCIDENCE_DEGREE") == 0) 01851 { 01852 return("Row Incidence Degree"); 01853 } 01854 else 01855 if(m_s_VertexOrderingVariant.compare("COLUMN_INCIDENCE_DEGREE") == 0) 01856 { 01857 return("Column Incidence Degree"); 01858 } 01859 else 01860 { 01861 return("Unknown"); 01862 } 01863 } 01864 01865 //Public Function 2364 01866 void BipartiteGraphPartialOrdering::GetOrderedVertices(vector<int> &output) 01867 { 01868 output = (m_vi_OrderedVertices); 01869 } 01870 01871 int BipartiteGraphPartialOrdering::OrderVertices(string s_OrderingVariant, string s_ColoringVariant) { 01872 s_ColoringVariant = toUpper(s_ColoringVariant); 01873 s_OrderingVariant = toUpper(s_OrderingVariant); 01874 01875 if(s_ColoringVariant == "ROW_PARTIAL_DISTANCE_TWO") 01876 { 01877 if((s_OrderingVariant.compare("NATURAL") == 0)) 01878 { 01879 return(RowNaturalOrdering()); 01880 } 01881 else 01882 if((s_OrderingVariant.compare("LARGEST_FIRST") == 0)) 01883 { 01884 return(RowLargestFirstOrdering()); 01885 } 01886 else 01887 if((s_OrderingVariant.compare("SMALLEST_LAST") == 0)) 01888 { 01889 return(RowSmallestLastOrdering()); 01890 } 01891 else 01892 if((s_OrderingVariant.compare("INCIDENCE_DEGREE") == 0)) 01893 { 01894 return(RowIncidenceDegreeOrdering()); 01895 } 01896 else 01897 if((s_OrderingVariant.compare("RANDOM") == 0)) 01898 { 01899 return(RowRandomOrdering()); 01900 } 01901 else 01902 { 01903 cerr<<endl; 01904 cerr<<"Unknown Ordering Method"; 01905 cerr<<endl; 01906 } 01907 } 01908 else 01909 if(s_ColoringVariant == "COLUMN_PARTIAL_DISTANCE_TWO") 01910 { 01911 if((s_OrderingVariant.compare("NATURAL") == 0)) 01912 { 01913 return(ColumnNaturalOrdering()); 01914 } 01915 else 01916 if((s_OrderingVariant.compare("LARGEST_FIRST") == 0)) 01917 { 01918 return(ColumnLargestFirstOrdering()); 01919 } 01920 else 01921 if((s_OrderingVariant.compare("SMALLEST_LAST") == 0)) 01922 { 01923 return(ColumnSmallestLastOrdering()); 01924 } 01925 else 01926 if((s_OrderingVariant.compare("INCIDENCE_DEGREE") == 0)) 01927 { 01928 return(ColumnIncidenceDegreeOrdering()); 01929 } 01930 else 01931 if((s_OrderingVariant.compare("RANDOM") == 0)) 01932 { 01933 return(ColumnRandomOrdering()); 01934 } 01935 else 01936 { 01937 cerr<<endl; 01938 cerr<<"Unknown Ordering Method: "<<s_OrderingVariant; 01939 cerr<<endl; 01940 } 01941 } 01942 else 01943 { 01944 cerr<<endl; 01945 cerr<<"Invalid s_ColoringVariant = \""<<s_ColoringVariant<<"\", must be either \"COLUMN_PARTIAL_DISTANCE_TWO\" or \"ROW_PARTIAL_DISTANCE_TWO\"."; 01946 cerr<<endl; 01947 } 01948 01949 return(_TRUE); 01950 } 01951 01952 01953 void BipartiteGraphPartialOrdering::PrintVertexOrdering() { 01954 cout<<"PrintVertexOrdering() "<<m_s_VertexOrderingVariant<<endl; 01955 for(unsigned int i=0; i<m_vi_OrderedVertices.size();i++) { 01956 //printf("\t [%d] %d \n", i, m_vi_OrderedVertices[i]); 01957 cout<<"\t["<<setw(5)<<i<<"] "<<setw(5)<<m_vi_OrderedVertices[i]<<endl; 01958 } 01959 cout<<endl; 01960 } 01961 01962 double BipartiteGraphPartialOrdering::GetVertexOrderingTime() { 01963 return m_d_OrderingTime; 01964 } 01965 01966 }