ColPack
GraphColoring/GraphInputOutput.cpp
Go to the documentation of this file.
00001 /************************************************************************************
00002     Copyright (C) 2005-2008 Assefaw H. Gebremedhin, Arijit Tarafdar, Duc Nguyen,
00003     Alex Pothen
00004 
00005     This file is part of ColPack.
00006 
00007     ColPack is free software: you can redistribute it and/or modify
00008     it under the terms of the GNU Lesser General Public License as published
00009     by the Free Software Foundation, either version 3 of the License, or
00010     (at your option) any later version.
00011 
00012     ColPack is distributed in the hope that it will be useful,
00013     but WITHOUT ANY WARRANTY; without even the implied warranty of
00014     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00015     GNU Lesser General Public License for more details.
00016 
00017     You should have received a copy of the GNU Lesser General Public License
00018     along with ColPack.  If not, see <http://www.gnu.org/licenses/>.
00019 ************************************************************************************/
00020 
00021 #include "ColPackHeaders.h"
00022 
00023 using namespace std;
00024 
00025 namespace ColPack
00026 {
00027         //Private Function 1201
00028         int GraphInputOutput::ParseWidth(string FortranFormat)
00029         {
00030                 int i;
00031 
00032                 int LetterCount;
00033 
00034                 char PresentLetter;
00035 
00036                 string FieldWidth;
00037 
00038                 boolean FOUND;
00039 
00040                 FieldWidth.clear();
00041 
00042                 FOUND = FALSE;
00043 
00044                 LetterCount = (signed) FortranFormat.size();
00045 
00046                 for(i=0; i<LetterCount; i++)
00047                 {
00048                         PresentLetter = FortranFormat[i];
00049 
00050                         if(FOUND == TRUE)
00051                         {
00052                           FieldWidth += PresentLetter;
00053                         }
00054 
00055                         if(PresentLetter == 'I' || PresentLetter == 'Z' || PresentLetter == 'F' || PresentLetter == 'E' || PresentLetter == 'G' || PresentLetter == 'D' || PresentLetter == 'L' || PresentLetter == 'A')
00056                         {
00057                                 FOUND = TRUE;
00058                         }
00059                         else
00060                         if(PresentLetter == '.' || PresentLetter == ')')
00061                         {
00062                                 FOUND = FALSE;
00063 
00064                                  break;
00065                         }
00066                 }
00067 
00068                 return(atoi(FieldWidth.c_str()));
00069         }
00070 
00071 
00072         //Private Function 1202
00073         void GraphInputOutput::CalculateVertexDegrees()
00074         {
00075                 int i_VertexCount = STEP_DOWN((signed) m_vi_Vertices.size());
00076 
00077                 for(int i = 0; i < i_VertexCount; i++)
00078                 {
00079                         int i_VertexDegree = m_vi_Vertices[i + 1] - m_vi_Vertices[i];
00080 
00081                         if(m_i_MaximumVertexDegree < i_VertexDegree)
00082                         {
00083                                 m_i_MaximumVertexDegree = i_VertexDegree;
00084                         }
00085 
00086                         if(m_i_MinimumVertexDegree == _UNKNOWN)
00087                         {
00088                                 m_i_MinimumVertexDegree = i_VertexDegree;
00089                         }
00090                         else
00091                         if(m_i_MinimumVertexDegree > i_VertexDegree)
00092                         {
00093                                 m_i_MinimumVertexDegree = i_VertexDegree;
00094                         }
00095                 }
00096 
00097                 m_d_AverageVertexDegree = (double)m_vi_Edges.size()/i_VertexCount;
00098 
00099                 return;
00100 
00101         }
00102 
00103 
00104         //Public Constructor 1251
00105         GraphInputOutput::GraphInputOutput()
00106         {
00107                 Clear();
00108 
00109                 GraphCore::Clear();
00110         }
00111 
00112 
00113         //Public Destructor 1252
00114         GraphInputOutput::~GraphInputOutput()
00115         {
00116                 Clear();
00117         }
00118 
00119         //Virtual Function 1254
00120         void GraphInputOutput::Clear()
00121         {
00122                 GraphCore::Clear();
00123 
00124                 return;
00125         }
00126 
00127         //Public Function 1255
00128         string GraphInputOutput::GetInputFile()
00129         {
00130                 return(m_s_InputFile);
00131         }
00132         
00133         int GraphInputOutput::WriteMatrixMarket(string s_OutputFile, bool b_getStructureOnly) {
00134                 ofstream out (s_OutputFile.c_str());
00135                 if(!out) {
00136                         cout<<"Error creating file: \""<<s_OutputFile<<"\""<<endl;
00137                         exit(1);
00138                 }
00139                 
00140                 bool b_printValue = ( (!b_getStructureOnly) && (m_vd_Values.size()==m_vi_Edges.size()) );
00141                 int i_NumOfLines = 0;
00142                 int max = m_vi_Vertices.size()-1;
00143                 
00144                 out<<"%%MatrixMarket matrix coordinate real symmetric"<<endl;
00145                 
00146                 //Count i_NumOfLines
00147                 for(int i = 1; i<max;i++) {
00148                   for(int j = m_vi_Vertices[i]; j < m_vi_Vertices[i+1]; j++) {
00149                     //Only print out the entries in the lower triangular portion of the matrix
00150                     if(i>m_vi_Edges[j]) {
00151                       i_NumOfLines++;
00152                     }
00153                   }
00154                 }
00155                 
00156                 out<<m_vi_Vertices.size()-1<<" "<<m_vi_Vertices.size()-1<<" "<< i_NumOfLines<<endl;
00157                 
00158                 out<<setprecision(10)<<scientific<<showpoint;
00159                 for(int i = 1; i<max;i++) {
00160                   for(int j = m_vi_Vertices[i]; j < m_vi_Vertices[i+1]; j++) {
00161                     //Only print out the entries in the lower triangular portion of the matrix
00162                     if(i>m_vi_Edges[j]) {
00163                       out<<i+1<<" "<<m_vi_Edges[j]+1;
00164                       if (b_printValue) out<<" "<<m_vd_Values[j];
00165                       out<<endl;
00166                     }
00167                   }
00168                 }
00169                 
00170                 return 0;
00171         }
00172 
00173 
00174         //Public Function 1257
00175         int GraphInputOutput::ReadMatrixMarketAdjacencyGraph(string s_InputFile, bool b_getStructureOnly)
00176         {
00177                 //Pause();
00178                 Clear();
00179 
00180                 m_s_InputFile=s_InputFile;
00181 
00182                 //initialize local data
00183                 int col=0, row=0, rowIndex=0, colIndex=0;
00184                 int entry_counter = 0, num_of_entries = 0;
00185                 bool value_not_specified = false;
00186                 //int num=0, numCount=0;
00187                 double value;
00188                 bool b_getValue = !b_getStructureOnly, b_symmetric;
00189                 istringstream in2;
00190                 string line="";
00191                 map<int,vector<int> > nodeList;
00192                 map<int,vector<double> > valueList;
00193 
00194                 //READ IN BANNER
00195                 MM_typecode matcode;
00196                 FILE *f;
00197                 if ((f = fopen(m_s_InputFile.c_str(), "r")) == NULL)  {
00198                   cout<<m_s_InputFile<<" not Found!"<<endl;
00199                   exit(1);
00200                 }
00201                 else cout<<"Found file "<<m_s_InputFile<<endl;
00202 
00203                 if (mm_read_banner(f, &matcode) != 0)
00204                 {
00205                     printf("Could not process Matrix Market banner.\n");
00206                     exit(1);
00207                 }
00208 
00209 
00210                 if( mm_is_pattern(matcode) ) {
00211                   b_getValue = false;
00212                 }
00213                 if(mm_is_symmetric(matcode)) {
00214                   b_symmetric = true;
00215                 }
00216                 else b_symmetric = false;
00217                 //Check and make sure that the input file is supported
00218                 char * result = mm_typecode_to_str(matcode);
00219                 printf("Graph of Market Market type: [%s]\n", result);
00220                 free(result);
00221                 if (b_getValue) printf("\t Graph structure and VALUES will be read\n");
00222                 else printf("\t Read graph struture only. Values will NOT be read\n");
00223                 if( !( mm_is_coordinate(matcode) && (mm_is_symmetric(matcode) || mm_is_general(matcode) ) && ( mm_is_real(matcode) || mm_is_pattern(matcode) || mm_is_integer(matcode) ) ) ) {
00224                   printf("Sorry, this application does not support this type.");
00225                   exit(-1);
00226                 }
00227 
00228                 fclose(f);
00229                 //DONE - READ IN BANNER
00230 
00231                 // FIND OUT THE SIZE OF THE MATRIX
00232                 ifstream in (m_s_InputFile.c_str());
00233                 if(!in) {
00234                         cout<<m_s_InputFile<<" not Found!"<<endl;
00235                         exit(1);
00236                 }
00237                 else {
00238                   //cout<<"Found file "<<m_s_InputFile<<endl;
00239                 }
00240 
00241                 getline(in,line);
00242                 while(line.size()>0&&line[0]=='%') {//ignore comment line
00243                         getline(in,line);
00244                 }
00245                 in2.str(line);
00246                 in2>>row>>col>>num_of_entries;
00247                 //cout<<"row="<<row<<"; col="<<col<<"; num_of_entries="<<num_of_entries<<endl;
00248 
00249                 if(row!=col) {
00250                         cout<<"* WARNING: GraphInputOutput::ReadMatrixMarketAdjacencyGraph()"<<endl;
00251                         cout<<"*\t row!=col. This is not a square matrix. Can't process."<<endl;
00252                         return _FALSE;
00253                 }
00254                 // DONE - FIND OUT THE SIZE OF THE MATRIX
00255 
00256                 while(!in.eof() && entry_counter<num_of_entries) //there should be (num_of_entries+1) lines in the input file (excluding the comments)
00257                 {
00258                         getline(in,line);
00259                         entry_counter++;
00260                         //cout<<"Line "<<entry_counter<<"="<<line<<endl;
00261                         if(line!="")
00262                         {
00263                                 in2.clear();
00264                                 in2.str(line);
00265                                 
00266                                 //value =-999999999;
00267                                 //value_not_specified=false;
00268                                 
00269                                 in2 >> rowIndex >> colIndex >> value;
00270                                 
00271                                 /*if(value == -999999999 && in2.eof()) {                
00272                                   // "value" entry is not specified
00273                                   value_not_specified = true;
00274                                   if(b_getValue) {
00275                                     cerr<<"ERROR: GraphInputOutput::ReadMatrixMarketAdjacencyGraph()"<<endl;
00276                                     cerr<<"\t entry #"<<entry_counter<<": \""<< line <<"\""<<endl;
00277                                     cerr<<"\t \"value\" entry is not specified"<<endl;
00278                                     exit(1);
00279                                   }
00280                                 }
00281                                 else if (value == 0) {
00282                                   continue;
00283                                 }
00284                                 */
00285 
00286                                 rowIndex--;
00287                                 colIndex--;
00288                                 
00289                                 if(b_symmetric) {
00290                                         if(rowIndex>colIndex) {
00291                                                 //cout<<"\t"<<setw(4)<<rowIndex<<setw(4)<<colIndex<<setw(4)<<nz_counter<<endl;
00292                                                 nodeList[rowIndex].push_back(colIndex);
00293                                                 nodeList[colIndex].push_back(rowIndex);
00294                                                 
00295                                                 if(b_getValue) {
00296                                                         //cout<<"Value = "<<value<<endl;
00297                                                         valueList[rowIndex].push_back(value);
00298                                                         valueList[colIndex].push_back(value);
00299                                                 }
00300 
00301                                         }
00302                                         else if (rowIndex == colIndex) {
00303                                           continue;
00304                                         }
00305                                         else {
00306                                                 cerr<<"* WARNING: GraphInputOutput::ReadMatrixMarketAdjacencyGraph()"<<endl;
00307                                                 cerr<<"\t Found a nonzero in the upper triangular. A symmetric Matrix Market file format should only specify the nonzeros in the lower triangular."<<endl;
00308                                                 exit( -1);
00309                                         }
00310                                 }
00311                                 else { // !b_symmetric
00312                                         if(rowIndex!=colIndex) {
00313                                                 //cout<<"\t"<<setw(4)<<rowIndex<<setw(4)<<colIndex<<setw(4)<<nz_counter<<endl;
00314                                                 nodeList[rowIndex].push_back(colIndex);
00315                                                 
00316                                                 if(b_getValue) {
00317                                                         //cout<<"Value = "<<value<<endl;
00318                                                         valueList[rowIndex].push_back(value);
00319                                                 }
00320                                         }
00321                                         else { //rowIndex == colIndex
00322                                                 continue;
00323                                         }
00324                                 }
00325                                 
00326                         }
00327                         else
00328                         {
00329                                 cerr<<"* WARNING: GraphInputOutput::ReadMatrixMarketAdjacencyGraph()"<<endl;
00330                                 cerr<<"*\t Entry == \"\" at row "<<entry_counter<<". Empty line. Wrong input format. Can't process."<<endl;
00331                                 cerr<<"\t # entries so far: "<<entry_counter<<"/"<<num_of_entries<<endl;
00332                                 return _FALSE;
00333                         }
00334                 }
00335                 if(entry_counter<num_of_entries) { //entry_counter should be == num_of_entries
00336                                 cerr<<"* WARNING: GraphInputOutput::ReadMatrixMarketAdjacencyGraph()"<<endl;
00337                                 cerr<<"*\t entry_counter<num_of_entries. Wrong input format. Can't process."<<endl;
00338                                 cerr<<"\t # entries so far: "<<entry_counter<<"/"<<num_of_entries<<endl;
00339                                 return _FALSE;
00340                 }
00341 
00342 
00343                 //now construct the graph
00344                 m_vi_Vertices.push_back(m_vi_Edges.size()); //m_vi_Edges.size() == 0 at this point
00345 
00346                 for(int i=0;i<row; i++) {
00347                         m_vi_Edges.insert(m_vi_Edges.end(),nodeList[i].begin(),nodeList[i].end());
00348                         m_vi_Vertices.push_back(m_vi_Edges.size());
00349                 }
00350                 if(b_getValue) {
00351                         for(int i=0;i<row; i++) {
00352                                 m_vd_Values.insert(m_vd_Values.end(),valueList[i].begin(),valueList[i].end());
00353                         }
00354                 }
00355 
00356 //PrintGraph();
00357 //cout<<"DONE PrintGraph();"<<endl;
00358 
00359                 CalculateVertexDegrees();
00360 
00361                 return(_TRUE);
00362         }
00363 
00364 
00365         int GraphInputOutput::ReadHarwellBoeingAdjacencyGraph(string s_InputFile) {
00366                 Clear();
00367 
00368                 m_s_InputFile=s_InputFile;
00369                 ifstream in (m_s_InputFile.c_str());
00370 
00371                 if(!in)
00372                 {
00373                         cout<<"File "<<m_s_InputFile<<" Not Found"<<endl;
00374                         return _FALSE;
00375                 }
00376                 else
00377                 {
00378                         cout<<"Found File "<<m_s_InputFile<<endl;
00379                 }
00380                 int i_Dummy, i, j;
00381                 int num, counter;
00382                 double d;
00383                 int nnz;
00384                 string line, num_string;
00385                 istringstream iin;
00386                 vector< vector<int> > vvi_VertexAdjacency;
00387                 vector< vector<double> > vvd_Values;
00388                 vector<int> vi_ColumnStartPointers, vi_RowIndices;
00389                 vector<double> vd_Values;
00390                 
00391                 //ignore the first line, which is the tittle and key
00392                 getline(in, line);
00393                 
00394                 // Get line 2
00395                 int TOTCRD; // (ignored) Total number of lines excluding header
00396                 int PTRCRD; // (ignored) Number of lines for pointers
00397                 int INDCRD; // (ignored) Number of lines for row (or variable) indices
00398                 int VALCRD; // Number of lines for numerical values. VALCRD == 0 if no values is presented
00399                 int RHSCRD; // (ignored) Number of lines for right-hand sides. RHSCRD == 0 if no right-hand side data is presented
00400                 
00401                 getline(in, line);
00402                 iin.clear();
00403                 iin.str(line);
00404                 iin >> TOTCRD >> PTRCRD >> INDCRD >> VALCRD >> RHSCRD;
00405                 
00406                 // Get line 3
00407                 string MXTYPE; //Matrix type. We only accept: (R | P) (S | U) (A)
00408                 int NROW; // Number of rows (or left vertices)
00409                 int NCOL; // Number of columns (or  right vertices)
00410                 int NNZERO; // Number of nonzeros 
00411                             // in case of symmetric matrix, it is the number of nonzeros IN THE UPPER TRIANGULAR including the diagonal
00412                 int NELTVL; // (ignored) Number of elemental matrix entries (zero in the case of assembled matrices) 
00413                 bool b_symmetric; // true if this matrix is symmetric (MXTYPE[1] == 'S'), false otherwise.
00414                 
00415                 getline(in, line);
00416                 iin.clear();
00417                 iin.str(line);
00418                 iin >> MXTYPE >> NROW >> NCOL >> NNZERO >> NELTVL;
00419                 // We only accept MXTYPE = (R|P)(S|U)A
00420                 if(MXTYPE[0] == 'C') { //Complex matrix 
00421                   cerr<<"ERR: Complex matrix format is not supported"<<endl;
00422                   exit(-1);
00423                 }
00424                 if(MXTYPE[1] == 'S') {
00425                   b_symmetric = true;
00426                 }
00427                 else {
00428                   b_symmetric = false;
00429                   if(MXTYPE[1] != 'U') { //H, Z, R types are not supported
00430                     cerr<<"ERR: Matrix format is not supported. MXTYPE[1] != 'S' && MXTYPE[1] != 'U'"<<endl;
00431                     exit(-1);
00432                   }
00433                 }
00434                 if(MXTYPE[2] == 'E') { //Elemental matrices (unassembled) 
00435                   cerr<<"ERR: Elemental matrices (unassembled) format is not supported"<<endl;
00436                   exit(-1);
00437                 }
00438 
00439                 if(NROW != NCOL) {
00440                         cout<<"* WARNING: GraphInputOutput::ReadHarwellBoeingAdjacencyGraph()"<<endl;
00441                         cout<<"*\t row!=col. This is not a square matrix. Can't process."<<endl;
00442                         return _FALSE;
00443                 }
00444 
00445                 // Ignore line 4 for now
00446                 getline(in, line);
00447                 
00448                 //If the right-hand side data is presented, ignore the 5th header line
00449                 if(RHSCRD) getline(in, line);
00450                 
00451                 //Initialize data structures
00452                 m_vi_Vertices.clear();
00453                 m_vi_Vertices.resize(NROW+1, _UNKNOWN);
00454                 vvi_VertexAdjacency.clear();
00455                 vvi_VertexAdjacency.resize(NROW);
00456                 vvd_Values.clear();
00457                 vvd_Values.resize(NROW);
00458                 vi_ColumnStartPointers.clear();
00459                 vi_ColumnStartPointers.resize(NCOL+1);
00460                 vi_RowIndices.clear();
00461                 vi_RowIndices.resize(NNZERO);
00462                 vd_Values.clear();
00463                 vd_Values.resize(NNZERO);
00464                 
00465                 // get the 2nd data block: column start pointers                
00466                 for(int i=0; i<NCOL+1; i++) {
00467                   in>> vi_ColumnStartPointers[i];
00468                 }
00469                 
00470                 // get the 3rd data block: row (or variable) indices,
00471                 for(i=0; i<NNZERO; i++) {
00472                   in >> num;
00473                   vi_RowIndices[i] = num-1;
00474                 }
00475 
00476                 // get the 4th data block: numerical values 
00477                 if(VALCRD !=0) {
00478                   for(i=0; i<NNZERO; i++) {
00479                     in >> num_string;
00480                     ConvertHarwellBoeingDouble(num_string);
00481                     iin.clear();
00482                     iin.str(num_string);
00483                     iin >> d;
00484                     vd_Values[i] = d;
00485                   }
00486                 }
00487 
00488                 //populate vvi_VertexAdjacency & vvd_Values
00489                 nnz = 0;
00490                 counter = 0;
00491                 for(i=0; i<NCOL; i++) {
00492                   for(j=vi_ColumnStartPointers[i]; j< vi_ColumnStartPointers[i+1]; j++) {
00493                     num = vi_RowIndices[counter];
00494                     d = vd_Values[counter];
00495                     
00496                     if(num != i) {
00497                       if(b_symmetric) {
00498                         vvi_VertexAdjacency[i].push_back(num);                  
00499                         vvi_VertexAdjacency[num].push_back(i);
00500                         
00501                         if(VALCRD !=0) {
00502                           vvd_Values[i].push_back(d);
00503                           vvd_Values[num].push_back(d);
00504                         }
00505                         
00506                         nnz+=2;
00507                       }
00508                       else { // !b_symmetric
00509                         vvi_VertexAdjacency[i].push_back(num);
00510                         if(VALCRD !=0) vvd_Values[i].push_back(d);
00511                         nnz++;
00512                       }
00513                     }
00514                     counter++;
00515                   }
00516                 }
00517                 
00518                 m_vi_Edges.clear();
00519                 m_vi_Edges.resize(nnz, _UNKNOWN);
00520                 if(VALCRD !=0) {
00521                   m_vd_Values.clear();
00522                   m_vd_Values.resize(nnz, _UNKNOWN);
00523                 }
00524                 //populate the m_vi_Vertices, their Edges and Values at the same time
00525                 m_vi_Vertices[0]=0;
00526                 for(i=0; i<NROW; i++) {
00527                   for(j=0; j<vvi_VertexAdjacency[i].size();j++) {
00528                     m_vi_Edges[m_vi_Vertices[i]+j] = vvi_VertexAdjacency[i][j];
00529                     if(VALCRD !=0) m_vd_Values[m_vi_Vertices[i]+j] = vvd_Values[i][j];
00530                   }
00531                   
00532                   m_vi_Vertices[i+1] = m_vi_Vertices[i]+vvi_VertexAdjacency[i].size();
00533                 }
00534                 
00535                 PrintGraph();
00536                 Pause();
00537                                 
00538           return 0;
00539         }
00540 
00541         int GraphInputOutput::ReadMeTiSAdjacencyGraph(string s_InputFile)
00542         {
00543                 istringstream in2;
00544                 int i, j;
00545 
00546                 int i_LineCount, i_TokenCount;
00547 
00548                 int i_Vertex;
00549 
00550                 int i_VertexCount, i_VertexDegree;
00551 
00552                 int i_EdgeCount;
00553 
00554                 int i_VertexWeights, i_EdgeWeights;
00555 
00556                 string _GAP(" ");
00557 
00558                 string s_InputLine;
00559 
00560                 ifstream InputStream;
00561 
00562                 vector<string> vs_InputTokens;
00563 
00564                 vector<double> vi_EdgeWeights;
00565                 vector<double> vi_VertexWeights;
00566 
00567                 vector< vector<int> > v2i_VertexAdjacency;
00568 
00569                 vector< vector<double> > v2i_VertexWeights;
00570 
00571                 Clear();
00572 
00573                 m_s_InputFile = s_InputFile;
00574 
00575                 InputStream.open(m_s_InputFile.c_str());
00576                 if(!InputStream) {
00577                         cout<<m_s_InputFile<<" not Found!"<<endl;
00578                         return (_FALSE);
00579                 }
00580                 else cout<<"Found file "<<m_s_InputFile<<endl;
00581 
00582                 vi_EdgeWeights.clear();
00583 
00584                 v2i_VertexWeights.clear();
00585 
00586                 i_VertexWeights = i_EdgeWeights = _FALSE;
00587 
00588                 i_LineCount = _FALSE;
00589 
00590                 do
00591                 {
00592                         getline(InputStream, s_InputLine);
00593 
00594                         if(!InputStream)
00595                         {
00596                                 break;
00597                         }
00598 
00599                         if(s_InputLine[0] == '%')
00600                         {
00601                                 continue;
00602                         }
00603 
00604                         if(i_LineCount == _FALSE)
00605                         {
00606                                 in2.clear();
00607                                 in2.str(s_InputLine);
00608                                 in2>>i_VertexCount>>i_EdgeCount;
00609 
00610                                 i_VertexWeights = _FALSE;
00611                                 i_EdgeWeights = _FALSE;
00612 
00613                                 if(!in2.eof())
00614                                 {
00615                                   int Weights;
00616                                   in2>>Weights;
00617                                         if(Weights == 1)
00618                                         {
00619                                                 i_EdgeWeights = _TRUE;
00620                                         }
00621                                         else
00622                                         if(Weights == 10)
00623                                         {
00624                                                 i_VertexWeights = _TRUE;
00625                                         }
00626                                         else
00627                                         if(Weights == 11)
00628                                         {
00629                                                 i_EdgeWeights = _TRUE;
00630                                                 i_VertexWeights = _TRUE;
00631                                         }
00632                            }
00633 
00634                                 if(!in2.eof())
00635                                 {
00636                                         in2>>i_VertexWeights ;
00637                                 }
00638 
00639                                 v2i_VertexAdjacency.clear();
00640                                 v2i_VertexAdjacency.resize((unsigned) i_VertexCount);
00641 
00642                                 i_LineCount++;
00643 
00644                         }
00645                         else
00646                         {
00647                                 in2.clear();
00648                                 //remove trailing space or tab in s_InputLine
00649                                 int input_end= s_InputLine.size() - 1;
00650                                 if(input_end>=0) {
00651                                   while(s_InputLine[input_end] == ' ' || s_InputLine[input_end] == '\t') input_end--;
00652                                 }
00653                                 if(input_end<0) s_InputLine = "";
00654                                 else s_InputLine = s_InputLine.substr(0, input_end+1);
00655                                 
00656                                 in2.str(s_InputLine);
00657                                 string tokens;
00658 
00659                                 vs_InputTokens.clear();
00660 
00661                                 while( !in2.eof() )
00662                                 {
00663                                         in2>>tokens;
00664                                         vs_InputTokens.push_back(tokens);
00665                                 }
00666 
00667                                 i_TokenCount = (signed) vs_InputTokens.size();
00668 
00669                                 vi_VertexWeights.clear();
00670 
00671                                 for(i=0; i<i_VertexWeights; i++)
00672                                 {
00673                                         vi_VertexWeights.push_back(atoi(vs_InputTokens[i].c_str()));
00674                                 }
00675 
00676                                 if(i_VertexWeights != _FALSE)
00677                                 {
00678                                         v2i_VertexWeights.push_back(vi_VertexWeights);
00679                                 }
00680 
00681                                 if(i_EdgeWeights == _FALSE)
00682                                 {
00683                                         for(i=i_VertexWeights; i<i_TokenCount; i++)
00684                                         {
00685                                                 if(vs_InputTokens[i] != "") {
00686                                                         i_Vertex = STEP_DOWN(atoi(vs_InputTokens[i].c_str()));
00687                                                         
00688                                                         //if(i_Vertex == -1) {
00689                                                         //  cout<<"i_Vertex == -1, i = "<<i<<", vs_InputTokens[i] = "<<vs_InputTokens[i]<<endl;
00690                                                         //  Pause();
00691                                                         //}
00692 
00693                                                         if(i_Vertex != STEP_DOWN(i_LineCount))
00694                                                         {
00695                                                                 v2i_VertexAdjacency[STEP_DOWN(i_LineCount)].push_back(i_Vertex);
00696                                                         }
00697                                                 }
00698                                         }
00699                                 }
00700                                 else
00701                                 {
00702                                         for(i=i_VertexWeights; i<i_TokenCount; i=i+2)
00703                                         {
00704                                                 i_Vertex = STEP_DOWN(atoi(vs_InputTokens[i].c_str()));
00705 
00706                                                 if(i_Vertex != STEP_DOWN(i_LineCount))
00707                                                 {
00708                                                         v2i_VertexAdjacency[STEP_DOWN(i_LineCount)].push_back(i_Vertex);
00709 
00710                                                         vi_EdgeWeights.push_back(STEP_DOWN(atof(vs_InputTokens[STEP_UP(i)].c_str())));
00711                                                 }
00712                                         }
00713                                 }
00714 
00715                                 i_LineCount++;
00716                         }
00717 
00718                 }
00719                 while(InputStream);
00720 
00721                 InputStream.close();
00722 
00723                 i_VertexCount = (signed) v2i_VertexAdjacency.size();
00724 
00725                 for(i=0; i<i_VertexCount; i++)
00726                 {
00727                         m_vi_Vertices.push_back((signed) m_vi_Edges.size());
00728 
00729                         i_VertexDegree = (signed) v2i_VertexAdjacency[i].size();
00730 
00731                         for(j=0; j<i_VertexDegree; j++)
00732                         {
00733                                 m_vi_Edges.push_back(v2i_VertexAdjacency[i][j]);
00734                         }
00735                 }
00736 
00737                 m_vi_Vertices.push_back((signed) m_vi_Edges.size());
00738 
00739                 CalculateVertexDegrees();
00740 
00741 #if DEBUG == 1259
00742 
00743                 cout<<endl;
00744                 cout<<"DEBUG 1259 | Graph Coloring | Vertex Adjacency | "<<m_s_InputFile<<endl;
00745                 cout<<endl;
00746 
00747                 i_EdgeCount = _FALSE;
00748 
00749                 for(i=0; i<i_VertexCount; i++)
00750                 {
00751                         cout<<"Vertex "<<STEP_UP(i)<<"\t"<<" : ";
00752 
00753                         i_VertexDegree = (signed) v2i_VertexAdjacency[i].size();
00754 
00755                         for(j=0; j<i_VertexDegree; j++)
00756                         {
00757                                 if(j == STEP_DOWN(i_VertexDegree))
00758                                 {
00759                                         cout<<STEP_UP(v2i_VertexAdjacency[i][j])<<" ("<<i_VertexDegree<<")";
00760 
00761                                         i_EdgeCount++;
00762                                 }
00763                                 else
00764                                 {
00765                                         cout<<STEP_UP(v2i_VertexAdjacency[i][j])<<", ";
00766 
00767                                         i_EdgeCount++;
00768                                 }
00769                         }
00770 
00771                         cout<<endl;
00772                 }
00773 
00774                 cout<<endl;
00775                 cout<<"[Vertices = "<<i_VertexCount<<"; Edges = "<<i_EdgeCount/2<<"]"<<endl;
00776                 cout<<endl;
00777 
00778 #endif
00779 
00780                 return(_TRUE);
00781 
00782         }
00783 
00784 
00785         //Public Function 1259
00786         int GraphInputOutput::ReadMeTiSAdjacencyGraph2(string s_InputFile)
00787         {
00788                 int i, j;
00789 
00790                 int i_LineCount, i_TokenCount;
00791 
00792                 int i_Vertex;
00793 
00794                 int i_VertexCount, i_VertexDegree;
00795 
00796                 int i_EdgeCount;
00797 
00798                 int i_VertexWeights, i_EdgeWeights;
00799 
00800                 string _GAP(" ");
00801 
00802                 string s_InputLine;
00803 
00804                 ifstream InputStream;
00805 
00806                 vector<string> vs_InputTokens;
00807 
00808                 vector<double> vi_EdgeWeights;
00809                 vector<double> vi_VertexWeights;
00810 
00811                 vector< vector<int> > v2i_VertexAdjacency;
00812 
00813                 vector< vector<double> > v2i_VertexWeights;
00814 
00815                 Clear();
00816 
00817                 m_s_InputFile = s_InputFile;
00818 
00819                 InputStream.open(m_s_InputFile.c_str());
00820                 if(!InputStream) {
00821                         cout<<m_s_InputFile<<" not Found!"<<endl;
00822                         return (_FALSE);
00823                 }
00824                 else cout<<"Found file "<<m_s_InputFile<<endl;
00825 
00826                 vi_EdgeWeights.clear();
00827 
00828                 v2i_VertexWeights.clear();
00829 
00830                 i_VertexWeights = i_EdgeWeights = _FALSE;
00831 
00832                 i_LineCount = _FALSE;
00833 
00834                 do
00835                 {
00836                         getline(InputStream, s_InputLine);
00837 
00838                         if(!InputStream)
00839                         {
00840                                 break;
00841                         }
00842 
00843                         if(s_InputLine[0] == '%')
00844                         {
00845                                 continue;
00846                         }
00847 
00848                         if(i_LineCount == _FALSE)
00849                         {
00850                                 StringTokenizer GapTokenizer(s_InputLine, _GAP);
00851 
00852                                 vs_InputTokens.clear();
00853 
00854                                 while(GapTokenizer.HasMoreTokens())
00855                                 {
00856                                         vs_InputTokens.push_back(GapTokenizer.GetNextToken());
00857                                 }
00858 
00859                                 i_VertexCount = atoi(vs_InputTokens[0].c_str());
00860                                 i_EdgeCount = atoi(vs_InputTokens[1].c_str());
00861 
00862                                 i_VertexWeights = _FALSE;
00863                                 i_EdgeWeights = _FALSE;
00864 
00865                                 if(vs_InputTokens.size() > 2)
00866                                 {
00867                                         if(atoi(vs_InputTokens[2].c_str()) == 1)
00868                                         {
00869                                                 i_EdgeWeights = _TRUE;
00870                                         }
00871                                         else
00872                                         if(atoi(vs_InputTokens[2].c_str()) == 10)
00873                                         {
00874                                                 i_VertexWeights = _TRUE;
00875                                         }
00876                                         else
00877                                         if(atoi(vs_InputTokens[2].c_str()) == 11)
00878                                         {
00879                                                 i_EdgeWeights = _TRUE;
00880                                                 i_VertexWeights = _TRUE;
00881                                         }
00882                            }
00883 
00884                                 if(vs_InputTokens.size() > 3)
00885                                 {
00886                                         i_VertexWeights = atoi(vs_InputTokens[3].c_str());
00887                                 }
00888 
00889                                 v2i_VertexAdjacency.clear();
00890                                 v2i_VertexAdjacency.resize((unsigned) i_VertexCount);
00891 
00892                                 i_LineCount++;
00893 
00894                         }
00895                         else
00896                         {
00897                                 StringTokenizer GapTokenizer(s_InputLine, _GAP);
00898 
00899                                 vs_InputTokens.clear();
00900 
00901                                 while(GapTokenizer.HasMoreTokens())
00902                                 {
00903                                         vs_InputTokens.push_back(GapTokenizer.GetNextToken());
00904                                 }
00905 
00906                                 i_TokenCount = (signed) vs_InputTokens.size();
00907 
00908                                 vi_VertexWeights.clear();
00909 
00910                                 for(i=0; i<i_VertexWeights; i++)
00911                                 {
00912                                         vi_VertexWeights.push_back(atoi(vs_InputTokens[i].c_str()));
00913                                 }
00914 
00915                                 if(i_VertexWeights != _FALSE)
00916                                 {
00917                                         v2i_VertexWeights.push_back(vi_VertexWeights);
00918                                 }
00919 
00920                                 if(i_EdgeWeights == _FALSE)
00921                                 {
00922                                         for(i=i_VertexWeights; i<i_TokenCount; i++)
00923                                         {
00924                                                 i_Vertex = STEP_DOWN(atoi(vs_InputTokens[i].c_str()));
00925 
00926                                                 if(i_Vertex != STEP_DOWN(i_LineCount))
00927                                                 {
00928                                                         v2i_VertexAdjacency[STEP_DOWN(i_LineCount)].push_back(i_Vertex);
00929                                                 }
00930                                         }
00931                                 }
00932                                 else
00933                                 {
00934                                         for(i=i_VertexWeights; i<i_TokenCount; i=i+2)
00935                                         {
00936                                                 i_Vertex = STEP_DOWN(atoi(vs_InputTokens[i].c_str()));
00937 
00938                                                 if(i_Vertex != STEP_DOWN(i_LineCount))
00939                                                 {
00940                                                         v2i_VertexAdjacency[STEP_DOWN(i_LineCount)].push_back(i_Vertex);
00941 
00942                                                         vi_EdgeWeights.push_back(STEP_DOWN(atof(vs_InputTokens[STEP_UP(i)].c_str())));
00943                                                 }
00944                                         }
00945                                 }
00946 
00947                                 i_LineCount++;
00948                         }
00949 
00950                 }
00951                 while(InputStream);
00952 
00953                 InputStream.close();
00954 
00955                 i_VertexCount = (signed) v2i_VertexAdjacency.size();
00956 
00957                 for(i=0; i<i_VertexCount; i++)
00958                 {
00959                         m_vi_Vertices.push_back((signed) m_vi_Edges.size());
00960 
00961                         i_VertexDegree = (signed) v2i_VertexAdjacency[i].size();
00962 
00963                         for(j=0; j<i_VertexDegree; j++)
00964                         {
00965                                 m_vi_Edges.push_back(v2i_VertexAdjacency[i][j]);
00966                         }
00967                 }
00968 
00969                 m_vi_Vertices.push_back((signed) m_vi_Edges.size());
00970 
00971                 CalculateVertexDegrees();
00972 
00973 #if DEBUG == 1259
00974 
00975                 cout<<endl;
00976                 cout<<"DEBUG 1259 | Graph Coloring | Vertex Adjacency | "<<m_s_InputFile<<endl;
00977                 cout<<endl;
00978 
00979                 i_EdgeCount = _FALSE;
00980 
00981                 for(i=0; i<i_VertexCount; i++)
00982                 {
00983                         cout<<"Vertex "<<STEP_UP(i)<<"\t"<<" : ";
00984 
00985                         i_VertexDegree = (signed) v2i_VertexAdjacency[i].size();
00986 
00987                         for(j=0; j<i_VertexDegree; j++)
00988                         {
00989                                 if(j == STEP_DOWN(i_VertexDegree))
00990                                 {
00991                                         cout<<STEP_UP(v2i_VertexAdjacency[i][j])<<" ("<<i_VertexDegree<<")";
00992 
00993                                         i_EdgeCount++;
00994                                 }
00995                                 else
00996                                 {
00997                                         cout<<STEP_UP(v2i_VertexAdjacency[i][j])<<", ";
00998 
00999                                         i_EdgeCount++;
01000                                 }
01001                         }
01002 
01003                         cout<<endl;
01004                 }
01005 
01006                 cout<<endl;
01007                 cout<<"[Vertices = "<<i_VertexCount<<"; Edges = "<<i_EdgeCount/2<<"]"<<endl;
01008                 cout<<endl;
01009 
01010 #endif
01011 
01012                 return(_TRUE);
01013 
01014         }
01015 
01016 
01017         //Public Function 1260
01018         int GraphInputOutput::PrintGraph()
01019         {
01020                 int i;
01021 
01022                 int i_VertexCount, i_EdgeCount;
01023 
01024                 i_VertexCount = (signed) m_vi_Vertices.size();
01025 
01026                 cout<<endl;
01027                 cout<<"Graph Coloring | Vertex List | "<<m_s_InputFile<<endl;
01028                 cout<<endl;
01029 
01030                 for(i=0; i<i_VertexCount; i++)
01031                 {
01032                         if(i == STEP_DOWN(i_VertexCount))
01033                         {
01034                                 cout<<STEP_UP(m_vi_Vertices[i])<<" ("<<i_VertexCount<<")"<<endl;
01035                         }
01036                         else
01037                         {
01038                                 cout<<STEP_UP(m_vi_Vertices[i])<<", ";
01039                         }
01040                 }
01041 
01042                 i_EdgeCount = (signed) m_vi_Edges.size();
01043 
01044                 cout<<endl;
01045                 cout<<"Graph Coloring | Edge List | "<<m_s_InputFile<<endl;
01046                 cout<<endl;
01047 
01048                 for(i=0; i<i_EdgeCount; i++)
01049                 {
01050                         if(i == STEP_DOWN(i_EdgeCount))
01051                         {
01052                                 cout<<STEP_UP(m_vi_Edges[i])<<" ("<<i_EdgeCount<<")"<<endl;
01053                         }
01054                         else
01055                         {
01056                                 cout<<STEP_UP(m_vi_Edges[i])<<", ";
01057                         }
01058                 }
01059 
01060                 if(m_vd_Values.size() > _FALSE)
01061                 {
01062 
01063                         cout<<endl;
01064                         cout<<"Graph Coloring | Nonzero List | "<<m_s_InputFile<<endl;
01065                         cout<<endl;
01066 
01067                         for(i=0; i<i_EdgeCount; i++)
01068                         {
01069                                 if(i == STEP_DOWN(i_EdgeCount))
01070                                 {
01071                                         cout<<m_vd_Values[i]<<" ("<<i_EdgeCount<<")"<<endl;
01072                                 }
01073                                 else
01074                                 {
01075                                         cout<<m_vd_Values[i]<<", ";
01076                                 }
01077                         }
01078 
01079                         cout<<endl;
01080                         cout<<"[Vertices = "<<STEP_DOWN(i_VertexCount)<<"; Edges = "<<i_EdgeCount/2<<"; Nonzeros = "<<i_EdgeCount/2<<"]"<<endl;
01081                         cout<<endl;
01082                 }
01083                 else
01084                 {
01085                         cout<<endl;
01086                         cout<<"[Vertices = "<<STEP_DOWN(i_VertexCount)<<"; Edges = "<<i_EdgeCount/2<<"]"<<endl;
01087                         cout<<endl;
01088 
01089                 }
01090 
01091                 return(_TRUE);
01092         }
01093 
01094 
01095         //Public Function 1261
01096         int GraphInputOutput::PrintGraphStructure()
01097         {
01098                 int i;
01099 
01100                 int i_VertexCount, i_EdgeCount;
01101 
01102                 i_VertexCount = (signed) m_vi_Vertices.size();
01103 
01104                 cout<<endl;
01105                 cout<<"Graph Coloring | Vertex List | "<<m_s_InputFile<<endl;
01106                 cout<<endl;
01107 
01108                 for(i=0; i<i_VertexCount; i++)
01109                 {
01110                         if(i == STEP_DOWN(i_VertexCount))
01111                         {
01112                                 cout<<STEP_UP(m_vi_Vertices[i])<<" ("<<i_VertexCount<<")"<<endl;
01113                         }
01114                         else
01115                         {
01116                                 cout<<STEP_UP(m_vi_Vertices[i])<<", ";
01117                         }
01118                 }
01119 
01120                 i_EdgeCount = (signed) m_vi_Edges.size();
01121 
01122                 cout<<endl;
01123                 cout<<"Graph Coloring | Edge List | "<<m_s_InputFile<<endl;
01124                 cout<<endl;
01125 
01126                 for(i=0; i<i_EdgeCount; i++)
01127                 {
01128                         if(i == STEP_DOWN(i_EdgeCount))
01129                         {
01130                                 cout<<STEP_UP(m_vi_Edges[i])<<" ("<<i_EdgeCount<<")"<<endl;
01131                         }
01132                         else
01133                         {
01134                                 cout<<STEP_UP(m_vi_Edges[i])<<", ";
01135                         }
01136                 }
01137 
01138                 cout<<endl;
01139                 cout<<"[Vertices = "<<STEP_DOWN(i_VertexCount)<<"; Edges = "<<i_EdgeCount/2<<"]"<<endl;
01140                 cout<<endl;
01141 
01142                 return(_TRUE);
01143         }
01144 
01145         int GraphInputOutput::PrintGraphStructure2()
01146         {
01147                 int i;
01148 
01149                 int i_VertexCount, i_EdgeCount;
01150 
01151                 i_VertexCount = (signed) m_vi_Vertices.size();
01152 
01153                 cout<<endl;
01154                 cout<<"PrintGraphStructure2() for graph: "<<m_s_InputFile<<endl;
01155                 cout<<"Format: Vertex id (# of edges): D1 neighbor #1, D1 neighbor #2, ... (all vertices is displayed using 1-based index)"<<endl;
01156                 cout<<endl;
01157 
01158                 for(i=0; i<i_VertexCount-1; i++)
01159                 {
01160                         cout<<"Vertex "<<STEP_UP(i)<<" ("<<m_vi_Vertices[i+1] - m_vi_Vertices[i]<<"): ";
01161                         
01162                         for(int j=m_vi_Vertices[i]; j<m_vi_Vertices[i+1]; j++)
01163                         {
01164                                 cout<<STEP_UP(m_vi_Edges[j])<<", ";
01165                         }
01166                         cout<<endl;
01167                 }
01168 
01169                 
01170                 cout<<endl;
01171 
01172                 return(_TRUE);
01173         }
01174 
01175 
01176         //Public Function 1262
01177         int GraphInputOutput::PrintMatrix()
01178         {
01179                 int i, j;
01180 
01181                 int i_VertexCount;
01182 
01183                 cout<<endl;
01184                 cout<<"Graph Coloring | Matrix Elements | "<<m_s_InputFile<<endl;
01185                 cout<<endl;
01186 
01187                 i_VertexCount = (signed) m_vi_Vertices.size();
01188 
01189                 for(i=0; i<i_VertexCount-1; i++)
01190                 {
01191                         for(j=m_vi_Vertices[i]; j<m_vi_Vertices[STEP_UP(i)]; j++)
01192                         {
01193                                 cout<<"Element["<<STEP_UP(i)<<"]["<<STEP_UP(m_vi_Edges[j])<<"] = "<<m_vd_Values[j]<<endl;
01194                         }
01195                 }
01196 
01197                 cout<<endl;
01198 
01199                 return(_TRUE);
01200         }
01201 
01202 
01203         //Public Function 1263
01204         int GraphInputOutput::PrintMatrix(vector<int> & vi_Vertices, vector<int> & vi_Edges, vector<double> & vd_Values)
01205         {
01206                 int i, j;
01207 
01208                 int i_VertexCount;
01209 
01210                 cout<<endl;
01211                 cout<<"Graph Coloring | Matrix Elements | "<<m_s_InputFile<<endl;
01212                 cout<<endl;
01213                 i_VertexCount = (signed) vi_Vertices.size();
01214 
01215                 for(i=0; i<i_VertexCount-1; i++)
01216                 {
01217                         for(j=vi_Vertices[i]; j<vi_Vertices[STEP_UP(i)]; j++)
01218                         {
01219                                 cout<<"Element["<<STEP_UP(i)<<"]["<<STEP_UP(vi_Edges[j])<<"] = "<<vd_Values[j]<<endl;
01220                         }
01221                 }
01222 
01223                 cout<<endl;
01224 
01225                 return(_TRUE);
01226         }
01227 
01228 
01229         //Public Function 1264
01230         void GraphInputOutput::PrintVertexDegrees()
01231         {
01232                 cout<<endl;
01233                 cout<<"Graph | "<<m_s_InputFile<<" | Maximum Vertex Degree | "<<m_i_MaximumVertexDegree<<endl;
01234                 cout<<"Graph | "<<m_s_InputFile<<" | Minimum Vertex Degree | "<<m_i_MinimumVertexDegree<<endl;
01235                 cout<<"Graph | "<<m_s_InputFile<<" | Average Vertex Degree | "<<m_d_AverageVertexDegree<<endl;
01236                 cout<<endl;
01237 
01238                 return;
01239         }
01240 
01241         int GraphInputOutput::BuildGraphFromRowCompressedFormat(unsigned int ** uip2_HessianSparsityPattern, int i_RowCount) {
01242           int i, j;
01243 
01244           int i_ElementCount, i_PositionCount;
01245 
01246           int i_HighestDegree;
01247 
01248 #if DEBUG == 1
01249 
01250           cout<<endl;
01251           cout<<"DEBUG | Graph Coloring | Sparsity Pattern"<<endl;
01252           cout<<endl;
01253 
01254           for(i=0; i<i_RowCount; i++)
01255             {
01256               cout<<i<<"\t"<<" : ";
01257 
01258               i_PositionCount = uip2_HessianSparsityPattern[i][0];
01259 
01260               for(j=0; j<i_PositionCount; j++)
01261                 {
01262                   if(j == STEP_DOWN(i_PositionCount))
01263                     {
01264                       cout<<uip2_HessianSparsityPattern[i][STEP_UP(j)]<<" ("<<i_PositionCount<<")";
01265                     }
01266                   else
01267                     {
01268                       cout<<uip2_HessianSparsityPattern[i][STEP_UP(j)]<<", ";
01269                     }
01270 
01271                 }
01272 
01273               cout<<endl;
01274             }
01275 
01276           cout<<endl;
01277 
01278 #endif
01279 
01280           m_vi_Vertices.clear();
01281           m_vi_Vertices.push_back(_FALSE);
01282 
01283           m_vi_Edges.clear();
01284 
01285           i_HighestDegree = _UNKNOWN;
01286 
01287           for(i=0; i<i_RowCount; i++)
01288             {
01289               i_ElementCount = _FALSE;
01290 
01291               i_PositionCount = uip2_HessianSparsityPattern[i][0];
01292 
01293               if(i_HighestDegree < i_PositionCount)
01294                 {
01295                   i_HighestDegree = i_PositionCount;
01296                 }
01297 
01298               for(j=0; j<i_PositionCount; j++)
01299                 {
01300                   if((signed) uip2_HessianSparsityPattern[i][STEP_UP(j)] != i)
01301                     {
01302                       m_vi_Edges.push_back((signed) uip2_HessianSparsityPattern[i][STEP_UP(j)]);
01303 
01304                       i_ElementCount++;
01305                     }
01306 
01307                 }
01308 
01309               m_vi_Vertices.push_back(m_vi_Vertices.back() + i_ElementCount);
01310             }
01311 
01312 #if DEBUG == 1
01313 
01314           int i_VertexCount, i_EdgeCount;
01315 
01316           cout<<endl;
01317           cout<<"DEBUG | Graph Coloring | Graph Format"<<endl;
01318           cout<<endl;
01319 
01320           cout<<"Vertices"<<"\t"<<" : ";
01321 
01322           i_VertexCount = (signed) m_vi_Vertices.size();
01323 
01324           for(i=0; i<i_VertexCount; i++)
01325             {
01326               if(i == STEP_DOWN(i_VertexCount))
01327                 {
01328                   cout<<m_vi_Vertices[i]<<" ("<<i_VertexCount<<")";
01329                 }
01330               else
01331                 {
01332                   cout<<m_vi_Vertices[i]<<", ";
01333                 }
01334             }
01335 
01336           cout<<endl;
01337 
01338           cout<<"Edges"<<"\t"<<" : ";
01339 
01340           i_EdgeCount = (signed) m_vi_Edges.size();
01341 
01342           for(i=0; i<i_EdgeCount; i++)
01343             {
01344               if(i == STEP_DOWN(i_EdgeCount))
01345                 {
01346                   cout<<m_vi_Edges[i]<<" ("<<i_EdgeCount<<")";
01347                 }
01348               else
01349                 {
01350                   cout<<m_vi_Edges[i]<<", ";
01351                 }
01352             }
01353 
01354           cout<<endl;
01355 
01356 #endif
01357 
01358           CalculateVertexDegrees();
01359 
01360           return(i_HighestDegree);
01361         }
01362 
01363         int GraphInputOutput::ReadAdjacencyGraph(string s_InputFile, string s_fileFormat)
01364         {
01365                 if (s_fileFormat == "AUTO_DETECTED" || s_fileFormat == "") {
01366                         File file(s_InputFile);
01367                         string fileExtension = file.GetFileExtension();
01368                         if (isHarwellBoeingFormat(fileExtension)) {
01369                                 cout<<"ReadHarwellBoeingAdjacencyGraph"<<endl;
01370                                 return ReadHarwellBoeingAdjacencyGraph(s_InputFile);
01371                         }
01372                         else if (isMeTiSFormat(fileExtension)) {
01373                                 cout<<"ReadMeTiSAdjacencyGraph"<<endl;
01374                                 return ReadMeTiSAdjacencyGraph(s_InputFile);
01375                         }
01376                         else if (isMatrixMarketFormat(fileExtension)) {
01377                                 cout<<"ReadMatrixMarketAdjacencyGraph"<<endl;
01378                                 return ReadMatrixMarketAdjacencyGraph(s_InputFile);
01379                         }
01380                         else { //other extensions
01381                                 cout<<"unfamiliar extension \""<<fileExtension<<"\", use ReadMatrixMarketAdjacencyGraph"<<endl;
01382                                 return ReadMatrixMarketAdjacencyGraph(s_InputFile);
01383                         }
01384                 }
01385                 else if (s_fileFormat == "MM") {
01386                         cout<<"ReadMatrixMarketAdjacencyGraph"<<endl;
01387                         return ReadMatrixMarketAdjacencyGraph(s_InputFile);
01388                 }
01389                 else if (s_fileFormat == "HB") {
01390                         cout<<"ReadHarwellBoeingAdjacencyGraph"<<endl;
01391                         return ReadHarwellBoeingAdjacencyGraph(s_InputFile);
01392                 }
01393                 else if (s_fileFormat == "MeTiS") {
01394                         cout<<"ReadMeTiSAdjacencyGraph"<<endl;
01395                         return ReadMeTiSAdjacencyGraph(s_InputFile);
01396                 }
01397                 else {
01398                         cerr<<"GraphInputOutput::ReadAdjacencyGraph s_fileFormat is not recognized"<<endl;
01399                         exit(1);
01400                 }
01401 
01402                 return(_TRUE);
01403         }
01404 
01405 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines