• Main Page
  • Modules
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

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         //Public Function 3455
00088         int BipartiteGraphOrdering::NaturalOrdering()
00089         {
00090                 if(CheckVertexOrdering("NATURAL"))
00091                 {
00092                         return(_TRUE);
00093                 }
00094 
00095                 int i;
00096 
00097                 int i_LeftVertexCount, i_RightVertexCount;
00098 
00099                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00100                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00101 
00102                 m_vi_OrderedVertices.clear();
00103                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00104 
00105                 for(i=0; i<i_LeftVertexCount; i++)
00106                 {
00107                         m_vi_OrderedVertices.push_back(i);
00108                 }
00109 
00110                 for(i=0; i<i_RightVertexCount; i++)
00111                 {
00112                         m_vi_OrderedVertices.push_back(i + i_LeftVertexCount);
00113                 }
00114 
00115                 return(_TRUE);
00116         }
00117 
00118         int BipartiteGraphOrdering::RandomOrdering()
00119         {
00120                 if(CheckVertexOrdering("RANDOM"))
00121                 {
00122                         return(_TRUE);
00123                 }
00124 
00125                 m_s_VertexOrderingVariant = "RANDOM";
00126 
00127                 int i;
00128 
00129                 int i_LeftVertexCount, i_RightVertexCount;
00130 
00131                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00132                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00133 
00134                 m_vi_OrderedVertices.clear();
00135 
00136                 //Order left vertices
00137                 m_vi_OrderedVertices.resize((unsigned) i_LeftVertexCount);
00138 
00139                 for(unsigned int i = 0; i<i_LeftVertexCount; i++) {
00140                         m_vi_OrderedVertices[i] = i;
00141                 }
00142 
00143                 randomOrdering(m_vi_OrderedVertices);
00144 
00145                 //Order right vertices
00146                 vector<int> tempOrdering;
00147 
00148                 tempOrdering.resize((unsigned) i_RightVertexCount);
00149 
00150                 for(unsigned int i = 0; i<i_RightVertexCount; i++) {
00151                         tempOrdering[i] = i + i_LeftVertexCount;
00152                 }
00153 
00154                 randomOrdering(tempOrdering);
00155                 
00156                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00157 
00158                 //Now, populate vector m_vi_OrderedVertices with the right vertices
00159                 for(unsigned int i = 0; i<i_RightVertexCount; i++) {
00160                         m_vi_OrderedVertices.push_back(tempOrdering[i]);
00161                 }
00162                 
00163                 return(_TRUE);
00164         }
00165 
00166         //Public Function 3456
00167         int BipartiteGraphOrdering::LargestFirstOrdering()
00168         {
00169                 if(CheckVertexOrdering("LARGEST_FIRST"))
00170                 {
00171                         return(_TRUE);
00172                 }
00173 
00174                 int i, j;
00175 
00176                 int i_LeftVertexCount, i_RightVertexCount;
00177 
00178                 int i_HighestDegreeVertex;
00179 
00180                 int i_VertexDegree, i_VertexDegreeCount;
00181 
00182                 vector< vector< int > > v2i_GroupedVertexDegree;
00183 
00184                 m_i_MaximumVertexDegree = _FALSE;
00185 
00186                 i_HighestDegreeVertex = _UNKNOWN;
00187 
00188                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00189                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00190 
00191                 v2i_GroupedVertexDegree.clear();
00192                 v2i_GroupedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00193 
00194                 for(i=0; i<i_LeftVertexCount; i++)
00195                 {
00196                         i_VertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00197 
00198                         v2i_GroupedVertexDegree[i_VertexDegree].push_back(i);
00199 
00200                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00201                         {
00202                                 m_i_MaximumVertexDegree = i_VertexDegree;
00203 
00204                                 i_HighestDegreeVertex = i;
00205                         }
00206                 }
00207 
00208                 for(i=0; i<i_RightVertexCount; i++)
00209                 {
00210                         i_VertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00211 
00212                         v2i_GroupedVertexDegree[i_VertexDegree].push_back(i + i_LeftVertexCount);
00213 
00214                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00215                         {
00216                                 m_i_MaximumVertexDegree = i_VertexDegree;
00217 
00218                                 i_HighestDegreeVertex = i + i_LeftVertexCount;
00219                         }
00220                 }
00221 
00222                 m_vi_OrderedVertices.clear();
00223                 m_vi_OrderedVertices.reserve(i_LeftVertexCount + i_RightVertexCount);
00224 
00225                 if(i_HighestDegreeVertex < i_LeftVertexCount)
00226                 {
00227                         for(i=m_i_MaximumVertexDegree; i>=0; i--)
00228                         {
00229                                 i_VertexDegreeCount = (signed) v2i_GroupedVertexDegree[i].size();
00230 
00231                                 for(j=0; j<i_VertexDegreeCount; j++)
00232                                 {
00233                                         m_vi_OrderedVertices.push_back(v2i_GroupedVertexDegree[i][j]);
00234                                 }
00235                         }
00236                 }
00237                 else
00238                 {
00239                         for(i=m_i_MaximumVertexDegree; i>=0; i--)
00240                         {
00241                                 i_VertexDegreeCount = (signed) v2i_GroupedVertexDegree[i].size();
00242 
00243                                 for(j=STEP_DOWN(i_VertexDegreeCount); j>=0; j--)
00244                                 {
00245                                         m_vi_OrderedVertices.push_back(v2i_GroupedVertexDegree[i][j]);
00246                                 }
00247                         }
00248                 }
00249 
00250                 v2i_GroupedVertexDegree.clear();
00251 
00252                 return(_TRUE);
00253         }
00254 
00255         int BipartiteGraphOrdering::SmallestLastOrdering()
00256         {
00257                 if(CheckVertexOrdering("SMALLEST_LAST"))
00258                 {
00259                         return(_TRUE);
00260                 }
00261 
00262                 int i, u, l, v;
00263 
00264                 int _FOUND;
00265 
00266                 int i_HighestInducedVertexDegree, i_HighestInducedDegreeVertex;
00267 
00268                 int i_LeftVertexCount, i_RightVertexCount;
00269 
00270                 int i_VertexCountMinus1; // = i_LeftVertexCount + i_RightVertexCount - 1, used when inserting selected vertices into m_vi_OrderedVertices
00271 
00272                 int i_InducedVertexDegree;
00273 
00274                 int i_InducedVertexDegreeCount;
00275 
00276                 int i_SelectedVertex, i_SelectedVertexCount;
00277 
00278                 vector <int> vi_InducedVertexDegree;
00279 
00280                 vector < vector < int > > vvi_GroupedInducedVertexDegree;
00281 
00282                 vector <  int > vi_VertexLocation;
00283                 
00284                 vector <  int > vi_LeftSidedVertexinBucket;
00285 
00286 
00287                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00288                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00289 
00290                 i_VertexCountMinus1 = i_LeftVertexCount + i_RightVertexCount - 1;
00291 
00292                 vi_InducedVertexDegree.clear();
00293                 vi_InducedVertexDegree.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00294 
00295                 vvi_GroupedInducedVertexDegree.clear();
00296                 vvi_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00297 
00298                 vi_VertexLocation.clear();
00299                 vi_VertexLocation.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00300 
00301                 vi_LeftSidedVertexinBucket.clear();
00302                 vi_LeftSidedVertexinBucket.reserve((unsigned) i_LeftVertexCount + i_RightVertexCount);
00303 
00304                 i_HighestInducedVertexDegree = _FALSE;
00305 
00306                 i_HighestInducedDegreeVertex = _UNKNOWN;
00307 
00308                 i_SelectedVertex = _UNKNOWN;
00309 
00310                 for(i=0; i<i_LeftVertexCount; i++)
00311                 {
00312                         i_InducedVertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00313 
00314                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00315 
00316                         vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].push_back(i);
00317 
00318                         vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1);
00319 
00320                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00321                         {
00322                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00323 
00324                                 i_HighestInducedDegreeVertex = i;
00325                         }
00326                 }
00327 
00328                 
00329                 // get the bucket division positions now
00330                 for(i= 0; i < i_LeftVertexCount + i_RightVertexCount; i++)
00331                         vi_LeftSidedVertexinBucket.push_back(vvi_GroupedInducedVertexDegree[i].size());
00332         
00333 
00334                 for(i=0; i<i_RightVertexCount; i++)
00335                 {
00336                         i_InducedVertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00337 
00338                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00339 
00340                         vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].push_back(i + i_LeftVertexCount);
00341 
00342                         vi_VertexLocation.push_back(vvi_GroupedInducedVertexDegree[i_InducedVertexDegree].size() - 1);
00343 
00344                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00345                         {
00346                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00347 
00348                                 i_HighestInducedDegreeVertex = i + i_LeftVertexCount;
00349                         }
00350                 }
00351 
00352                 m_vi_OrderedVertices.clear();
00353                 m_vi_OrderedVertices.resize(i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
00354 
00355                 i_SelectedVertexCount = _FALSE;
00356 
00357                 int iMin = 1;
00358 
00359                 while(i_SelectedVertexCount < i_LeftVertexCount + i_RightVertexCount)
00360                 {
00361                         if(iMin != 0 && vvi_GroupedInducedVertexDegree[iMin -1].size() != _FALSE)
00362                                 iMin--;
00363 
00364                         for(i=iMin; i<STEP_UP(i_HighestInducedVertexDegree); i++)
00365                         {
00366                                 i_InducedVertexDegreeCount = (signed) vvi_GroupedInducedVertexDegree[i].size();
00367 
00368                                 if(i_InducedVertexDegreeCount == _FALSE)
00369                                 {
00370                                         iMin++;
00371                                         continue;
00372                                 }
00373 
00374                                 if(i_HighestInducedDegreeVertex < i_LeftVertexCount)
00375                                 {
00376                                         _FOUND = _FALSE;
00377 
00378                                         /*
00379                                         if(vi_LeftSidedVertexinBucket[i] > 0)
00380                                         {
00381                                                 vi_LeftSidedVertexinBucket[i]--;
00382                                                 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i][vi_LeftSidedVertexinBucket[i]];
00383 
00384                                                 vvi_GroupedInducedVertexDegree[i][vi_LeftSidedVertexinBucket[i]] = vvi_GroupedInducedVertexDegree[i].back();
00385                                                 vi_VertexLocation[vvi_GroupedInducedVertexDegree[i].back()] = vi_VertexLocation[u];
00386 
00387                                                 _FOUND = _TRUE;                                         
00388                                         }
00389                                         */
00390                                         if(vi_LeftSidedVertexinBucket[i] > 0)
00391                                         for(int j  = 0; j < vvi_GroupedInducedVertexDegree[i].size(); j++)
00392                                         {
00393                                                 u = vvi_GroupedInducedVertexDegree[i][j];
00394                                                 if(u < i_LeftVertexCount)
00395                                                 {
00396                                                         i_SelectedVertex = u;
00397 
00398                                                         if(vvi_GroupedInducedVertexDegree[i].size() > 1)
00399                                                         {
00400                                                                 // swap this node with the last node
00401                                                                 vvi_GroupedInducedVertexDegree[i][j] = vvi_GroupedInducedVertexDegree[i].back();
00402                                                                 vi_VertexLocation[vvi_GroupedInducedVertexDegree[i].back()] = vi_VertexLocation[u];
00403                                                         }
00404                                                         _FOUND = _TRUE;
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                                                 i_SelectedVertex = vvi_GroupedInducedVertexDegree[i].back();
00439                                 }
00440 
00441                                 break;
00442                         }
00443 
00444                         vvi_GroupedInducedVertexDegree[i].pop_back(); // remove the selected vertex from the bucket
00445 
00446                         if(i_SelectedVertex < i_LeftVertexCount)
00447                         {
00448                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
00449                                 {
00450                                         u = m_vi_Edges[i] + i_LeftVertexCount; // neighbour are always right sided
00451 
00452                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00453                                         {
00454                                                 continue;
00455                                         }
00456                                         
00457                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00458                                         if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1)
00459                                         {
00460                                                 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back();
00461 
00462                                                 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l;
00463 
00464                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00465                                         }
00466                                         // remove last element from this bucket
00467                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back();
00468         
00469                                         // reduce degree of u by 1
00470                                         vi_InducedVertexDegree[u]--;
00471 
00472                                         // move u to appropriate bucket
00473                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u);
00474 
00475                                         // update vi_VertexLocation[u] since it has now been changed
00476                                          vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00477                                         
00478                                         /*
00479                                         if(u < i_LeftVertexCount)
00480                                         {
00481                                                 // swap this vertex and location                                                
00482                                                 v = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00483                                                 if(v > 0)
00484                                                 {
00485                                                         l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v];
00486 
00487                                                         swap(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]], vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v]);
00488                                                         swap(vi_VertexLocation[u], vi_VertexLocation[l]);
00489                                                 }
00490                                                 vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]++;
00491                                         }*/
00492                                 }
00493                         }
00494                         else
00495                         {
00496                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
00497                                 {
00498                                         u = m_vi_Edges[i]; // neighbour are always left sided
00499 
00500                                         if(vi_InducedVertexDegree[u] == _UNKNOWN)
00501                                         {
00502                                                 continue;
00503                                         }
00504                                         
00505                                         // move the last element in this bucket to u's position to get rid of expensive erase operation
00506                                         if(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() > 1)
00507                                         {
00508                                                 l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].back();
00509 
00510                                                 vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_VertexLocation[u]] = l;
00511 
00512                                                 vi_VertexLocation[l] = vi_VertexLocation[u];
00513                                         }
00514                                         
00515                                         // remove last element from this bucket
00516                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].pop_back();
00517         
00518                                         // reduce degree of u by 1
00519                                         vi_InducedVertexDegree[u]--;
00520 
00521                                         // move u to appropriate bucket
00522                                         vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].push_back(u);
00523 
00524                                         // update vi_VertexLocation[u] since it has now been changed
00525                                         vi_VertexLocation[u] = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00526 
00527                                         /*
00528                                         if(u < i_LeftVertexCount)
00529                                         {
00530                                                 // swap this vertex and location                                                
00531                                                 v = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]].size() - 1;
00532                                                 if(v > 0)
00533                                                 {
00534                                                         l = vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v];
00535 
00536                                                         swap(vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]], vvi_GroupedInducedVertexDegree[vi_InducedVertexDegree[u]][v]);
00537                                                         swap(vi_VertexLocation[u], vi_VertexLocation[l]);
00538                                                 }
00539                                                 vi_LeftSidedVertexinBucket[vi_InducedVertexDegree[u]]++;
00540                                         }*/
00541 
00542                                 }
00543                         }
00544 
00545                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
00546 
00547                         m_vi_OrderedVertices[i_VertexCountMinus1 - i_SelectedVertexCount] = i_SelectedVertex;
00548 
00549                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
00550                 }
00551 
00552                 vi_InducedVertexDegree.clear();
00553                 vvi_GroupedInducedVertexDegree.clear();
00554                 vi_VertexLocation.clear();
00555                 vi_LeftSidedVertexinBucket.clear();
00556 
00557                 return(_TRUE);
00558         }
00559 
00560 /*
00561         //Public Function 3457
00562         int BipartiteGraphOrdering::SmallestLastOrdering()
00563         {
00564                 if(CheckVertexOrdering("SMALLEST_LAST"))
00565                 {
00566                         return(_TRUE);
00567                 }
00568 
00569                 int i;
00570 
00571                 int _FOUND;
00572 
00573                 int i_HighestInducedVertexDegree, i_HighestInducedDegreeVertex;
00574 
00575                 int i_LeftVertexCount, i_RightVertexCount;
00576 
00577                 int i_VertexCountMinus1; // = i_LeftVertexCount + i_RightVertexCount - 1, used when inserting selected vertices into m_vi_OrderedVertices
00578 
00579                 int i_InducedVertexDegree;
00580 
00581                 int i_InducedVertexDegreeCount;
00582 
00583                 int i_SelectedVertex, i_SelectedVertexCount;
00584 
00585                 vector<int> vi_InducedVertexDegree;
00586 
00587                 vector< list<int> > vli_GroupedInducedVertexDegree;
00588 
00589                 vector< list<int>::iterator > vlit_VertexLocation;
00590 
00591                 list<int>::iterator lit_ListIterator;
00592 
00593                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00594                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00595 
00596                 i_VertexCountMinus1 = i_LeftVertexCount + i_RightVertexCount - 1;
00597 
00598                 vi_InducedVertexDegree.clear();
00599 
00600                 vli_GroupedInducedVertexDegree.clear();
00601                 vli_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00602 
00603                 vlit_VertexLocation.clear();
00604 
00605                 i_HighestInducedVertexDegree = _FALSE;
00606 
00607                 i_HighestInducedDegreeVertex = _UNKNOWN;
00608 
00609                 i_SelectedVertex = _UNKNOWN;
00610 
00611                 for(i=0; i<i_LeftVertexCount; i++)
00612                 {
00613                         i_InducedVertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00614 
00615                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00616 
00617                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i);
00618 
00619                         vlit_VertexLocation.push_back(vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin());
00620 
00621                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00622                         {
00623                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00624 
00625                                 i_HighestInducedDegreeVertex = i;
00626                         }
00627                 }
00628 
00629                 for(i=0; i<i_RightVertexCount; i++)
00630                 {
00631                         i_InducedVertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00632 
00633                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
00634 
00635                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i + i_LeftVertexCount);
00636 
00637                         vlit_VertexLocation.push_back(vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin());
00638 
00639                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
00640                         {
00641                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
00642 
00643                                 i_HighestInducedDegreeVertex = i + i_LeftVertexCount;
00644                         }
00645                 }
00646 
00647 
00648 #if DEBUG == 3457
00649 
00650                 int j;
00651 
00652                 cout<<endl;
00653                 cout<<"DEBUG 3457 | Bipartite Graph Coloring | Bipartite Smallest Last Ordering"<<endl;
00654                 cout<<endl;
00655 
00656                 for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
00657                 {
00658                         cout<<"Degree "<<i<<"\t"<<" : ";
00659 
00660                         i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
00661 
00662                         j = _FALSE;
00663 
00664                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
00665                         {
00666                                 if(j==STEP_DOWN(i_InducedVertexDegreeCount))
00667                                 {
00668                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_InducedVertexDegreeCount<<")";
00669                                 }
00670                                 else
00671                                 {
00672                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
00673                                 }
00674 
00675                                 j++;
00676                         }
00677 
00678                         cout<<endl;
00679                 }
00680 
00681                 cout<<endl;
00682                 cout<<"[Highest Degree Vertex = "<<STEP_UP(i_HighestInducedDegreeVertex)<<"; Highest Induced Vertex Degree = "<<i_HighestInducedVertexDegree<<"]"<<endl;
00683                 cout<<endl;
00684 
00685 #endif
00686 
00687                 m_vi_OrderedVertices.clear();
00688                 m_vi_OrderedVertices.resize(i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
00689 
00690                 i_SelectedVertexCount = _FALSE;
00691 
00692                 while(i_SelectedVertexCount < i_LeftVertexCount + i_RightVertexCount)
00693                 {
00694                         for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
00695                         {
00696                                 i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
00697 
00698                                 if(i_InducedVertexDegreeCount == _FALSE)
00699                                 {
00700                                         continue;
00701                                 }
00702 
00703                                 if(i_HighestInducedDegreeVertex < i_LeftVertexCount)
00704                                 {
00705                                         _FOUND = _FALSE;
00706 
00707                                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
00708                                         {
00709                                                 if(*lit_ListIterator < i_LeftVertexCount)
00710                                                 {
00711                                                         i_SelectedVertex = *lit_ListIterator;
00712 
00713                                                         _FOUND = _TRUE;
00714 
00715                                                         break;
00716                                                 }
00717                                         }
00718 
00719                                         if(!_FOUND)
00720                                         {
00721                                                 i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
00722                                         }
00723 
00724                                         break;
00725                                 }
00726                                 else
00727                                 {
00728                                         _FOUND = _FALSE;
00729 
00730                                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
00731                                         {
00732                                                 if(*lit_ListIterator >= i_LeftVertexCount)
00733                                                 {
00734                                                         i_SelectedVertex = *lit_ListIterator;
00735 
00736                                                         _FOUND = _TRUE;
00737 
00738                                                         break;
00739                                                 }
00740                                         }
00741 
00742                                         if(!_FOUND)
00743                                         {
00744                                                 i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
00745                                         }
00746                                 }
00747 
00748                                 break;
00749                         }
00750 
00751                         if(i_SelectedVertex < i_LeftVertexCount)
00752                         {
00753                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
00754                                 {
00755                                         if(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
00756                                         {
00757                                                 continue;
00758                                         }
00759 
00760                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
00761 
00762                                         vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
00763 
00764                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
00765 
00766                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
00767                                 }
00768                         }
00769                         else
00770                         {
00771                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
00772                                 {
00773                                         if(vi_InducedVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
00774                                         {
00775                                                 continue;
00776                                         }
00777 
00778                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
00779 
00780                                         vi_InducedVertexDegree[m_vi_Edges[i]] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i]]);
00781 
00782                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
00783 
00784                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].begin();
00785                                 }
00786                         }
00787 
00788                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
00789 
00790                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
00791 
00792                         //insert at the beginning (push_front)? BAD!!!! due to moving
00793                         //m_vi_OrderedVertices.insert(m_vi_OrderedVertices.begin(), i_SelectedVertex);
00794 
00795                         //Alternative method for the above (commented) line
00796 //printf("m_vi_OrderedVertices[i_VertexCountMinus1 (%d) - i_SelectedVertexCount (%d)] (%d) = i_SelectedVertex (%d) \n", i_VertexCountMinus1, i_SelectedVertexCount, i_VertexCountMinus1 - i_SelectedVertexCount, i_SelectedVertex);
00797                         m_vi_OrderedVertices[i_VertexCountMinus1 - i_SelectedVertexCount] = i_SelectedVertex;
00798 
00799                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
00800                 }
00801 
00802 #if DEBUG == 3457
00803 
00804                 int i_OrderedVertexCount;
00805 
00806                 cout<<endl;
00807                 cout<<"DEBUG 3457 | Bipartite Graph Coloring | Bipartite Smallest Last Ordering"<<endl;
00808                 cout<<endl;
00809 
00810                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
00811 
00812                 for(i=0; i<i_OrderedVertexCount; i++)
00813                 {
00814                         if(i == STEP_DOWN(i_OrderedVertexCount))
00815                         {
00816                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_OrderedVertexCount<<")"<<endl;
00817                         }
00818                         else
00819                         {
00820                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
00821                         }
00822                 }
00823 
00824                 cout<<endl;
00825                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
00826                 cout<<endl;
00827 
00828 #endif
00829 
00830                 return(_TRUE);
00831         }
00832         */
00833 
00834         //Public Function 3458
00835         int BipartiteGraphOrdering::IncidenceDegreeOrdering()
00836         {
00837                 if(CheckVertexOrdering("INCIDENCE_DEGREE"))
00838                 {
00839                         return(_TRUE);
00840                 }
00841 
00842                 int i;
00843 
00844                 int _FOUND;
00845 
00846                 int i_HighestIncidenceDegreeVertex, i_HighestIncidenceVertexDegree;
00847 
00848                 int i_LeftVertexCount, i_RightVertexCount;
00849 
00850                 int i_VertexDegree;
00851 
00852                 int i_IncidenceVertexDegree, i_IncidenceVertexDegreeCount;
00853 
00854                 int i_SelectedVertex, i_SelectedVertexCount;
00855 
00856                 vector<int> vi_IncidenceVertexDegree;
00857 
00858                 vector< list<int> > vli_GroupedIncidenceVertexDegree;
00859 
00860                 vector< list<int>::iterator > vlit_VertexLocation;
00861 
00862                 list<int>::iterator lit_ListIterator;
00863 
00864                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
00865                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
00866 
00867                 vi_IncidenceVertexDegree.clear();
00868 
00869                 vli_GroupedIncidenceVertexDegree.clear();
00870                 vli_GroupedIncidenceVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
00871 
00872                 vlit_VertexLocation.clear();
00873 
00874                 i_HighestIncidenceDegreeVertex = i_HighestIncidenceVertexDegree = _UNKNOWN;
00875 
00876                 i_IncidenceVertexDegree = _FALSE;
00877 
00878                 i_SelectedVertex = _UNKNOWN;
00879 
00880                 for(i=0; i<i_LeftVertexCount; i++)
00881                 {
00882                         vi_IncidenceVertexDegree.push_back(i_IncidenceVertexDegree);
00883 
00884                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i);
00885 
00886                         vlit_VertexLocation.push_back(vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin());
00887 
00888                         i_VertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
00889 
00890                         if(i_HighestIncidenceVertexDegree < i_VertexDegree)
00891                         {
00892                                 i_HighestIncidenceVertexDegree = i_VertexDegree;
00893 
00894                                 i_HighestIncidenceDegreeVertex = i;
00895                         }
00896                 }
00897 
00898                 for(i=0; i<i_RightVertexCount; i++)
00899                 {
00900                         vi_IncidenceVertexDegree.push_back(i_IncidenceVertexDegree);
00901 
00902                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i + i_LeftVertexCount);
00903 
00904                         vlit_VertexLocation.push_back(vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin());
00905 
00906                         i_VertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
00907 
00908                         if(i_HighestIncidenceVertexDegree < i_VertexDegree)
00909                         {
00910                                 i_HighestIncidenceVertexDegree = i_VertexDegree;
00911 
00912                                 i_HighestIncidenceDegreeVertex = i + i_LeftVertexCount;
00913                         }
00914                 }
00915 
00916 #if DEBUG == 3458
00917 
00918                 int j;
00919 
00920                 cout<<endl;
00921                 cout<<"DEBUG 3458 | Bipartite Graph Coloring | Bipartite Incidence Degree Ordering"<<endl;
00922                 cout<<endl;
00923 
00924                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
00925                 {
00926                         cout<<"Degree "<<i<<"\t"<<" : ";
00927 
00928                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
00929 
00930                         j = _FALSE;
00931 
00932                         for(lit_ListIterator = vli_GroupedIncidenceVertexDegree[i].begin(); lit_ListIterator != vli_GroupedIncidenceVertexDegree[i].end(); lit_ListIterator++)
00933                         {
00934                                 if(j==STEP_DOWN(i_IncidenceVertexDegreeCount))
00935                                 {
00936                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_IncidenceVertexDegreeCount<<")";
00937                                 }
00938                                 else
00939                                 {
00940                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
00941                                 }
00942 
00943                                 j++;
00944                         }
00945 
00946                         cout<<endl;
00947                 }
00948 
00949                 cout<<endl;
00950                 cout<<"[Highest Degree Vertex = "<<STEP_UP(i_HighestIncidenceDegreeVertex)<<"; Highest Vertex Degree = "<<i_HighestIncidenceVertexDegree<<"]"<<endl;
00951                 cout<<endl;
00952 
00953 #endif
00954 
00955                 m_vi_OrderedVertices.clear();
00956 
00957                 i_SelectedVertexCount = _FALSE;
00958 
00959                 while(i_SelectedVertexCount < i_LeftVertexCount + i_RightVertexCount)
00960                 {
00961                         if(i_SelectedVertexCount == _FALSE)
00962                         {
00963                                 i_SelectedVertex = i_HighestIncidenceDegreeVertex;
00964                         }
00965                         else
00966                         {
00967                                 for(i=i_HighestIncidenceVertexDegree; i>=0; i--)
00968                                 {
00969                                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
00970 
00971                                         if(i_IncidenceVertexDegreeCount == _FALSE)
00972                                         {
00973                                                 continue;
00974                                         }
00975 
00976                                         if(i_HighestIncidenceDegreeVertex < i_LeftVertexCount)
00977                                         {
00978                                                 _FOUND = _FALSE;
00979 
00980                                                 for(lit_ListIterator = vli_GroupedIncidenceVertexDegree[i].begin(); lit_ListIterator != vli_GroupedIncidenceVertexDegree[i].end(); lit_ListIterator++)
00981                                                 {
00982                                                         if(*lit_ListIterator < i_LeftVertexCount)
00983                                                         {
00984                                                                 i_SelectedVertex = *lit_ListIterator;
00985 
00986                                                                 _FOUND = _TRUE;
00987 
00988                                                                 break;
00989                                                         }
00990                                                 }
00991 
00992                                                 if(!_FOUND)
00993                                                 {
00994                                                         i_SelectedVertex = vli_GroupedIncidenceVertexDegree[i].front();
00995                                                 }
00996 
00997                                                 break;
00998                                         }
00999                                         else
01000                                         {
01001                                                 _FOUND = _FALSE;
01002 
01003                                                 for(lit_ListIterator = vli_GroupedIncidenceVertexDegree[i].begin(); lit_ListIterator != vli_GroupedIncidenceVertexDegree[i].end(); lit_ListIterator++)
01004                                                 {
01005                                                         if(*lit_ListIterator >= i_LeftVertexCount)
01006                                                         {
01007                                                                 i_SelectedVertex = *lit_ListIterator;
01008 
01009                                                                 _FOUND = _TRUE;
01010 
01011                                                                 break;
01012                                                         }
01013                                                 }
01014 
01015                                                 if(!_FOUND)
01016                                                 {
01017                                                         i_SelectedVertex = vli_GroupedIncidenceVertexDegree[i].front();
01018                                                 }
01019                                         }
01020 
01021                                         break;
01022                                 }
01023                         }
01024 
01025                         if(i_SelectedVertex < i_LeftVertexCount)
01026                         {
01027                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01028                                 {
01029                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01030                                         {
01031                                                 continue;
01032                                         }
01033 
01034                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01035 
01036                                         vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01037 
01038                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01039 
01040                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01041 
01042                                 }
01043                         }
01044                         else
01045                         {
01046                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01047                                 {
01048                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01049                                         {
01050                                                 continue;
01051                                         }
01052 
01053                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01054 
01055                                         vi_IncidenceVertexDegree[m_vi_Edges[i]] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i]]);
01056 
01057                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01058 
01059                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].begin();
01060 
01061                                 }
01062 
01063                         }
01064 
01065                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01066 
01067                         vi_IncidenceVertexDegree[i_SelectedVertex] = _UNKNOWN;
01068 
01069                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
01070 
01071                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01072                 }
01073 
01074 #if DEBUG == 3458
01075 
01076                 int i_OrderedVertexCount;
01077 
01078                 cout<<endl;
01079                 cout<<"DEBUG 3458 | Bipartite Graph Coloring | Bipartite Incidence Degree Ordering"<<endl;
01080                 cout<<endl;
01081 
01082                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01083 
01084                 for(i=0; i<i_OrderedVertexCount; i++)
01085                 {
01086                         if(i == STEP_DOWN(i_OrderedVertexCount))
01087                         {
01088                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_OrderedVertexCount<<")"<<endl;
01089                         }
01090                         else
01091                         {
01092                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
01093                         }
01094                 }
01095 
01096                 cout<<endl;
01097                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
01098                 cout<<endl;
01099 
01100 #endif
01101 
01102                 return(_TRUE);
01103         }
01104 
01105         //Public Function 3459
01106         int BipartiteGraphOrdering::SelectiveLargestFirstOrdering()
01107         {
01108                 if(CheckVertexOrdering("SELECTVE_LARGEST_FIRST"))
01109                 {
01110                         return(_TRUE);
01111                 }
01112 
01113                 int i, j;
01114 
01115                 int i_LeftVertexCount, i_RightVertexCount;
01116 
01117                 int i_VertexDegree, i_VertexDegreeCount;
01118 
01119                 vector< vector<int> > v2i_GroupedVertexDegree;
01120 
01121                 m_i_MaximumVertexDegree = _FALSE;
01122 
01123                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01124                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01125 
01126                 v2i_GroupedVertexDegree.clear();
01127                 v2i_GroupedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01128 
01129                 for(i=0; i<i_LeftVertexCount; i++)
01130                 {
01131                         if(m_vi_IncludedLeftVertices[i] == _FALSE)
01132                         {
01133                                 continue;
01134                         }
01135 
01136                         i_VertexDegree = _FALSE;
01137 
01138                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01139                         {
01140                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01141                                 {
01142                                         continue;
01143                                 }
01144 
01145                                 i_VertexDegree++;
01146                         }
01147 
01148                         v2i_GroupedVertexDegree[i_VertexDegree].push_back(i);
01149 
01150                         if(m_i_MaximumVertexDegree < i_VertexDegree)
01151                         {
01152                                 m_i_MaximumVertexDegree = i_VertexDegree;
01153                         }
01154                 }
01155 
01156                 for(i=0; i<i_RightVertexCount; i++)
01157                 {
01158                         if(m_vi_IncludedRightVertices[i] == _FALSE)
01159                         {
01160                                 continue;
01161                         }
01162 
01163                         i_VertexDegree = _FALSE;
01164 
01165                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01166                         {
01167                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01168                                 {
01169                                         continue;
01170                                 }
01171 
01172                                 i_VertexDegree++;
01173                         }
01174 
01175                         v2i_GroupedVertexDegree[i_VertexDegree].push_back(i + i_LeftVertexCount);
01176 
01177                         if(m_i_MaximumVertexDegree < i_VertexDegree)
01178                         {
01179                                 m_i_MaximumVertexDegree = i_VertexDegree;
01180                         }
01181                 }
01182 
01183                 m_vi_OrderedVertices.clear();
01184 
01185                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01186                 {
01187                         i_VertexDegreeCount = (signed) v2i_GroupedVertexDegree[i].size();
01188 
01189                         for(j=0; j<i_VertexDegreeCount; j++)
01190                         {
01191                                 m_vi_OrderedVertices.push_back(v2i_GroupedVertexDegree[i][j]);
01192                         }
01193                 }
01194 
01195 #if DEBUG == 3459
01196 
01197                 int i_VertexCount;
01198 
01199                 cout<<endl;
01200                 cout<<"DEBUG 3459 | Bipartite Graph Bicoloring | Largest First Ordering"<<endl;
01201                 cout<<endl;
01202 
01203                 i_VertexCount = (signed) m_vi_OrderedVertices.size();
01204 
01205                 for(i=0; i<i_VertexCount; i++)
01206                 {
01207                         if(i == STEP_DOWN(i_VertexCount))
01208                         {
01209                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_VertexCount<<")"<<endl;
01210                         }
01211                         else
01212                         {
01213                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
01214                         }
01215                 }
01216 
01217                 cout<<endl;
01218                 cout<<"[Highest Vertex Degree = "<<m_i_MaximumVertexDegree<<"]"<<endl;
01219                 cout<<endl;
01220 
01221 #endif
01222 
01223                 return(_TRUE);
01224         }
01225 
01226         //Public Function 3460
01227         int BipartiteGraphOrdering::SelectiveSmallestLastOrdering()
01228         {
01229                 if(CheckVertexOrdering("SELECTIVE_SMALLEST_LAST"))
01230                 {
01231                         return(_TRUE);
01232                 }
01233 
01234                 int i, j;
01235 
01236                 int i_HighestInducedVertexDegree;
01237 
01238                 int i_LeftVertexCount, i_RightVertexCount;
01239 
01240                 int i_InducedVertexDegree;
01241 
01242                 int i_InducedVertexDegreeCount;
01243 
01244                 int i_IncludedVertexCount;
01245 
01246                 int i_SelectedVertex, i_SelectedVertexCount;
01247 
01248                 vector<int> vi_InducedVertexDegree;
01249 
01250                 vector< list<int> > vli_GroupedInducedVertexDegree;
01251 
01252                 vector< list<int>::iterator > vlit_VertexLocation;
01253 
01254                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01255                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01256 
01257                 vi_InducedVertexDegree.clear();
01258                 vi_InducedVertexDegree.resize((signed) i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
01259 
01260                 vli_GroupedInducedVertexDegree.clear();
01261                 vli_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01262 
01263                 vlit_VertexLocation.clear();
01264                 vlit_VertexLocation.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01265 
01266                 i_IncludedVertexCount = _FALSE;
01267 
01268                 i_HighestInducedVertexDegree = _FALSE;
01269 
01270                 i_SelectedVertex = _UNKNOWN;
01271 
01272                 for(i=0; i<i_LeftVertexCount; i++)
01273                 {
01274                 if(m_vi_IncludedLeftVertices[i] == _FALSE)
01275                         {
01276                                 continue;
01277                         }
01278 
01279                         i_IncludedVertexCount++;
01280 
01281                         i_InducedVertexDegree = _FALSE;
01282 
01283                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01284                         {
01285                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01286                                 {
01287                                         continue;
01288                                 }
01289 
01290                                 i_InducedVertexDegree++;
01291                         }
01292 
01293                         vi_InducedVertexDegree[i] = i_InducedVertexDegree;
01294 
01295                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i);
01296 
01297                         vlit_VertexLocation[vli_GroupedInducedVertexDegree[i_InducedVertexDegree].front()] = vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin();
01298 
01299                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01300                         {
01301                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01302                         }
01303                 }
01304 
01305                 for(i=0; i<i_RightVertexCount; i++)
01306                 {
01307                 if(m_vi_IncludedRightVertices[i] == _FALSE)
01308                         {
01309                                 continue;
01310                         }
01311 
01312                         i_IncludedVertexCount++;
01313 
01314                         i_InducedVertexDegree = _FALSE;
01315 
01316                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01317                         {
01318                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01319                                 {
01320                                         continue;
01321                                 }
01322 
01323                                 i_InducedVertexDegree++;
01324                         }
01325 
01326                         vi_InducedVertexDegree[i + i_LeftVertexCount] = i_InducedVertexDegree;
01327 
01328                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i + i_LeftVertexCount);
01329 
01330                         vlit_VertexLocation[vli_GroupedInducedVertexDegree[i_InducedVertexDegree].front()] = vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin();
01331 
01332                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01333                         {
01334                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01335                         }
01336                 }
01337 
01338 
01339 #if DEBUG == 3460
01340 
01341                 list<int>::iterator lit_ListIterator;
01342 
01343                 cout<<endl;
01344                 cout<<"DEBUG 3460 | Vertex Ordering | Vertex Degree"<<endl;
01345                 cout<<endl;
01346 
01347                 for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
01348                 {
01349                         cout<<"Degree "<<i<<"\t"<<" : ";
01350 
01351                         i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01352 
01353                         j = _FALSE;
01354 
01355                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
01356                         {
01357                                 if(j==STEP_DOWN(i_InducedVertexDegreeCount))
01358                                 {
01359                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_InducedVertexDegreeCount<<")";
01360                                 }
01361                                 else
01362                                 {
01363                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
01364                                 }
01365 
01366                                 j++;
01367                         }
01368 
01369                         cout<<endl;
01370                 }
01371 
01372                 cout<<endl;
01373 
01374 #endif
01375 
01376                 m_vi_OrderedVertices.clear();
01377 
01378                 i_SelectedVertexCount = _FALSE;
01379 
01380                 while(i_SelectedVertexCount < i_IncludedVertexCount)
01381                 {
01382                         for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
01383                         {
01384                                 i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01385 
01386                                 if(i_InducedVertexDegreeCount != _FALSE)
01387                                 {
01388                                         i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
01389 
01390                                         break;
01391                                 }
01392                         }
01393 
01394                         if(i_SelectedVertex < i_LeftVertexCount)
01395                         {
01396                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01397                                 {
01398                                         if(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01399                                         {
01400                                                 continue;
01401                                         }
01402 
01403                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01404 
01405                                         vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01406 
01407                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01408 
01409                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01410                                 }
01411                         }
01412                         else
01413                         {
01414                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01415                                 {
01416                                         if(vi_InducedVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01417                                         {
01418                                                 continue;
01419                                         }
01420 
01421                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01422 
01423                                         vi_InducedVertexDegree[m_vi_Edges[i]] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i]]);
01424 
01425                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01426 
01427                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].begin();
01428                                 }
01429                         }
01430 
01431                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01432 
01433                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
01434 
01435                         m_vi_OrderedVertices.insert(m_vi_OrderedVertices.begin(), i_SelectedVertex);
01436 
01437                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01438                 }
01439 
01440 
01441 #if DEBUG == 3460
01442 
01443                 int i_OrderedVertexCount;
01444 
01445                 cout<<endl;
01446                 cout<<"DEBUG 3460 | Vertex Ordering | Smallest Last"<<endl;
01447                 cout<<endl;
01448 
01449                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01450 
01451                 for(i=0; i<i_OrderedVertexCount; i++)
01452                 {
01453                         cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_OrderedVertices[i])<<endl;
01454                 }
01455 
01456                 cout<<endl;
01457                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
01458                 cout<<endl;
01459 
01460 #endif
01461 
01462                 return(_TRUE);
01463         }
01464 
01465 
01466         //Public Function 3461
01467         int BipartiteGraphOrdering::SelectiveIncidenceDegreeOrdering()
01468         {
01469                 if(CheckVertexOrdering("SELECTIVE_INCIDENCE_DEGREE"))
01470                 {
01471                         return(_TRUE);
01472                 }
01473 
01474                 int i, j;
01475 
01476                 int i_HighestDegreeVertex, m_i_MaximumVertexDegree;
01477 
01478                 int i_LeftVertexCount, i_RightVertexCount;
01479 
01480                 int i_IncidenceVertexDegree, i_IncidenceVertexDegreeCount;
01481 
01482                 int i_IncludedVertexCount;
01483 
01484                 int i_SelectedVertex, i_SelectedVertexCount;
01485 
01486                 vector<int> vi_IncidenceVertexDegree;
01487 
01488                 vector< list<int> > vli_GroupedIncidenceVertexDegree;
01489 
01490                 vector< list<int>::iterator > vlit_VertexLocation;
01491 
01492                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01493                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01494 
01495                 vi_IncidenceVertexDegree.clear();
01496                 vi_IncidenceVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount, _UNKNOWN);
01497 
01498                 vli_GroupedIncidenceVertexDegree.clear();
01499                 vli_GroupedIncidenceVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01500 
01501                 vlit_VertexLocation.clear();
01502                 vlit_VertexLocation.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01503 
01504                 i_SelectedVertex = _UNKNOWN;
01505 
01506                 i_IncludedVertexCount = _FALSE;
01507 
01508                 i_HighestDegreeVertex = m_i_MaximumVertexDegree = _UNKNOWN;
01509 
01510                 for(i=0; i<i_LeftVertexCount; i++)
01511                 {
01512                         if(m_vi_IncludedLeftVertices[i] == _FALSE)
01513                         {
01514                                 continue;
01515                         }
01516 
01517                         i_IncludedVertexCount++;
01518 
01519                         i_IncidenceVertexDegree = _FALSE;
01520 
01521                         vi_IncidenceVertexDegree[i] = i_IncidenceVertexDegree;
01522 
01523                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i);
01524 
01525                         vlit_VertexLocation[vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].front()] = vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin();
01526 
01527                         for(j=m_vi_LeftVertices[i]; j<m_vi_LeftVertices[STEP_UP(i)]; j++)
01528                         {
01529                                 if(m_vi_IncludedRightVertices[m_vi_Edges[j]] == _FALSE)
01530                                 {
01531                                         continue;
01532                                 }
01533 
01534                                 i_IncidenceVertexDegree++;
01535                         }
01536 
01537                         if(m_i_MaximumVertexDegree < i_IncidenceVertexDegree)
01538                         {
01539                                 m_i_MaximumVertexDegree = i_IncidenceVertexDegree;
01540 
01541                                 i_HighestDegreeVertex = i;
01542                         }
01543                 }
01544 
01545                 for(i=0; i<i_RightVertexCount; i++)
01546                 {
01547                 if(m_vi_IncludedRightVertices[i] == _FALSE)
01548                         {
01549                                 continue;
01550                         }
01551 
01552                         i_IncludedVertexCount++;
01553 
01554                         i_IncidenceVertexDegree = _FALSE;
01555 
01556                         vi_IncidenceVertexDegree[i + i_LeftVertexCount] = i_IncidenceVertexDegree;
01557 
01558                         vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].push_front(i + i_LeftVertexCount);
01559 
01560                         vlit_VertexLocation[vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].front()] = vli_GroupedIncidenceVertexDegree[i_IncidenceVertexDegree].begin();
01561 
01562                         for(j=m_vi_RightVertices[i]; j<m_vi_RightVertices[STEP_UP(i)]; j++)
01563                         {
01564                                 if(m_vi_IncludedLeftVertices[m_vi_Edges[j]] == _FALSE)
01565                                 {
01566                                         continue;
01567                                 }
01568 
01569                                 i_IncidenceVertexDegree++;
01570                         }
01571 
01572                         if(m_i_MaximumVertexDegree < i_IncidenceVertexDegree)
01573                         {
01574                                 m_i_MaximumVertexDegree = i_IncidenceVertexDegree;
01575 
01576                                 i_HighestDegreeVertex = i + i_LeftVertexCount;
01577                         }
01578                 }
01579 
01580 #if DEBUG == 3461
01581 
01582                 list<int>::iterator lit_ListIterator;
01583 
01584                 cout<<endl;
01585                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Vertex Degrees"<<endl;
01586                 cout<<endl;
01587 
01588                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01589                 {
01590                         cout<<"Degree "<<i<<"\t"<<" : ";
01591 
01592                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
01593 
01594                         j = _FALSE;
01595 
01596                         for(lit_ListIterator = vli_GroupedIncidenceVertexDegree[i].begin(); lit_ListIterator != vli_GroupedIncidenceVertexDegree[i].end(); lit_ListIterator++)
01597                         {
01598                                 if(j==STEP_DOWN(i_IncidenceVertexDegreeCount))
01599                                 {
01600                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_IncidenceVertexDegreeCount<<")";
01601                                 }
01602                                 else
01603                                 {
01604                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
01605                                 }
01606 
01607                                 j++;
01608                         }
01609 
01610                  cout<<endl;
01611                 }
01612 
01613                 cout<<endl;
01614                 cout<<"[Highest Degree Vertex = "<<STEP_UP(i_HighestDegreeVertex)<<"; Highest Vertex Degree = "<<m_i_MaximumVertexDegree<<"; Candidate Vertex Count = "<<i_IncludedVertexCount<<"]"<<endl;
01615                 cout<<endl;
01616 
01617 #endif
01618 
01619                 m_vi_OrderedVertices.clear();
01620 
01621                 i_SelectedVertexCount = _FALSE;
01622 
01623                 while(i_SelectedVertexCount < i_IncludedVertexCount)
01624                 {
01625                         if(i_SelectedVertexCount == _FALSE)
01626                         {
01627                                 i_SelectedVertex = i_HighestDegreeVertex;
01628                         }
01629                         else
01630                         {
01631                                 for(i=m_i_MaximumVertexDegree; i>=0; i--)
01632                                 {
01633                                         i_IncidenceVertexDegreeCount = (signed) vli_GroupedIncidenceVertexDegree[i].size();
01634 
01635                                         if(i_IncidenceVertexDegreeCount != _FALSE)
01636                                         {
01637                                                 i_SelectedVertex = vli_GroupedIncidenceVertexDegree[i].front();
01638 
01639                                                 break;
01640                                         }
01641                                 }
01642                         }
01643 
01644                         if(i_SelectedVertex < i_LeftVertexCount)
01645                         {
01646 
01647 #if DEBUG == 3461
01648 
01649                                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Selected Left Vertex | "<<STEP_UP(i_SelectedVertex)<<" [Selection "<<STEP_UP(i_SelectedVertexCount)<<"]"<<endl;
01650 
01651 #endif
01652 
01653                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01654                                 {
01655                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01656                                         {
01657                                                 continue;
01658                                         }
01659 
01660                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01661 
01662                                         vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01663 
01664                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01665 
01666                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01667 
01668 #if DEBUG == 3461
01669 
01670                                         cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Repositioned Right Vertex | "<<STEP_UP(m_vi_Edges[i] + i_LeftVertexCount)<<endl;
01671 
01672 #endif
01673 
01674                                 }
01675                         }
01676                         else
01677                         {
01678 
01679 #if DEBUG == 3461
01680 
01681                                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Selected Right Vertex | "<<STEP_UP(i_SelectedVertex)<<" [Selection "<<STEP_UP(i_SelectedVertexCount)<<"]"<<endl;
01682 
01683 #endif
01684 
01685                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01686                                 {
01687                                         if(vi_IncidenceVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01688                                         {
01689                                                 continue;
01690                                         }
01691 
01692                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01693 
01694                                         vi_IncidenceVertexDegree[m_vi_Edges[i]] = STEP_UP(vi_IncidenceVertexDegree[m_vi_Edges[i]]);
01695 
01696                                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01697 
01698                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[m_vi_Edges[i]]].begin();
01699 
01700 #if DEBUG == 3461
01701 
01702                                         cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree | Repositioned Left Vertex | "<<STEP_UP(m_vi_Edges[i])<<endl;
01703 
01704 #endif
01705 
01706                                 }
01707                         }
01708 
01709                         vli_GroupedIncidenceVertexDegree[vi_IncidenceVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01710 
01711                         vi_IncidenceVertexDegree[i_SelectedVertex] = _UNKNOWN;
01712 
01713                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
01714 
01715                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01716 
01717                 }
01718 
01719 #if DEBUG == 3461
01720 
01721                 int i_OrderedVertexCount;
01722 
01723                 cout<<endl;
01724                 cout<<"DEBUG 3461 | Vertex Ordering | Incidence Degree"<<endl;
01725                 cout<<endl;
01726 
01727                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01728 
01729                 for(i=0; i<i_OrderedVertexCount; i++)
01730                 {
01731                         cout<<STEP_UP(i)<<"\t"<<" : "<<STEP_UP(m_vi_OrderedVertices[i])<<endl;
01732                 }
01733 
01734                 cout<<endl;
01735                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
01736                 cout<<endl;
01737 
01738 #endif
01739 
01740                 return(_TRUE);
01741         }
01742 
01743 
01744         //Public Function 3462
01745         int BipartiteGraphOrdering::DynamicLargestFirstOrdering()
01746         {
01747                 if(CheckVertexOrdering("DYNAMIC_LARGEST_FIRST"))
01748                 {
01749                         return(_TRUE);
01750                 }
01751 
01752                 int i;
01753 
01754                 int _FOUND;
01755 
01756                 int i_HighestInducedVertexDegree, i_HighestInducedDegreeVertex;
01757 
01758                 int i_LeftVertexCount, i_RightVertexCount;
01759 
01760                 int i_InducedVertexDegree;
01761 
01762                 int i_InducedVertexDegreeCount;
01763 
01764                 int i_SelectedVertex, i_SelectedVertexCount;
01765 
01766                 vector<int> vi_InducedVertexDegree;
01767 
01768                 vector< list<int> > vli_GroupedInducedVertexDegree;
01769 
01770                 vector< list<int>::iterator > vlit_VertexLocation;
01771 
01772                 list<int>::iterator lit_ListIterator;
01773 
01774                 i_LeftVertexCount = STEP_DOWN((signed) m_vi_LeftVertices.size());
01775                 i_RightVertexCount = STEP_DOWN((signed) m_vi_RightVertices.size());
01776 
01777                 vi_InducedVertexDegree.clear();
01778 
01779                 vli_GroupedInducedVertexDegree.clear();
01780                 vli_GroupedInducedVertexDegree.resize((unsigned) i_LeftVertexCount + i_RightVertexCount);
01781 
01782                 vlit_VertexLocation.clear();
01783 
01784                 i_SelectedVertex = _UNKNOWN;
01785 
01786                 i_HighestInducedVertexDegree = _FALSE;
01787 
01788                 i_HighestInducedDegreeVertex = _UNKNOWN;
01789 
01790                 for(i=0; i<i_LeftVertexCount; i++)
01791                 {
01792                         i_InducedVertexDegree = m_vi_LeftVertices[STEP_UP(i)] - m_vi_LeftVertices[i];
01793 
01794                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
01795 
01796                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i);
01797 
01798                         vlit_VertexLocation.push_back(vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin());
01799 
01800                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01801                         {
01802                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01803 
01804                                 i_HighestInducedDegreeVertex = i;
01805                         }
01806                 }
01807 
01808                 for(i=0; i<i_RightVertexCount; i++)
01809                 {
01810                         i_InducedVertexDegree = m_vi_RightVertices[STEP_UP(i)] - m_vi_RightVertices[i];
01811 
01812                         vi_InducedVertexDegree.push_back(i_InducedVertexDegree);
01813 
01814                         vli_GroupedInducedVertexDegree[i_InducedVertexDegree].push_front(i + i_LeftVertexCount);
01815 
01816                         vlit_VertexLocation.push_back(vli_GroupedInducedVertexDegree[i_InducedVertexDegree].begin());
01817 
01818                         if(i_HighestInducedVertexDegree < i_InducedVertexDegree)
01819                         {
01820                                 i_HighestInducedVertexDegree = i_InducedVertexDegree;
01821 
01822                                 i_HighestInducedDegreeVertex = i + i_LeftVertexCount;
01823                         }
01824                 }
01825 
01826 
01827 #if DEBUG == 3462
01828 
01829                 int j;
01830 
01831                 cout<<endl;
01832                 cout<<"DEBUG 3462 | Bipartite Graph Coloring | Bipartite Smallest Last Ordering"<<endl;
01833                 cout<<endl;
01834 
01835                 for(i=0; i<STEP_UP(i_HighestInducedVertexDegree); i++)
01836                 {
01837                         cout<<"Degree "<<i<<"\t"<<" : ";
01838 
01839                         i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01840 
01841                         j = _FALSE;
01842 
01843                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
01844                         {
01845                                 if(j==STEP_DOWN(i_InducedVertexDegreeCount))
01846                                 {
01847                                         cout<<STEP_UP(*lit_ListIterator)<<" ("<<i_InducedVertexDegreeCount<<")";
01848                                 }
01849                                 else
01850                                 {
01851                                         cout<<STEP_UP(*lit_ListIterator)<<", ";
01852                                 }
01853 
01854                                 j++;
01855                         }
01856 
01857                          cout<<endl;
01858                 }
01859 
01860                 cout<<endl;
01861                 cout<<"[Highest Degree Vertex = "<<STEP_UP(i_HighestInducedDegreeVertex)<<"; Highest Induced Vertex Degree = "<<i_HighestInducedVertexDegree<<"]"<<endl;
01862                 cout<<endl;
01863 
01864 #endif
01865 
01866                 m_vi_OrderedVertices.clear();
01867 
01868                 i_SelectedVertexCount = _FALSE;
01869 
01870                 while(i_SelectedVertexCount < i_LeftVertexCount + i_RightVertexCount)
01871                 {
01872                         for(i=i_HighestInducedVertexDegree; i>=0; i--)
01873                         {
01874                                 i_InducedVertexDegreeCount = (signed) vli_GroupedInducedVertexDegree[i].size();
01875 
01876                                 if(i_InducedVertexDegreeCount == _FALSE)
01877                                 {
01878                                         continue;
01879                                 }
01880 
01881                                 if(i_HighestInducedDegreeVertex < i_LeftVertexCount)
01882                                 {
01883                                         _FOUND = _FALSE;
01884 
01885                                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
01886                                         {
01887                                                 if(*lit_ListIterator < i_LeftVertexCount)
01888                                                 {
01889                                                         i_SelectedVertex = *lit_ListIterator;
01890 
01891                                                         _FOUND = _TRUE;
01892 
01893                                                         break;
01894                                                 }
01895                                         }
01896 
01897                                         if(!_FOUND)
01898                                         {
01899                                                 i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
01900                                         }
01901 
01902                                         break;
01903                                 }
01904                                 else
01905                                 {
01906                                         _FOUND = _FALSE;
01907 
01908                                         for(lit_ListIterator = vli_GroupedInducedVertexDegree[i].begin(); lit_ListIterator != vli_GroupedInducedVertexDegree[i].end(); lit_ListIterator++)
01909                                         {
01910                                                 if(*lit_ListIterator >= i_LeftVertexCount)
01911                                                 {
01912                                                         i_SelectedVertex = *lit_ListIterator;
01913 
01914                                                         _FOUND = _TRUE;
01915 
01916                                                         break;
01917                                                 }
01918                                         }
01919 
01920                                         if(!_FOUND)
01921                                         {
01922                                                 i_SelectedVertex = vli_GroupedInducedVertexDegree[i].front();
01923                                         }
01924                                 }
01925 
01926                                 break;
01927                         }
01928 
01929                         if(i_SelectedVertex < i_LeftVertexCount)
01930                         {
01931                                 for(i=m_vi_LeftVertices[i_SelectedVertex]; i<m_vi_LeftVertices[STEP_UP(i_SelectedVertex)]; i++)
01932                                 {
01933                                         if(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] == _UNKNOWN)
01934                                         {
01935                                                 continue;
01936                                         }
01937 
01938                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].erase(vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount]);
01939 
01940                                         vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]);
01941 
01942                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].push_front(m_vi_Edges[i] + i_LeftVertexCount);
01943 
01944                                         vlit_VertexLocation[m_vi_Edges[i] + i_LeftVertexCount] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i] + i_LeftVertexCount]].begin();
01945                                 }
01946                         }
01947                         else
01948                         {
01949                                 for(i=m_vi_RightVertices[i_SelectedVertex - i_LeftVertexCount]; i<m_vi_RightVertices[STEP_UP(i_SelectedVertex - i_LeftVertexCount)]; i++)
01950                                 {
01951                                         if(vi_InducedVertexDegree[m_vi_Edges[i]] == _UNKNOWN)
01952                                         {
01953                                                 continue;
01954                                         }
01955 
01956                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].erase(vlit_VertexLocation[m_vi_Edges[i]]);
01957 
01958                                         vi_InducedVertexDegree[m_vi_Edges[i]] = STEP_DOWN(vi_InducedVertexDegree[m_vi_Edges[i]]);
01959 
01960                                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].push_front(m_vi_Edges[i]);
01961 
01962                                         vlit_VertexLocation[m_vi_Edges[i]] = vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[m_vi_Edges[i]]].begin();
01963                                 }
01964                         }
01965 
01966                         vli_GroupedInducedVertexDegree[vi_InducedVertexDegree[i_SelectedVertex]].erase(vlit_VertexLocation[i_SelectedVertex]);
01967 
01968                         vi_InducedVertexDegree[i_SelectedVertex] = _UNKNOWN;
01969 
01970                         m_vi_OrderedVertices.push_back(i_SelectedVertex);
01971 
01972                         i_SelectedVertexCount = STEP_UP(i_SelectedVertexCount);
01973                 }
01974 
01975 
01976 #if DEBUG == 3462
01977 
01978                 int i_OrderedVertexCount;
01979 
01980                 cout<<endl;
01981                 cout<<"DEBUG 3462 | Bipartite Graph Coloring | Bipartite Dynamic Largest First Ordering"<<endl;
01982                 cout<<endl;
01983 
01984                 i_OrderedVertexCount = (signed) m_vi_OrderedVertices.size();
01985 
01986                 for(i=0; i<i_OrderedVertexCount; i++)
01987                 {
01988                         if(i == STEP_DOWN(i_OrderedVertexCount))
01989                         {
01990                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<" ("<<i_OrderedVertexCount<<")"<<endl;
01991                         }
01992                         else
01993                         {
01994                                 cout<<STEP_UP(m_vi_OrderedVertices[i])<<", ";
01995                         }
01996                 }
01997 
01998                 cout<<endl;
01999                 cout<<"[Ordered Vertex Count = "<<i_OrderedVertexCount<<"/"<<i_LeftVertexCount + i_RightVertexCount<<"]"<<endl;
02000                 cout<<endl;
02001 
02002 #endif
02003 
02004                 return(_TRUE);
02005         }
02006 
02007         //Public Function 3463
02008         string BipartiteGraphOrdering::GetVertexOrderingVariant()
02009         {
02010 
02011                 if(m_s_VertexOrderingVariant.compare("NATURAL") == 0)
02012                 {
02013                         return("Natural");
02014                 }
02015                 else
02016                 if(m_s_VertexOrderingVariant.compare("LARGEST_FIRST") == 0)
02017                 {
02018                         return("Largest First");
02019                 }
02020                 else
02021                 if(m_s_VertexOrderingVariant.compare("SMALLEST_LAST") == 0)
02022                 {
02023                         return("Smallest Last");
02024                 }
02025                 else
02026                 if(m_s_VertexOrderingVariant.compare("INCIDENCE_DEGREE") == 0)
02027                 {
02028                         return("Incidence Degree");
02029                 }
02030                 else
02031                 if(m_s_VertexOrderingVariant.compare("SELECTVE_LARGEST_FIRST") == 0)
02032                 {
02033                         return("Selective Largest First");
02034                 }
02035                 else
02036                 if(m_s_VertexOrderingVariant.compare("SELECTVE_SMALLEST_FIRST") == 0)
02037                 {
02038                         return("Selective Smallest Last");
02039                 }
02040                 else
02041                 if(m_s_VertexOrderingVariant.compare("SELECTIVE_INCIDENCE_DEGREE") == 0)
02042                 {
02043                         return("Selective Incidence Degree");
02044                 }
02045                 else
02046                 if(m_s_VertexOrderingVariant.compare("DYNAMIC_LARGEST_FIRST") == 0)
02047                 {
02048                         return("Dynamic Largest First");
02049                 }
02050                 else
02051                 {
02052                         return("Unknown");
02053                 }
02054         }
02055 
02056         //Public Function 3464
02057         void BipartiteGraphOrdering::GetOrderedVertices(vector<int> &output)
02058         {
02059                 output = (m_vi_OrderedVertices);
02060         }
02061 
02062         int BipartiteGraphOrdering::OrderVertices(string s_OrderingVariant) {
02063                 s_OrderingVariant = toUpper(s_OrderingVariant);
02064 
02065                 if((s_OrderingVariant.compare("NATURAL") == 0))
02066                 {
02067                         return(NaturalOrdering());
02068                 }
02069                 else
02070                 if((s_OrderingVariant.compare("LARGEST_FIRST") == 0))
02071                 {
02072                         return(LargestFirstOrdering());
02073                 }
02074                 else
02075                 if((s_OrderingVariant.compare("DYNAMIC_LARGEST_FIRST") == 0))
02076                 {
02077                         return(DynamicLargestFirstOrdering());
02078                 }
02079                 else
02080                 if((s_OrderingVariant.compare("SMALLEST_LAST") == 0))
02081                 {
02082                         return(SmallestLastOrdering());
02083                 }
02084                 else
02085                 if((s_OrderingVariant.compare("INCIDENCE_DEGREE") == 0))
02086                 {
02087                         return(IncidenceDegreeOrdering());
02088                 }
02089                 else
02090                 if((s_OrderingVariant.compare("RANDOM") == 0))
02091                 {
02092                         return(RandomOrdering());
02093                 }
02094                 else
02095                 {
02096                         cerr<<endl;
02097                         cerr<<"Unknown Ordering Method: "<<s_OrderingVariant;
02098                         cerr<<endl;
02099                 }
02100 
02101                 return(_TRUE);
02102         }
02103 
02104         void BipartiteGraphOrdering::PrintVertexOrdering() {
02105                 cout<<"PrintVertexOrdering() "<<m_s_VertexOrderingVariant<<endl;
02106                 for(unsigned int i=0; i<m_vi_OrderedVertices.size();i++) {
02107                         //printf("\t [%d] %d \n", i, m_vi_OrderedVertices[i]);
02108                         cout<<"\t["<<setw(5)<<i<<"] "<<setw(5)<<m_vi_OrderedVertices[i]<<endl;
02109                 }
02110                 cout<<endl;
02111         }
02112         
02113         double BipartiteGraphOrdering::GetVertexOrderingTime() {
02114           return m_d_OrderingTime;
02115         }
02116 
02117 }

Generated on Tue Sep 7 2010 15:28:10 for ColPack by  doxygen 1.7.1