ColPack
BipartiteGraphBicoloring/BipartiteGraphOrdering.cpp
Go to the documentation of this file.
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 3401
00028         int BipartiteGraphOrdering::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 3451
00045         BipartiteGraphOrdering::BipartiteGraphOrdering()
00046         {
00047                 Clear();
00048         }
00049 
00050 
00051         //Public Destructor 3452
00052         BipartiteGraphOrdering::~BipartiteGraphOrdering()
00053         {
00054                 Clear();
00055         }
00056 
00057 
00058         //Virtual Function 3453
00059         void BipartiteGraphOrdering::Clear()
00060         {
00061                 BipartiteGraphVertexCover::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 3454
00074         void BipartiteGraphOrdering::Reset()
00075         {
00076                 BipartiteGraphVertexCover::Reset();
00077 
00078                 m_d_OrderingTime = _UNKNOWN;
00079 
00080                 m_s_VertexOrderingVariant.clear();
00081 
00082                 m_vi_OrderedVertices.clear();
00083 
00084                 return;
00085         }
00086 
00087         int BipartiteGraphOrdering::NaturalOrdering()
00088         {
00089                 if(CheckVertexOrdering("NATURAL"))
00090                 {
00091                         return(_TRUE);
00092                 }
00093 
00094                 int i;
00095 
00096                 int i_LeftVertexCount, i_RightVertexCount;
00097 
00098                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00099                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00100 
00101                 m_vi_OrderedVertices.clear();
00102                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00103 
00104                 for(i=0; i<i_LeftVertexCount; i++)
00105                 {
00106                         m_vi_OrderedVertices.push_back(i);
00107                 }
00108 
00109                 for(i=0; i<i_RightVertexCount; i++)
00110                 {
00111                         m_vi_OrderedVertices.push_back(i + i_LeftVertexCount);
00112                 }
00113 
00114                 return(_TRUE);
00115         }
00116 
00117         int BipartiteGraphOrdering::RandomOrdering()
00118         {
00119                 if(CheckVertexOrdering("RANDOM"))
00120                 {
00121                         return(_TRUE);
00122                 }
00123 
00124                 m_s_VertexOrderingVariant = "RANDOM";
00125 
00126                 int i;
00127 
00128                 int i_LeftVertexCount, i_RightVertexCount;
00129 
00130                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00131                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00132 
00133                 m_vi_OrderedVertices.clear();
00134 
00135                 //Order left vertices
00136                 m_vi_OrderedVertices.resize((unsigned) i_LeftVertexCount);
00137 
00138                 for(unsigned int i = 0; i<i_LeftVertexCount; i++) {
00139                         m_vi_OrderedVertices[i] = i;
00140                 }
00141 
00142                 randomOrdering(m_vi_OrderedVertices);
00143 
00144                 //Order right vertices
00145                 vector<int> tempOrdering;
00146 
00147                 tempOrdering.resize((unsigned) i_RightVertexCount);
00148 
00149                 for(unsigned int i = 0; i<i_RightVertexCount; i++) {
00150                         tempOrdering[i] = i + i_LeftVertexCount;
00151                 }
00152 
00153                 randomOrdering(tempOrdering);
00154                 
00155                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00156 
00157                 //Now, populate vector m_vi_OrderedVertices with the right vertices
00158                 for(unsigned int i = 0; i<i_RightVertexCount; i++) {
00159                         m_vi_OrderedVertices.push_back(tempOrdering[i]);
00160                 }
00161                 
00162                 return(_TRUE);
00163         }
00164 
00165         int BipartiteGraphOrdering::LargestFirstOrdering()
00166         {
00167                 if(CheckVertexOrdering("LARGEST_FIRST"))
00168                 {
00169                         return(_TRUE);
00170                 }
00171 
00172                 int i, j;
00173 
00174                 int i_LeftVertexCount, i_RightVertexCount;
00175 
00176                 int i_HighestDegreeVertex;
00177 
00178                 int i_VertexDegree, i_VertexDegreeCount;
00179 
00180                 vector< vector< int > > vvi_GroupedVertexDegree;
00181 
00182                 m_i_MaximumVertexDegree = _FALSE;
00183 
00184                 i_HighestDegreeVertex = _UNKNOWN;
00185 
00186                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00187                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00188 
00189                 vvi_GroupedVertexDegree.clear();
00190                 vvi_GroupedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00191 
00192                 for(i=0; i<i_LeftVertexCount; i++)
00193                 {
00194                         i_VertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00195 
00196                         vvi_GroupedVertexDegree[i_VertexDegree].push_back(i);
00197 
00198                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00199                         {
00200                                 m_i_MaximumVertexDegree = i_VertexDegree;
00201 
00202                                 i_HighestDegreeVertex = i;
00203                         }
00204                 }
00205 
00206                 for(i=0; i<i_RightVertexCount; i++)
00207                 {
00208                         i_VertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00209 
00210                         vvi_GroupedVertexDegree[i_VertexDegree].push_back(i + i_LeftVertexCount);
00211 
00212                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00213                         {
00214                                 m_i_MaximumVertexDegree = i_VertexDegree;
00215 
00216                                 i_HighestDegreeVertex = i + i_LeftVertexCount;
00217                         }
00218                 }
00219 
00220                 m_vi_OrderedVertices.clear();
00221                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00222 
00223                 if(i_HighestDegreeVertex < i_LeftVertexCount)
00224                 {
00225                         for(i=m_i_MaximumVertexDegree; i>=0; i--)
00226                         {
00227                                 i_VertexDegreeCount = (signed) vvi_GroupedVertexDegree[i].size();
00228 
00229                                 for(j=0; j<i_VertexDegreeCount; j++)
00230                                 {
00231                                         m_vi_OrderedVertices.push_back(vvi_GroupedVertexDegree[i][j]);
00232                                 }
00233                         }
00234                 }
00235                 else
00236                 {
00237                         for(i=m_i_MaximumVertexDegree; i>=0; i--)
00238                         {
00239                                 i_VertexDegreeCount = (signed) vvi_GroupedVertexDegree[i].size();
00240 
00241                                 for(j=STEP_DOWN(i_VertexDegreeCount); j>=0; j--)
00242                                 {
00243                                         m_vi_OrderedVertices.push_back(vvi_GroupedVertexDegree[i][j]);
00244                                 }
00245                         }
00246                 }
00247 
00248                 vvi_GroupedVertexDegree.clear();
00249 
00250                 return(_TRUE);
00251         }
00252 
00253         int BipartiteGraphOrdering::SmallestLastOrdering()
00254         {
00255                 if(CheckVertexOrdering("SMALLEST_LAST"))
00256                 {
00257                         return(_TRUE);
00258                 }
00259 
00260                 int i, u, l, v;
00261 
00262                 int _FOUND;
00263 
00264                 int i_HighestInducedVertexDegree, i_HighestInducedDegreeVertex;
00265 
00266                 int i_LeftVertexCount, i_RightVertexCount;
00267 
00268                 int i_VertexCountMinus1; // = i_LeftVertexCount + i_RightVertexCount - 1, used when inserting selected vertices into m_vi_OrderedVertices
00269 
00270                 int i_InducedVertexDegree;
00271 
00272                 int i_InducedVertexDegreeCount;
00273 
00274                 int i_SelectedVertex, i_SelectedVertexCount;
00275 
00276                 vector <int> vi_InducedVertexDegree;
00277 
00278                 vector < vector < int > > vvi_GroupedInducedVertexDegree;
00279 
00280                 vector <  int > vi_VertexLocation;
00281                 
00282                 vector <  int > vi_LeftSidedVertexinBucket;
00283 
00284 
00285                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00286                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00287 
00288                 i_VertexCountMinus1 = i_LeftVertexCount + i_RightVertexCount - 1;
00289 
00290                 vi_InducedVertexDegree.clear();
00291                 vi_InducedVertexDegree.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00292 
00293                 vvi_GroupedInducedVertexDegree.clear();
00294                 vvi_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00295 
00296                 vi_VertexLocation.clear();
00297                 vi_VertexLocation.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00298 
00299                 vi_LeftSidedVertexinBucket.clear();
00300                 vi_LeftSidedVertexinBucket.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00301 
00302                 i_HighestInducedVertexDegree = _FALSE;
00303 
00304                 i_HighestInducedDegreeVertex = _UNKNOWN;
00305 
00306                 i_SelectedVertex = _UNKNOWN;
00307 
00308                 for(i=0; i<i_LeftVertexCount; i++)
00309                 {
00310                         i_InducedVertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00311 
00312                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00313 
00314                         vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].push_back(i);
00315 
00316                         vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1);
00317 
00318                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00319                         {
00320                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00321 
00322                                 i_HighestInducedDegreeVertex = i;
00323                         }
00324                 }
00325 
00326                 
00327                 // get the bucket division positions now
00328                 for(i= 0; i < i_LeftVertexCount + i_RightVertexCount; i++)
00329                         vi_LeftSidedVertexinBucket.push_back(vvi_GroupedInducedVertexDegree[i].size());
00330         
00331 
00332                 for(i=0; i<i_RightVertexCount; i++)
00333                 {
00334                         i_InducedVertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00335 
00336                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00337 
00338                         vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].push_back(i + i_LeftVertexCount);
00339 
00340                         vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1);
00341 
00342                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00343                         {
00344                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00345 
00346                                 i_HighestInducedDegreeVertex = i + i_LeftVertexCount;
00347                         }
00348                 }
00349 
00350                 m_vi_OrderedVertices.clear();
00351                 m_vi_OrderedVertices.resize(i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
00352 
00353                 i_SelectedVertexCount = _FALSE;
00354 
00355                 int iMin = 1;
00356 
00357                 while(i_SelectedVertexCount < i_LeftVertexCount + i_RightVertexCount)
00358                 {
00359                         if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin -1].size() != _FALSE)
00360                                 iMin--;
00361 
00362                         for(i=iMin; i<STEP_UP(i_HighestInducedVertexDegree); i++)
00363                         {
00364                                 i_InducedVertexDegreeCount = (signed) vvi_GroupedInducedVertexDegree[i].size();
00365 
00366                                 if(i_InducedVertexDegreeCount == _FALSE)
00367                                 {
00368                                         iMin++;
00369                                         continue;
00370                                 }
00371 
00372                                 if(i_HighestInducedDegreeVertex < i_LeftVertexCount)
00373                                 {
00374                                         _FOUND = _FALSE;
00375 
00376                                         /*
00377                                         if(vi_LeftSidedVertexinBucket[i] > 0)
00378                                         {
00379                                                 vi_LeftSidedVertexinBucket[i]--;
00380                                                 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i][vi_LeftSidedVertexinBucket[i]];
00381 
00382                                                 vvi_GroupedInducedVertexDegree[i][vi_LeftSidedVertexinBucket[i]] = vvi_GroupedInducedVertexDegree[i].back();
00383                                                 vi_VertexLocation[vvi_GroupedInducedVertexDegree[i].back()] = vi_VertexLocation[u];
00384 
00385                                                 _FOUND = _TRUE;                                         
00386                                         }
00387                                         */
00388                                         if(vi_LeftSidedVertexinBucket[i] > 0)
00389                                         for(int j  = 0; j < vvi_GroupedInducedVertexDegree[i].size(); j++)
00390                                         {
00391                                                 u = vvi_GroupedInducedVertexDegree[i][j];
00392                                                 if(u < i_LeftVertexCount)
00393                                                 {
00394                                                         i_SelectedVertex = u;
00395 
00396                                                         if(vvi_GroupedInducedVertexDegree[i].size() > 1)
00397                                                         {
00398                                                                 // swap this node with the last node
00399                                                                 vvi_GroupedInducedVertexDegree[i][j] = vvi_GroupedInducedVertexDegree[i].back();
00400                                                                 vi_VertexLocation[vvi_GroupedInducedVertexDegree[i].back()] = vi_VertexLocation[u];
00401                                                         }
00402                                                         _FOUND = _TRUE;
00403                                                         vi_LeftSidedVertexinBucket[i]--;
00404 
00405                                                         break;
00406                                                 }
00407                                         }
00408 
00409                                         if(!_FOUND)
00410                                                 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back();
00411                                         
00412                                         break;
00413                                 }
00414                                 else
00415                                 {
00416                                         _FOUND = _FALSE;
00417 
00418                                         if((i_InducedVertexDegreeCount - vi_LeftSidedVertexinBucket[i]) > 0)
00419                                         for(int j = 0; j < vvi_GroupedInducedVertexDegree[i].size(); j++)
00420                                         {
00421                                                 u = vvi_GroupedInducedVertexDegree[i][j];
00422 
00423                                                 if(u >= i_LeftVertexCount)
00424                                                 {
00425                                                         i_SelectedVertex = u;
00426                                                         if(vvi_GroupedInducedVertexDegree[i].size() > 1)
00427                                                         {
00428                                                                 vvi_GroupedInducedVertexDegree[i][j] = vvi_GroupedInducedVertexDegree[i].back();
00429                                                                 vi_VertexLocation[vvi_GroupedInducedVertexDegree[i].back()] = vi_VertexLocation[u];
00430                                                         }
00431                                                         _FOUND = _TRUE;
00432 
00433                                                         break;
00434                                                 }
00435                                         }
00436 
00437                                         if(!_FOUND)
00438                                         {
00439                                                 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back();
00440                                                 vi_LeftSidedVertexinBucket[i]--;
00441                                         }
00442                                 }
00443 
00444                                 break;
00445                         }
00446 
00447                         vvi_GroupedInducedVertexDegree[i].pop_back(); // remove the selected vertex from the bucket
00448 
00449                         if(i_SelectedVertex < i_LeftVertexCount)
00450                         {
00451                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
00452                                 {
00453                                         u = m_vi_Edges[i] + i_LeftVertexCount; // neighbour are always right sided
00454 
00455                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00456                                         {
00457                                                 continue;
00458                                         }
00459                                         
00460                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00461                                         if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1)
00462                                         {
00463                                                 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back();
00464 
00465                                                 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l;
00466 
00467                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00468                                         }
00469                                         // remove last element from this bucket
00470                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back();
00471 
00472 
00473                                         // reduce degree of u by 1
00474                                         vi_InducedVertexDegree[u]--;
00475 
00476                                         // move u to appropriate bucket
00477                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u);
00478 
00479                                         // update vi_VertexLocation[u] since it has now been changed
00480                                          vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00481                                         
00482                                         /*
00483                                         if(u < i_LeftVertexCount)
00484                                         {
00485                                                 // swap this vertex and location                                                
00486                                                 v = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00487                                                 if(v > 0)
00488                                                 {
00489                                                         l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v];
00490 
00491                                                         swap(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]], vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v]);
00492                                                         swap(vi_VertexLocation[u], vi_VertexLocation[l]);
00493                                                 }
00494                                                 vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]++;
00495                                         }*/
00496                                 }
00497                         }
00498                         else
00499                         {
00500                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
00501                                 {
00502                                         u = m_vi_Edges[i]; // neighbour are always left sided
00503 
00504                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00505                                         {
00506                                                 continue;
00507                                         }
00508                                         
00509                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00510                                         if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1)
00511                                         {
00512                                                 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back();
00513 
00514                                                 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l;
00515 
00516                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00517                                         }
00518                                         
00519                                         // remove last element from this bucket
00520                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back();
00521                                         
00522                                         vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]--;
00523 
00524                                         // reduce degree of u by 1
00525                                         vi_InducedVertexDegree[u]--;
00526 
00527                                         vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]++;
00528 
00529                                         // move u to appropriate bucket
00530                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u);
00531 
00532                                         // update vi_VertexLocation[u] since it has now been changed
00533                                         vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00534 
00535                                         /*
00536                                         if(u < i_LeftVertexCount)
00537                                         {
00538                                                 // swap this vertex and location                                                
00539                                                 v = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00540                                                 if(v > 0)
00541                                                 {
00542                                                         l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v];
00543 
00544                                                         swap(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]], vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v]);
00545                                                         swap(vi_VertexLocation[u], vi_VertexLocation[l]);
00546                                                 }
00547                                                 vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]++;
00548                                         }*/
00549 
00550                                 }
00551                         }
00552 
00553                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
00554 
00555                         m_vi_OrderedVertices[i_VertexCountMinus1 - i_SelectedVertexCount] = i_SelectedVertex;
00556 
00557                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
00558                 }
00559 
00560                 vi_InducedVertexDegree.clear();
00561                 vvi_GroupedInducedVertexDegree.clear();
00562                 vi_VertexLocation.clear();
00563                 vi_LeftSidedVertexinBucket.clear();
00564 
00565                 return(_TRUE);
00566         }
00567 
00568         int BipartiteGraphOrdering::IncidenceDegreeOrdering()
00569         {
00570                 if(CheckVertexOrdering("INCIDENCE_DEGREE"))
00571                 {
00572                         return(_TRUE);
00573                 }
00574 
00575                 int i, u, l;
00576                 
00577                 int i_HighestIncidenceVertexDegree;
00578 
00579                 int i_LeftVertexCount, i_RightVertexCount, i_VertexCount;
00580 
00581                 int i_VertexDegree;
00582 
00583                 int i_IncidenceVertexDegree, i_IncidenceVertexDegreeCount;
00584 
00585                 int i_SelectedVertex, i_SelectedVertexCount;
00586 
00587                 vector<int> vi_IncidenceVertexDegree;
00588 
00589                 //Vertices of the same IncidenceDegree are differenciated into 
00590                 //  LeftVertices (vpvi_GroupedIncidenceVertexDegree.first) and
00591                 //  RightVertices (vpvi_GroupedIncidenceVertexDegree.second)
00592                 vector< pair<vector<int>, vector<int> > > vpvi_GroupedIncidenceVertexDegree;
00593 
00594                 vector< int > vi_VertexLocation;
00595 
00596                 list<int>::iterator lit_ListIterator; //???
00597 
00598                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00599                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00600                 i_VertexCount = i_LeftVertexCount + i_RightVertexCount;
00601 
00602                 vi_IncidenceVertexDegree.clear();
00603                 vi_IncidenceVertexDegree.reserve((unsigned) (i_VertexCount));
00604 
00605                 vpvi_GroupedIncidenceVertexDegree.clear();
00606                 vpvi_GroupedIncidenceVertexDegree.resize((unsigned) (i_VertexCount));
00607 
00608                 vi_VertexLocation.clear();
00609                 vi_VertexLocation.reserve((unsigned) (i_VertexCount));
00610 
00611                 i_HighestIncidenceVertexDegree = _UNKNOWN;
00612 
00613                 i_IncidenceVertexDegree = _FALSE;
00614 
00615                 i_SelectedVertex = _UNKNOWN;
00616 
00617                 for(i=0; i<i_LeftVertexCount; i++)
00618                 {
00619                         vi_IncidenceVertexDegree.push_back(i_IncidenceVertexDegree);
00620 
00621                         vpvi_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].first.push_back(i);
00622 
00623                         vi_VertexLocation.push_back(vpvi_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].first.size() - 1);
00624 
00625                         i_VertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00626 
00627                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00628                         {
00629                                 m_i_MaximumVertexDegree = i_VertexDegree;
00630                         }
00631                 }
00632 
00633                 for(i=0; i<i_RightVertexCount; i++)
00634                 {
00635                         vi_IncidenceVertexDegree.push_back(i_IncidenceVertexDegree);
00636 
00637                         vpvi_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].second.push_back(i + i_LeftVertexCount);
00638 
00639                         vi_VertexLocation.push_back(vpvi_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].second.size() - 1);
00640 
00641                         i_VertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00642 
00643                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00644                         {
00645                                 m_i_MaximumVertexDegree = i_VertexDegree;
00646                         }
00647                 }
00648                 
00649                 i_HighestIncidenceVertexDegree = 0;
00650 
00651                 m_vi_OrderedVertices.clear();
00652                 m_vi_OrderedVertices.reserve((unsigned) (i_VertexCount));
00653 
00654                 i_SelectedVertexCount = _FALSE;
00655 
00656                 while(i_SelectedVertexCount < i_VertexCount)
00657                 {
00658                         if(i_HighestIncidenceVertexDegree != m_i_MaximumVertexDegree &&
00659                             vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree+1].first.size() + 
00660                             vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree+1].second.size()  != 0) {
00661                           //We need to update the value of i_HighestIncidenceVertexDegree
00662                           i_HighestIncidenceVertexDegree++;
00663                         
00664                         }
00665                         else {
00666                           while(vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].first.size() + 
00667                             vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].second.size()  == 0) {
00668                             i_HighestIncidenceVertexDegree--;
00669                           }
00670                         }
00671                         
00672                         if(vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].first.size() != 0) {
00673                           //vertex with i_HighestIncidenceVertexDegree is a LeftVertex
00674                           i_SelectedVertex = vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].first.back();
00675                           vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].first.pop_back();
00676                         }
00677                         else {
00678                           //vertex with i_HighestIncidenceVertexDegree is a RightVertex
00679                           i_SelectedVertex = vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].second.back();
00680                           vpvi_GroupedIncidenceVertexDegree[i_HighestIncidenceVertexDegree].second.pop_back();
00681                         }
00682 
00683                         // Increase the IncidenceDegree of all the unvisited neighbor vertices by 1 and move them to the correct buckets
00684                         if(i_SelectedVertex < i_LeftVertexCount) // i_SelectedVertex is a LeftVertex
00685                         {
00686                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
00687                                 {
00688                                         u = m_vi_Edges[i] + i_LeftVertexCount;
00689                                         if(vi_IncidenceVertexDegree[u] == _UNKNOWN)
00690                                         {
00691                                                 continue;
00692                                         }
00693                                         
00694                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00695                                         if(vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second.size() > 1) {
00696                                                 l = vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second.back();
00697                                                 vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second[vi_VertexLocation[u]] = l;
00698                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00699                                         }
00700                                         
00701                                         //remove the last element from vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second
00702                                         vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second.pop_back();
00703                                         
00704                                         // increase incidence degree of u
00705                                         vi_IncidenceVertexDegree[u]++;
00706 
00707                                         // insert u into appropriate bucket
00708                                         vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second.push_back(u);
00709                                         
00710                                         // update location of u
00711                                         vi_VertexLocation[u] = vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].second.size() - 1;
00712                                 }
00713                         }
00714                         else
00715                         {
00716                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
00717                                 {
00718                                         u = m_vi_Edges[i];
00719                                         if(vi_IncidenceVertexDegree[u] == _UNKNOWN)
00720                                         {
00721                                                 continue;
00722                                         }
00723                                         
00724                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00725                                         if(vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first.size() > 1) {
00726                                                 l = vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first.back();
00727                                                 vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first[vi_VertexLocation[u]] = l;
00728                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00729                                         }
00730                                         
00731                                         //remove the last element from vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first
00732                                         vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first.pop_back();
00733                                         
00734                                         // increase incidence degree of u
00735                                         vi_IncidenceVertexDegree[u]++;
00736 
00737                                         // insert u into appropriate bucket
00738                                         vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first.push_back(u);
00739                                         
00740                                         // update location of u
00741                                         vi_VertexLocation[u] = vpvi_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[u]].first.size() - 1;
00742                                 }
00743 
00744                         }
00745         
00746                         // Mark that this vertex has been visited
00747                         vi_IncidenceVertexDegree[i_SelectedVertex] = _UNKNOWN;
00748 
00749                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
00750 
00751                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
00752                 }
00753 
00754 #if DEBUG == 3458
00755 
00756                 int i_OrderedVertexCount;
00757 
00758                 cout<<endl;
00759                 cout<<"DEBUG 3458 | Bipartite Graph Coloring | Bipartite Incidence Degree Ordering"<<endl;
00760                 cout<<endl;
00761 
00762                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
00763 
00764                 for(i=0; i<i_OrderedVertexCount; i++)
00765                 {
00766                         if(i == STEP_DOWN(i_OrderedVertexCount))
00767                         {
00768                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_OrderedVertexCount<<")"<<endl;
00769                         }
00770                         else
00771                         {
00772                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
00773                         }
00774                 }
00775 
00776                 cout<<endl;
00777                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_VertexCount<<"]"<<endl;
00778                 cout<<endl;
00779 
00780 #endif
00781 
00782                 return(_TRUE);
00783         }
00784 
00785 
00786         int BipartiteGraphOrdering::DynamicLargestFirstOrdering()
00787         {
00788                 if(CheckVertexOrdering("DYNAMIC_LARGEST_FIRST"))
00789                 {
00790                         return(_TRUE);
00791                 }
00792 
00793                 int i, u, l;
00794 
00795                 int i_HighestInducedVertexDegree;
00796 
00797                 int i_LeftVertexCount, i_RightVertexCount, i_VertexCount;
00798 
00799                 int i_InducedVertexDegree;
00800 
00801                 int i_SelectedVertex, i_SelectedVertexCount;
00802 
00803                 vector<int> vi_InducedVertexDegree;
00804 
00805                 vector< pair<vector<int>, vector<int> > > vpvi_GroupedInducedVertexDegree;
00806 
00807                 vector< int > vi_VertexLocation;
00808                 
00809 
00810                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00811                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00812                 i_VertexCount = i_LeftVertexCount + i_RightVertexCount;
00813 
00814                 vi_InducedVertexDegree.clear();
00815                 vi_InducedVertexDegree.reserve((unsigned) i_VertexCount);
00816 
00817                 vpvi_GroupedInducedVertexDegree.clear();
00818                 vpvi_GroupedInducedVertexDegree.resize((unsigned) i_VertexCount);
00819 
00820                 vi_VertexLocation.clear();
00821                 vi_VertexLocation.reserve((unsigned) i_VertexCount);
00822 
00823                 i_SelectedVertex = _UNKNOWN;
00824 
00825                 for(i=0; i<i_LeftVertexCount; i++)
00826                 {
00827                         i_InducedVertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00828 
00829                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00830 
00831                         vpvi_GroupedInducedVertexDegree[i_InducedVertexDegree].first.push_back(i);
00832 
00833                         vi_VertexLocation.push_back(vpvi_GroupedInducedVertexDegree[i_InducedVertexDegree].first.size() - 1);
00834 
00835                         if(m_i_MaximumVertexDegree < i_InducedVertexDegree)
00836                         {
00837                                 m_i_MaximumVertexDegree = i_InducedVertexDegree;
00838                         }
00839                 }
00840 
00841                 for(i=0; i<i_RightVertexCount; i++)
00842                 {
00843                         i_InducedVertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00844 
00845                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00846 
00847                         vpvi_GroupedInducedVertexDegree[i_InducedVertexDegree].second.push_back(i + i_LeftVertexCount);
00848 
00849                         vi_VertexLocation.push_back(vpvi_GroupedInducedVertexDegree[i_InducedVertexDegree].second.size() - 1);
00850 
00851                         if(m_i_MaximumVertexDegree < i_InducedVertexDegree)
00852                         {
00853                                 m_i_MaximumVertexDegree = i_InducedVertexDegree;
00854                         }
00855                 }
00856                 
00857                 i_HighestInducedVertexDegree = m_i_MaximumVertexDegree;
00858 
00859                 m_vi_OrderedVertices.clear();
00860                 m_vi_OrderedVertices.reserve((unsigned) i_VertexCount);
00861 
00862                 i_SelectedVertexCount = _FALSE;
00863 
00864                 // just counting the number of vertices that we have worked with,
00865                 // stop when i_SelectedVertexCount == i_VertexCount, i.e. we have looked through all the vertices
00866                 while(i_SelectedVertexCount < i_VertexCount)
00867                 {
00868                         while(vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].first.size() + 
00869                             vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].second.size()  == 0) {
00870                           i_HighestInducedVertexDegree--;
00871                         }
00872                         
00873                         if(vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].first.size() != 0) {
00874                           //vertex with i_HighestInducedVertexDegree is a LeftVertex
00875                           i_SelectedVertex = vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].first.back();
00876                           vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].first.pop_back();
00877                         }
00878                         else {
00879                           //vertex with i_HighestInducedVertexDegree is a RightVertex
00880                           i_SelectedVertex = vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].second.back();
00881                           vpvi_GroupedInducedVertexDegree[i_HighestInducedVertexDegree].second.pop_back();
00882                         }
00883 
00884                         // Decrease the InducedVertexDegree of all the unvisited neighbor vertices by 1 and move them to the correct buckets
00885                         if(i_SelectedVertex < i_LeftVertexCount) // i_SelectedVertex is a LeftVertex
00886                         {
00887                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
00888                                 {
00889                                         u = m_vi_Edges[i] + i_LeftVertexCount;
00890                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00891                                         {
00892                                                 continue;
00893                                         }
00894 
00895                                         
00896                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00897                                         if(vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second.size() > 1) {
00898                                                 l = vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second.back();
00899                                                 vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second[vi_VertexLocation[u]] = l;
00900                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00901                                         }
00902                                         
00903                                         //remove the last element from vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second
00904                                         vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second.pop_back();
00905                                         
00906                                         // increase incidence degree of u
00907                                         vi_InducedVertexDegree[u]--;
00908 
00909                                         // insert u into appropriate bucket
00910                                         vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second.push_back(u);
00911                                         
00912                                         // update location of u
00913                                         vi_VertexLocation[u] = vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].second.size() - 1;
00914                                 }
00915                         }
00916                         else
00917                         {
00918                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
00919                                 {
00920                                         u = m_vi_Edges[i];
00921                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00922                                         {
00923                                                 continue;
00924                                         }
00925 
00926                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00927                                         if(vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first.size() > 1) {
00928                                                 l = vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first.back();
00929                                                 vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first[vi_VertexLocation[u]] = l;
00930                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00931                                         }
00932                                         
00933                                         //remove the last element from vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first
00934                                         vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first.pop_back();
00935                                         
00936                                         // increase incidence degree of u
00937                                         vi_InducedVertexDegree[u]--;
00938 
00939                                         // insert u into appropriate bucket
00940                                         vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first.push_back(u);
00941                                         
00942                                         // update location of u
00943                                         vi_VertexLocation[u] = vpvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].first.size() - 1;
00944                                 }
00945                         }
00946 
00947                         // Mark that this vertex has been visited
00948                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
00949 
00950                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
00951 
00952                         i_SelectedVertexCount++;
00953                 }
00954 
00955 
00956 #if DEBUG == 3462
00957 
00958                 int i_OrderedVertexCount;
00959 
00960                 cout<<endl;
00961                 cout<<"DEBUG 3462 | Bipartite Graph Coloring | Bipartite Dynamic Largest First Ordering"<<endl;
00962                 cout<<endl;
00963 
00964                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
00965 
00966                 for(i=0; i<i_OrderedVertexCount; i++)
00967                 {
00968                         if(i == STEP_DOWN(i_OrderedVertexCount))
00969                         {
00970                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_OrderedVertexCount<<")"<<endl;
00971                         }
00972                         else
00973                         {
00974                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
00975                         }
00976                 }
00977 
00978                 cout<<endl;
00979                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_VertexCount<<"]"<<endl;
00980                 cout<<endl;
00981 
00982 #endif
00983 
00984                 return(_TRUE);
00985         }
00986 
00987         int BipartiteGraphOrdering::SelectiveLargestFirstOrdering()
00988         {
00989                 if(CheckVertexOrdering("SELECTVE_LARGEST_FIRST"))
00990                 {
00991                         return(_TRUE);
00992                 }
00993 
00994                 int i, j;
00995 
00996                 int i_LeftVertexCount, i_RightVertexCount;
00997 
00998                 int i_VertexDegree, i_VertexDegreeCount;
00999 
01000                 vector< vector<int> > vvi_GroupedVertexDegree;
01001 
01002                 m_i_MaximumVertexDegree = _FALSE;
01003 
01004                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01005                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01006 
01007                 vvi_GroupedVertexDegree.clear();
01008                 vvi_GroupedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01009 
01010                 for(i=0; i<i_LeftVertexCount; i++)
01011                 {
01012                         if(m_vi_IncludedLeftVertices[i] == _FALSE)
01013                         {
01014                                 continue;
01015                         }
01016 
01017                         i_VertexDegree = _FALSE;
01018 
01019                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01020                         {
01021                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01022                                 {
01023                                         continue;
01024                                 }
01025 
01026                                 i_VertexDegree++;
01027                         }
01028 
01029                         vvi_GroupedVertexDegree[i_VertexDegree].push_back(i);
01030 
01031                         if(m_i_MaximumVertexDegree < i_VertexDegree)
01032                         {
01033                                 m_i_MaximumVertexDegree = i_VertexDegree;
01034                         }
01035                 }
01036 
01037                 for(i=0; i<i_RightVertexCount; i++)
01038                 {
01039                         if(m_vi_IncludedRightVertices[i] == _FALSE)
01040                         {
01041                                 continue;
01042                         }
01043 
01044                         i_VertexDegree = _FALSE;
01045 
01046                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01047                         {
01048                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01049                                 {
01050                                         continue;
01051                                 }
01052 
01053                                 i_VertexDegree++;
01054                         }
01055 
01056                         vvi_GroupedVertexDegree[i_VertexDegree].push_back(i + i_LeftVertexCount);
01057 
01058                         if(m_i_MaximumVertexDegree < i_VertexDegree)
01059                         {
01060                                 m_i_MaximumVertexDegree = i_VertexDegree;
01061                         }
01062                 }
01063 
01064                 m_vi_OrderedVertices.clear();
01065 
01066                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01067                 {
01068                         i_VertexDegreeCount = (signed) vvi_GroupedVertexDegree[i].size();
01069 
01070                         for(j=0; j<i_VertexDegreeCount; j++)
01071                         {
01072                                 m_vi_OrderedVertices.push_back(vvi_GroupedVertexDegree[i][j]);
01073                         }
01074                 }
01075 
01076 #if DEBUG == 3459
01077 
01078                 int i_VertexCount;
01079 
01080                 cout<<endl;
01081                 cout<<"DEBUG 3459 | Bipartite Graph Bicoloring | Largest First Ordering"<<endl;
01082                 cout<<endl;
01083 
01084                 i_VertexCount = (signed) m_vi_OrderedVertices.size();
01085 
01086                 for(i=0; i<i_VertexCount; i++)
01087                 {
01088                         if(i == STEP_DOWN(i_VertexCount))
01089                         {
01090                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_VertexCount<<")"<<endl;
01091                         }
01092                         else
01093                         {
01094                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
01095                         }
01096                 }
01097 
01098                 cout<<endl;
01099                 cout<<"[Highest Vertex Degree = "<<m_i_MaximumVertexDegree<<"]"<<endl;
01100                 cout<<endl;
01101 
01102 #endif
01103 
01104                 return(_TRUE);
01105         }
01106 
01107         int BipartiteGraphOrdering::SelectiveSmallestLastOrdering()
01108         {
01109                 if(CheckVertexOrdering("SELECTIVE_SMALLEST_LAST"))
01110                 {
01111                         return(_TRUE);
01112                 }
01113 
01114                 int i, j;
01115 
01116                 int i_HighestInducedVertexDegree;
01117 
01118                 int i_LeftVertexCount, i_RightVertexCount;
01119 
01120                 int i_InducedVertexDegree;
01121 
01122                 int i_InducedVertexDegreeCount;
01123 
01124                 int i_IncludedVertexCount;
01125 
01126                 int i_SelectedVertex, i_SelectedVertexCount;
01127 
01128                 vector<int> vi_InducedVertexDegree;
01129 
01130                 vector< list<int> > vli_GroupedInducedVertexDegree;
01131 
01132                 vector< list<int>::iterator > vlit_VertexLocation;
01133 
01134                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01135                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01136 
01137                 vi_InducedVertexDegree.clear();
01138                 vi_InducedVertexDegree.resize((signed) i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
01139 
01140                 vli_GroupedInducedVertexDegree.clear();
01141                 vli_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01142 
01143                 vlit_VertexLocation.clear();
01144                 vlit_VertexLocation.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01145 
01146                 i_IncludedVertexCount = _FALSE;
01147 
01148                 i_HighestInducedVertexDegree = _FALSE;
01149 
01150                 i_SelectedVertex = _UNKNOWN;
01151 
01152                 for(i=0; i<i_LeftVertexCount; i++)
01153                 {
01154                 if(m_vi_IncludedLeftVertices[i] == _FALSE)
01155                         {
01156                                 continue;
01157                         }
01158 
01159                         i_IncludedVertexCount++;
01160 
01161                         i_InducedVertexDegree = _FALSE;
01162 
01163                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01164                         {
01165                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01166                                 {
01167                                         continue;
01168                                 }
01169 
01170                                 i_InducedVertexDegree++;
01171                         }
01172 
01173                         vi_InducedVertexDegree[i] = i_InducedVertexDegree;
01174 
01175                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i);
01176 
01177                         vlit_VertexLocation[vli_GroupedInducedVertexDegree[i_InducedVertexDegree].front()] = vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin();
01178 
01179                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01180                         {
01181                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01182                         }
01183                 }
01184 
01185                 for(i=0; i<i_RightVertexCount; i++)
01186                 {
01187                 if(m_vi_IncludedRightVertices[i] == _FALSE)
01188                         {
01189                                 continue;
01190                         }
01191 
01192                         i_IncludedVertexCount++;
01193 
01194                         i_InducedVertexDegree = _FALSE;
01195 
01196                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01197                         {
01198                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01199                                 {
01200                                         continue;
01201                                 }
01202 
01203                                 i_InducedVertexDegree++;
01204                         }
01205 
01206                         vi_InducedVertexDegree[i + i_LeftVertexCount] = i_InducedVertexDegree;
01207 
01208                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i + i_LeftVertexCount);
01209 
01210                         vlit_VertexLocation[vli_GroupedInducedVertexDegree[i_InducedVertexDegree].front()] = vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin();
01211 
01212                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01213                         {
01214                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01215                         }
01216                 }
01217 
01218 
01219 #if DEBUG == 3460
01220 
01221                 list<int>::iterator lit_ListIterator;
01222 
01223                 cout<<endl;
01224                 cout<<"DEBUG 3460 | Vertex Ordering | Vertex Degree"<<endl;
01225                 cout<<endl;
01226 
01227                 for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
01228                 {
01229                         cout<<"Degree "<<i<<"\t"<<" : ";
01230 
01231                         i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01232 
01233                         j = _FALSE;
01234 
01235                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
01236                         {
01237                                 if(j==STEP_DOWN(i_InducedVertexDegreeCount))
01238                                 {
01239                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_InducedVertexDegreeCount<<")";
01240                                 }
01241                                 else
01242                                 {
01243                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
01244                                 }
01245 
01246                                 j++;
01247                         }
01248 
01249                         cout<<endl;
01250                 }
01251 
01252                 cout<<endl;
01253 
01254 #endif
01255 
01256                 m_vi_OrderedVertices.clear();
01257 
01258                 i_SelectedVertexCount = _FALSE;
01259 
01260                 while(i_SelectedVertexCount < i_IncludedVertexCount)
01261                 {
01262                         for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
01263                         {
01264                                 i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01265 
01266                                 if(i_InducedVertexDegreeCount != _FALSE)
01267                                 {
01268                                         i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
01269 
01270                                         break;
01271                                 }
01272                         }
01273 
01274                         if(i_SelectedVertex < i_LeftVertexCount)
01275                         {
01276                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01277                                 {
01278                                         if(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01279                                         {
01280                                                 continue;
01281                                         }
01282 
01283                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01284 
01285                                         vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01286 
01287                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01288 
01289                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01290                                 }
01291                         }
01292                         else
01293                         {
01294                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01295                                 {
01296                                         if(vi_InducedVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01297                                         {
01298                                                 continue;
01299                                         }
01300 
01301                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01302 
01303                                         vi_InducedVertexDegree[m_vi_Edges[i]] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i]]);
01304 
01305                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01306 
01307                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].begin();
01308                                 }
01309                         }
01310 
01311                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01312 
01313                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
01314 
01315                         m_vi_OrderedVertices.insert(m_vi_OrderedVertices.begin(), i_SelectedVertex);
01316 
01317                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01318                 }
01319 
01320 
01321 #if DEBUG == 3460
01322 
01323                 int i_OrderedVertexCount;
01324 
01325                 cout<<endl;
01326                 cout<<"DEBUG 3460 | Vertex Ordering | Smallest Last"<<endl;
01327                 cout<<endl;
01328 
01329                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01330 
01331                 for(i=0; i<i_OrderedVertexCount; i++)
01332                 {
01333                         cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_OrderedVertices[i])<<endl;
01334                 }
01335 
01336                 cout<<endl;
01337                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
01338                 cout<<endl;
01339 
01340 #endif
01341 
01342                 return(_TRUE);
01343         }
01344 
01345 
01346         int BipartiteGraphOrdering::SelectiveIncidenceDegreeOrdering()
01347         {
01348                 if(CheckVertexOrdering("SELECTIVE_INCIDENCE_DEGREE"))
01349                 {
01350                         return(_TRUE);
01351                 }
01352 
01353                 int i, j;
01354 
01355                 int i_HighestDegreeVertex, m_i_MaximumVertexDegree;
01356 
01357                 int i_LeftVertexCount, i_RightVertexCount;
01358 
01359                 int i_IncidenceVertexDegree, i_IncidenceVertexDegreeCount;
01360 
01361                 int i_IncludedVertexCount;
01362 
01363                 int i_SelectedVertex, i_SelectedVertexCount;
01364 
01365                 vector<int> vi_IncidenceVertexDegree;
01366 
01367                 vector< list<int> > vli_GroupedIncidenceVertexDegree;
01368 
01369                 vector< list<int>::iterator > vlit_VertexLocation;
01370 
01371                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01372                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01373 
01374                 vi_IncidenceVertexDegree.clear();
01375                 vi_IncidenceVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
01376 
01377                 vli_GroupedIncidenceVertexDegree.clear();
01378                 vli_GroupedIncidenceVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01379 
01380                 vlit_VertexLocation.clear();
01381                 vlit_VertexLocation.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01382 
01383                 i_SelectedVertex = _UNKNOWN;
01384 
01385                 i_IncludedVertexCount = _FALSE;
01386 
01387                 i_HighestDegreeVertex = m_i_MaximumVertexDegree = _UNKNOWN;
01388 
01389                 for(i=0; i<i_LeftVertexCount; i++)
01390                 {
01391                         if(m_vi_IncludedLeftVertices[i] == _FALSE)
01392                         {
01393                                 continue;
01394                         }
01395 
01396                         i_IncludedVertexCount++;
01397 
01398                         i_IncidenceVertexDegree = _FALSE;
01399 
01400                         vi_IncidenceVertexDegree[i] = i_IncidenceVertexDegree;
01401 
01402                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i);
01403 
01404                         vlit_VertexLocation[vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].front()] = vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin();
01405 
01406                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01407                         {
01408                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01409                                 {
01410                                         continue;
01411                                 }
01412 
01413                                 i_IncidenceVertexDegree++;
01414                         }
01415 
01416                         if(m_i_MaximumVertexDegree < i_IncidenceVertexDegree)
01417                         {
01418                                 m_i_MaximumVertexDegree = i_IncidenceVertexDegree;
01419 
01420                                 i_HighestDegreeVertex = i;
01421                         }
01422                 }
01423 
01424                 for(i=0; i<i_RightVertexCount; i++)
01425                 {
01426                 if(m_vi_IncludedRightVertices[i] == _FALSE)
01427                         {
01428                                 continue;
01429                         }
01430 
01431                         i_IncludedVertexCount++;
01432 
01433                         i_IncidenceVertexDegree = _FALSE;
01434 
01435                         vi_IncidenceVertexDegree[i + i_LeftVertexCount] = i_IncidenceVertexDegree;
01436 
01437                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i + i_LeftVertexCount);
01438 
01439                         vlit_VertexLocation[vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].front()] = vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin();
01440 
01441                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01442                         {
01443                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01444                                 {
01445                                         continue;
01446                                 }
01447 
01448                                 i_IncidenceVertexDegree++;
01449                         }
01450 
01451                         if(m_i_MaximumVertexDegree < i_IncidenceVertexDegree)
01452                         {
01453                                 m_i_MaximumVertexDegree = i_IncidenceVertexDegree;
01454 
01455                                 i_HighestDegreeVertex = i + i_LeftVertexCount;
01456                         }
01457                 }
01458 
01459 #if DEBUG == 3461
01460 
01461                 list<int>::iterator lit_ListIterator;
01462 
01463                 cout<<endl;
01464                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Vertex Degrees"<<endl;
01465                 cout<<endl;
01466 
01467                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01468                 {
01469                         cout<<"Degree "<<i<<"\t"<<" : ";
01470 
01471                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
01472 
01473                         j = _FALSE;
01474 
01475                         for(lit_ListIterator = vli_GroupedIncidenceVertexDegree[i].begin(); lit_ListIterator != vli_GroupedIncidenceVertexDegree[i].end(); lit_ListIterator++)
01476                         {
01477                                 if(j==STEP_DOWN(i_IncidenceVertexDegreeCount))
01478                                 {
01479                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_IncidenceVertexDegreeCount<<")";
01480                                 }
01481                                 else
01482                                 {
01483                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
01484                                 }
01485 
01486                                 j++;
01487                         }
01488 
01489                  cout<<endl;
01490                 }
01491 
01492                 cout<<endl;
01493                 cout<<"[Highest Degree Vertex = "<<STEP_UP(i_HighestDegreeVertex)<<"; Highest Vertex Degree = "<<m_i_MaximumVertexDegree<<"; Candidate Vertex Count = "<<i_IncludedVertexCount<<"]"<<endl;
01494                 cout<<endl;
01495 
01496 #endif
01497 
01498                 m_vi_OrderedVertices.clear();
01499 
01500                 i_SelectedVertexCount = _FALSE;
01501 
01502                 while(i_SelectedVertexCount < i_IncludedVertexCount)
01503                 {
01504                         if(i_SelectedVertexCount == _FALSE)
01505                         {
01506                                 i_SelectedVertex = i_HighestDegreeVertex;
01507                         }
01508                         else
01509                         {
01510                                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01511                                 {
01512                                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
01513 
01514                                         if(i_IncidenceVertexDegreeCount != _FALSE)
01515                                         {
01516                                                 i_SelectedVertex = vli_GroupedIncidenceVertexDegree[i].front();
01517 
01518                                                 break;
01519                                         }
01520                                 }
01521                         }
01522 
01523                         if(i_SelectedVertex < i_LeftVertexCount)
01524                         {
01525 
01526 #if DEBUG == 3461
01527 
01528                                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Selected Left Vertex | "<<STEP_UP(i_SelectedVertex)<<" [Selection "<<STEP_UP(i_SelectedVertexCount)<<"]"<<endl;
01529 
01530 #endif
01531 
01532                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01533                                 {
01534                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01535                                         {
01536                                                 continue;
01537                                         }
01538 
01539                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01540 
01541                                         vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01542 
01543                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01544 
01545                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01546 
01547 #if DEBUG == 3461
01548 
01549                                         cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Repositioned Right Vertex | "<<STEP_UP(m_vi_Edges[i] + i_LeftVertexCount)<<endl;
01550 
01551 #endif
01552 
01553                                 }
01554                         }
01555                         else
01556                         {
01557 
01558 #if DEBUG == 3461
01559 
01560                                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Selected Right Vertex | "<<STEP_UP(i_SelectedVertex)<<" [Selection "<<STEP_UP(i_SelectedVertexCount)<<"]"<<endl;
01561 
01562 #endif
01563 
01564                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01565                                 {
01566                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01567                                         {
01568                                                 continue;
01569                                         }
01570 
01571                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01572 
01573                                         vi_IncidenceVertexDegree[m_vi_Edges[i]] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i]]);
01574 
01575                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01576 
01577                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].begin();
01578 
01579 #if DEBUG == 3461
01580 
01581                                         cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Repositioned Left Vertex | "<<STEP_UP(m_vi_Edges[i])<<endl;
01582 
01583 #endif
01584 
01585                                 }
01586                         }
01587 
01588                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01589 
01590                         vi_IncidenceVertexDegree[i_SelectedVertex] = _UNKNOWN;
01591 
01592                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
01593 
01594                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01595 
01596                 }
01597 
01598 #if DEBUG == 3461
01599 
01600                 int i_OrderedVertexCount;
01601 
01602                 cout<<endl;
01603                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree"<<endl;
01604                 cout<<endl;
01605 
01606                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01607 
01608                 for(i=0; i<i_OrderedVertexCount; i++)
01609                 {
01610                         cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_OrderedVertices[i])<<endl;
01611                 }
01612 
01613                 cout<<endl;
01614                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
01615                 cout<<endl;
01616 
01617 #endif
01618 
01619                 return(_TRUE);
01620         }
01621 
01622 
01623         string BipartiteGraphOrdering::GetVertexOrderingVariant()
01624         {
01625 
01626                 if(m_s_VertexOrderingVariant.compare("NATURAL") == 0)
01627                 {
01628                         return("Natural");
01629                 }
01630                 else
01631                 if(m_s_VertexOrderingVariant.compare("LARGEST_FIRST") == 0)
01632                 {
01633                         return("Largest First");
01634                 }
01635                 else
01636                 if(m_s_VertexOrderingVariant.compare("SMALLEST_LAST") == 0)
01637                 {
01638                         return("Smallest Last");
01639                 }
01640                 else
01641                 if(m_s_VertexOrderingVariant.compare("INCIDENCE_DEGREE") == 0)
01642                 {
01643                         return("Incidence Degree");
01644                 }
01645                 else
01646                 if(m_s_VertexOrderingVariant.compare("SELECTVE_LARGEST_FIRST") == 0)
01647                 {
01648                         return("Selective Largest First");
01649                 }
01650                 else
01651                 if(m_s_VertexOrderingVariant.compare("SELECTVE_SMALLEST_FIRST") == 0)
01652                 {
01653                         return("Selective Smallest Last");
01654                 }
01655                 else
01656                 if(m_s_VertexOrderingVariant.compare("SELECTIVE_INCIDENCE_DEGREE") == 0)
01657                 {
01658                         return("Selective Incidence Degree");
01659                 }
01660                 else
01661                 if(m_s_VertexOrderingVariant.compare("DYNAMIC_LARGEST_FIRST") == 0)
01662                 {
01663                         return("Dynamic Largest First");
01664                 }
01665                 else
01666                 {
01667                         return("Unknown");
01668                 }
01669         }
01670 
01671         void BipartiteGraphOrdering::GetOrderedVertices(vector<int> &output)
01672         {
01673                 output = (m_vi_OrderedVertices);
01674         }
01675 
01676         int BipartiteGraphOrdering::OrderVertices(string s_OrderingVariant) {
01677                 s_OrderingVariant = toUpper(s_OrderingVariant);
01678 
01679                 if((s_OrderingVariant.compare("NATURAL") == 0))
01680                 {
01681                         return(NaturalOrdering());
01682                 }
01683                 else
01684                 if((s_OrderingVariant.compare("LARGEST_FIRST") == 0))
01685                 {
01686                         return(LargestFirstOrdering());
01687                 }
01688                 else
01689                 if((s_OrderingVariant.compare("DYNAMIC_LARGEST_FIRST") == 0))
01690                 {
01691                         return(DynamicLargestFirstOrdering());
01692                 }
01693                 else
01694                 if((s_OrderingVariant.compare("SMALLEST_LAST") == 0))
01695                 {
01696                         return(SmallestLastOrdering());
01697                 }
01698                 else
01699                 if((s_OrderingVariant.compare("INCIDENCE_DEGREE") == 0))
01700                 {
01701                         return(IncidenceDegreeOrdering());
01702                 }
01703                 else
01704                 if((s_OrderingVariant.compare("RANDOM") == 0))
01705                 {
01706                         return(RandomOrdering());
01707                 }
01708                 else
01709                 {
01710                         cerr<<endl;
01711                         cerr<<"Unknown Ordering Method: "<<s_OrderingVariant;
01712                         cerr<<endl;
01713                 }
01714 
01715                 return(_TRUE);
01716         }
01717 
01718         void BipartiteGraphOrdering::PrintVertexOrdering() {
01719                 cout<<"PrintVertexOrdering() "<<m_s_VertexOrderingVariant<<endl;
01720                 for(unsigned int i=0; i<m_vi_OrderedVertices.size();i++) {
01721                         //printf("\t [%d] %d \n", i, m_vi_OrderedVertices[i]);
01722                         cout<<"\t["<<setw(5)<<i<<"] "<<setw(5)<<m_vi_OrderedVertices[i]<<endl;
01723                 }
01724                 cout<<endl;
01725         }
01726         
01727         double BipartiteGraphOrdering::GetVertexOrderingTime() {
01728           return m_d_OrderingTime;
01729         }
01730 
01731 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines