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