EpetraExt  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
EpetraExt_CrsMatrixIn.cpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ***********************************************************************
00003 //
00004 //     EpetraExt: Epetra Extended - Linear Algebra Services Package
00005 //                 Copyright (2011) Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ***********************************************************************
00040 //@HEADER
00041 #include "Epetra_ConfigDefs.h"
00042 #include "EpetraExt_CrsMatrixIn.h"
00043 #include "Epetra_Comm.h"
00044 #include "Epetra_CrsMatrix.h"
00045 #include "Epetra_Map.h"
00046 #include "Epetra_IntVector.h"
00047 #include "Epetra_IntSerialDenseVector.h"
00048 #include "Epetra_Import.h"
00049 #include "Epetra_Time.h"
00050 #include "Epetra_Util.h"
00051 
00052 #if defined(__PGI)
00053 #include <sstream>
00054 #endif
00055 
00058 // Functions to read a MatrixMarket file and load it into a matrix.
00059 // Adapted from
00060 //    Trilinos/packages/epetraext/src/inout/EpetraExt_CrsMatrixIn.cpp
00061 // Modified by Jon Berry and Karen Devine to make matrix reallocations
00062 // more efficient; rather than insert each non-zero separately, we
00063 // collect rows of non-zeros for insertion.
00064 // Modified by Karen Devine and Steve Plimpton to prevent all processors
00065 // from reading the same file at the same time (ouch!); chunks of the file
00066 // are read and broadcast by processor zero; each processor then extracts
00067 // its nonzeros from the chunk, sorts them by row, and inserts them into
00068 // the matrix.
00069 // The variable "chunk_read" can be changed to modify the size of the chunks
00070 // read from the file.  Larger values of chunk_read lead to faster execution
00071 // and greater memory use.
00074 
00075 using namespace EpetraExt;
00076 namespace EpetraExt {
00077 
00078 template<typename int_type>
00079 static void sort_three(int_type* list, int_type *parlista, double *parlistb,
00080                        int start, int end);
00081 
00083 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00084 int MatrixMarketFileToCrsMatrix(const char *filename,
00085     const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
00086 {
00087   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A));
00088   return(0);
00089 }
00090 
00092 int MatrixMarketFileToCrsMatrix(const char *filename,
00093     const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose)
00094 {
00095   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A,
00096                                                    0, 0, 0, 0, transpose));
00097   return(0);
00098 }
00100 
00101 int MatrixMarketFileToCrsMatrix(const char *filename,
00102     const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose,
00103     const bool verbose)
00104 {
00105   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A,
00106                                                    0, 0, 0, 0,
00107                                                    transpose, verbose));
00108   return(0);
00109 }
00110 
00112 int MatrixMarketFileToCrsMatrix(const char *filename,
00113     const Epetra_Map & rowMap, const Epetra_Map& rangeMap,
00114     const Epetra_Map& domainMap, Epetra_CrsMatrix * & A, const bool transpose,
00115     const bool verbose)
00116 {
00117   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
00118                                                    rowMap.Comm(), A,
00119                                                    &rowMap, 0,
00120                                                    &rangeMap, &domainMap,
00121                                                    transpose, verbose));
00122   return(0);
00123 }
00124 
00126 int MatrixMarketFileToCrsMatrix(const char *filename,
00127     const Epetra_Map & rowMap, Epetra_CrsMatrix * & A, const bool transpose,
00128     const bool verbose)
00129 {
00130   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
00131                                                    rowMap.Comm(), A,
00132                                                    &rowMap, 0, 0, 0,
00133                                                    transpose, verbose));
00134   return(0);
00135 }
00136 
00138 int MatrixMarketFileToCrsMatrix(const char *filename,
00139     const Epetra_Map & rowMap, const Epetra_Map & colMap,
00140     Epetra_CrsMatrix * & A, const bool transpose,
00141     const bool verbose)
00142 {
00143   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
00144                                                    rowMap.Comm(), A,
00145                                                    &rowMap, &colMap, 0, 0,
00146                                                    transpose, verbose));
00147   return(0);
00148 }
00149 
00151 int MatrixMarketFileToCrsMatrix(const char *filename,
00152     const Epetra_Map & rowMap, const Epetra_Map & colMap,
00153     const Epetra_Map& rangeMap, const Epetra_Map& domainMap,
00154     Epetra_CrsMatrix * & A, const bool transpose,
00155     const bool verbose)
00156 {
00157   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
00158                                                    rowMap.Comm(), A,
00159                                                    &rowMap, &colMap,
00160                                                    &rangeMap, &domainMap,
00161                                                    transpose, verbose));
00162   return(0);
00163 }
00164 #endif
00165 
00166 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00167 int MatrixMarketFileToCrsMatrix64(const char *filename,
00168     const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
00169 {
00170   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A));
00171   return(0);
00172 }
00173 
00175 int MatrixMarketFileToCrsMatrix64(const char *filename,
00176     const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose)
00177 {
00178   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A,
00179                                                    0, 0, 0, 0, transpose));
00180   return(0);
00181 }
00183 
00184 int MatrixMarketFileToCrsMatrix64(const char *filename,
00185     const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose,
00186     const bool verbose)
00187 {
00188   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A,
00189                                                    0, 0, 0, 0,
00190                                                    transpose, verbose));
00191   return(0);
00192 }
00193 
00195 int MatrixMarketFileToCrsMatrix64(const char *filename,
00196     const Epetra_Map & rowMap, const Epetra_Map& rangeMap,
00197     const Epetra_Map& domainMap, Epetra_CrsMatrix * & A, const bool transpose,
00198     const bool verbose)
00199 {
00200   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
00201                                                    rowMap.Comm(), A,
00202                                                    &rowMap, 0,
00203                                                    &rangeMap, &domainMap,
00204                                                    transpose, verbose));
00205   return(0);
00206 }
00207 
00209 int MatrixMarketFileToCrsMatrix64(const char *filename,
00210     const Epetra_Map & rowMap, Epetra_CrsMatrix * & A, const bool transpose,
00211     const bool verbose)
00212 {
00213   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
00214                                                    rowMap.Comm(), A,
00215                                                    &rowMap, 0, 0, 0,
00216                                                    transpose, verbose));
00217   return(0);
00218 }
00219 
00221 int MatrixMarketFileToCrsMatrix64(const char *filename,
00222     const Epetra_Map & rowMap, const Epetra_Map & colMap,
00223     Epetra_CrsMatrix * & A, const bool transpose,
00224     const bool verbose)
00225 {
00226   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
00227                                                    rowMap.Comm(), A,
00228                                                    &rowMap, &colMap, 0, 0,
00229                                                    transpose, verbose));
00230   return(0);
00231 }
00232 
00234 int MatrixMarketFileToCrsMatrix64(const char *filename,
00235     const Epetra_Map & rowMap, const Epetra_Map & colMap,
00236     const Epetra_Map& rangeMap, const Epetra_Map& domainMap,
00237     Epetra_CrsMatrix * & A, const bool transpose,
00238     const bool verbose)
00239 {
00240   EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
00241                                                    rowMap.Comm(), A,
00242                                                    &rowMap, &colMap,
00243                                                    &rangeMap, &domainMap,
00244                                                    transpose, verbose));
00245   return(0);
00246 }
00247 #endif
00248 
00250 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00251 int MatrixMarketFileToCrsMatrixHandle(const char *filename,
00252   const Epetra_Comm & comm,
00253   Epetra_CrsMatrix * & A,
00254   const Epetra_Map * rowMap,
00255   const Epetra_Map * colMap,
00256   const Epetra_Map * rangeMap,
00257   const Epetra_Map * domainMap,
00258   const bool transpose,
00259   const bool verbose
00260 )
00261 {
00262   const int chunk_read = 500000;  //  Modify this variable to change the
00263                                   //  size of the chunks read from the file.
00264   const int headerlineLength = 257;
00265   const int lineLength = 81;
00266   const int tokenLength = 35;
00267   char line[headerlineLength];
00268   char token1[tokenLength];
00269   char token2[tokenLength];
00270   char token3[tokenLength];
00271   char token4[tokenLength];
00272   char token5[tokenLength];
00273   int M, N, NZ;      // Matrix dimensions
00274   int me = comm.MyPID();
00275 
00276   Epetra_Time timer(comm);
00277 
00278   // Make sure domain and range maps are either both defined or undefined
00279   if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
00280     EPETRA_CHK_ERR(-3);
00281   }
00282 
00283   // check maps to see if domain and range are 1-to-1
00284 
00285   if (domainMap!=0) {
00286     if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00287     if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00288   }
00289   else {
00290     // If domain and range are not specified,
00291     // rowMap becomes both and must be 1-to-1 if specified
00292     if (rowMap!=0) {
00293       if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00294     }
00295   }
00296 
00297   FILE * handle = 0;
00298   if (me == 0) {
00299     if (verbose) std::cout << "Reading MatrixMarket file " << filename << std::endl;
00300     handle = fopen(filename,"r");  // Open file
00301     if (handle == 0)
00302       EPETRA_CHK_ERR(-1); // file not found
00303 
00304     // Check first line, which should be
00305     // %%MatrixMarket matrix coordinate real general
00306     if(fgets(line, headerlineLength, handle)==0) {
00307       if (handle!=0) fclose(handle);
00308       EPETRA_CHK_ERR(-1);
00309     }
00310     if(sscanf(line, "%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
00311       if (handle!=0) fclose(handle);
00312       EPETRA_CHK_ERR(-1);
00313     }
00314     if (strcmp(token1, "%%MatrixMarket") ||
00315         strcmp(token2, "matrix") ||
00316         strcmp(token3, "coordinate") ||
00317         strcmp(token4, "real") ||
00318         strcmp(token5, "general")) {
00319       if (handle!=0) fclose(handle);
00320       EPETRA_CHK_ERR(-1);
00321     }
00322 
00323     // Next, strip off header lines (which start with "%")
00324     do {
00325       if(fgets(line, headerlineLength, handle)==0) {
00326         if (handle!=0) fclose(handle);
00327         EPETRA_CHK_ERR(-1);
00328       }
00329     } while (line[0] == '%');
00330 
00331     // Next get problem dimensions: M, N, NZ
00332     if(sscanf(line, "%d %d %d", &M, &N, &NZ)==0) {
00333       if (handle!=0) fclose(handle);
00334       EPETRA_CHK_ERR(-1);
00335     }
00336   }
00337   comm.Broadcast(&M, 1, 0);
00338   comm.Broadcast(&N, 1, 0);
00339   comm.Broadcast(&NZ, 1, 0);
00340 
00341   // Now create matrix using user maps if provided.
00342 
00343 
00344   // Now read in chunks of triplets and broadcast them to other processors.
00345   char *buffer = new char[chunk_read*lineLength];
00346   int nchunk;
00347   int nmillion = 0;
00348   int nread = 0;
00349   int rlen;
00350 
00351   // Storage for this processor's nonzeros.
00352   const int localblock = 100000;
00353   int localsize = NZ / comm.NumProc() + localblock;
00354   int *iv = (int *) malloc(localsize * sizeof(int));
00355   int *jv = (int *) malloc(localsize * sizeof(int));
00356   double *vv = (double *) malloc(localsize * sizeof(double));
00357   int lnz = 0;   //  Number of non-zeros on this processor.
00358 
00359   if (!iv || !jv || !vv)
00360     EPETRA_CHK_ERR(-1);
00361 
00362   Epetra_Map *rowMap1;
00363   bool allocatedHere=false;
00364   if (rowMap != 0)
00365     rowMap1 = const_cast<Epetra_Map *>(rowMap);
00366   else {
00367     rowMap1 = new Epetra_Map(M, 0, comm);
00368     allocatedHere = true;
00369   }
00370   int ioffset = rowMap1->IndexBase()-1;
00371   int joffset = (colMap != 0 ? colMap->IndexBase()-1 : ioffset);
00372 
00373   int rowmajor = 1;  // Assume non-zeros are listed in row-major order,
00374   int prevrow = -1;  // but test to detect otherwise.  If non-zeros
00375                      // are row-major, we can avoid the sort.
00376 
00377   while (nread < NZ) {
00378     if (NZ-nread > chunk_read) nchunk = chunk_read;
00379     else nchunk = NZ - nread;
00380 
00381     if (me == 0) {
00382       char *eof;
00383       rlen = 0;
00384       for (int i = 0; i < nchunk; i++) {
00385         eof = fgets(&buffer[rlen],lineLength,handle);
00386         if (eof == NULL) {
00387           fprintf(stderr, "%s", "Unexpected end of matrix file.");
00388           EPETRA_CHK_ERR(-1);
00389         }
00390         rlen += strlen(&buffer[rlen]);
00391       }
00392       buffer[rlen++]='\n';
00393     }
00394     comm.Broadcast(&rlen, 1, 0);
00395     comm.Broadcast(buffer, rlen, 0);
00396 
00397     buffer[rlen++] = '\0';
00398     nread += nchunk;
00399 
00400     // Process the received data, saving non-zeros belonging on this proc.
00401     char *lineptr = buffer;
00402     for (rlen = 0; rlen < nchunk; rlen++) {
00403       char *next = strchr(lineptr,'\n');
00404       int I = atoi(strtok(lineptr," \t\n")) + ioffset;
00405       int J = atoi(strtok(NULL," \t\n")) + joffset;
00406       double V = atof(strtok(NULL," \t\n"));
00407       lineptr = next + 1;
00408       if (transpose) {
00409         // swap I and J indices.
00410         int tmp = I;
00411         I = J;
00412         J = tmp;
00413       }
00414       if (rowMap1->MyGID(I)) {
00415         //  This processor keeps this non-zero.
00416         if (lnz >= localsize) {
00417           // Need more memory to store nonzeros.
00418           localsize += localblock;
00419           iv = (int *) realloc(iv, localsize * sizeof(int));
00420           jv = (int *) realloc(jv, localsize * sizeof(int));
00421           vv = (double *) realloc(vv, localsize * sizeof(double));
00422         }
00423         iv[lnz] = I;
00424         jv[lnz] = J;
00425         vv[lnz] = V;
00426         lnz++;
00427         if (I < prevrow) rowmajor = 0;
00428         prevrow = I;
00429       }
00430     }
00431 
00432     // Status check
00433     if (nread / 1000000 > nmillion) {
00434       nmillion++;
00435       if (verbose && me == 0) std::cout << nmillion << "M ";
00436     }
00437   }
00438 
00439   delete [] buffer;
00440 
00441   if (!rowmajor) {
00442     // Sort into row-major order (by iv) so can insert entire rows at once.
00443     // Reorder jv and vv to parallel iv.
00444     if (verbose && me == 0) std::cout << std::endl << "   Sorting local nonzeros" << std::endl;
00445     sort_three(iv, jv, vv, 0, lnz-1);
00446   }
00447 
00448   //  Compute number of entries per local row for use in constructor;
00449   //  saves reallocs in FillComplete.
00450 
00451   //  Now construct the matrix.
00452   //  Count number of entries in each row so can use StaticProfile=2.
00453   if (verbose && me == 0) std::cout << std::endl << "   Constructing the matrix" << std::endl;
00454   int numRows = rowMap1->NumMyElements();
00455   int *numNonzerosPerRow = new int[numRows];
00456   for (int i = 0; i < numRows; i++) numNonzerosPerRow[i] = 0;
00457   for (int i = 0; i < lnz; i++)
00458     numNonzerosPerRow[rowMap1->LID(iv[i])]++;
00459 
00460   if (rowMap!=0 && colMap !=0)
00461     A = new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
00462   else if (rowMap!=0) {
00463     // Construct with StaticProfile=true since we know numNonzerosPerRow.
00464     // Less memory will be needed in FillComplete.
00465     A = new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow, true);
00466   }
00467   else {
00468     // Construct with StaticProfile=true since we know numNonzerosPerRow.
00469     // Less memory will be needed in FillComplete.
00470     A = new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow, true);
00471   }
00472   A->SetTracebackMode(1);
00473 
00474   // Rows are inserted in ascending global number, and the insertion uses numNonzerosPerRow.
00475   // Therefore numNonzerosPerRow must be permuted, as it was created in ascending local order.
00476   int *gidList = new int[numRows];
00477   for (int i = 0; i < numRows; i++) gidList[i] = rowMap1->GID(i);
00478   Epetra_Util::Sort(true,numRows,gidList,0,0,1,&numNonzerosPerRow);
00479   delete [] gidList;
00480 
00481   //  Insert the global values into the matrix row-by-row.
00482   if (verbose && me == 0) std::cout << "   Inserting global values" << std::endl;
00483   {
00484     int i = 0;
00485     for (int sum = 0; i < numRows; i++) {
00486       if (numNonzerosPerRow[i]) {
00487         int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
00488                                          &vv[sum], &jv[sum]);
00489         if (ierr<0) EPETRA_CHK_ERR(ierr);
00490         sum += numNonzerosPerRow[i];
00491       }
00492     }
00493   }
00494 
00495   delete [] numNonzerosPerRow;
00496   free(iv);
00497   free(jv);
00498   free(vv);
00499 
00500   if (verbose && me == 0) std::cout << "   Completing matrix fill" << std::endl;
00501   if (rangeMap != 0 && domainMap != 0) {
00502     EPETRA_CHK_ERR(A->FillComplete(*domainMap, *rangeMap));
00503   }
00504   else if (M!=N) {
00505     Epetra_Map newDomainMap(N,rowMap1->IndexBase(), comm);
00506     EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
00507   }
00508   else {
00509     EPETRA_CHK_ERR(A->FillComplete());
00510   }
00511 
00512   if (allocatedHere) delete rowMap1;
00513 
00514   if (handle!=0) fclose(handle);
00515   double dt = timer.ElapsedTime();
00516   if (verbose && me == 0) std::cout << "File Read time (secs):  " << dt << std::endl;
00517   return(0);
00518 }
00519 #endif
00520 
00521 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00522 int MatrixMarketFileToCrsMatrixHandle64(const char *filename,
00523   const Epetra_Comm & comm,
00524   Epetra_CrsMatrix * & A,
00525   const Epetra_Map * rowMap,
00526   const Epetra_Map * colMap,
00527   const Epetra_Map * rangeMap,
00528   const Epetra_Map * domainMap,
00529   const bool transpose,
00530   const bool verbose
00531 )
00532 {
00533   const int chunk_read = 500000;  //  Modify this variable to change the
00534                                   //  size of the chunks read from the file.
00535   const int headerlineLength = 257;
00536   const int lineLength = 81;
00537   const int tokenLength = 35;
00538   char line[headerlineLength];
00539   char token1[tokenLength];
00540   char token2[tokenLength];
00541   char token3[tokenLength];
00542   char token4[tokenLength];
00543   char token5[tokenLength];
00544   long long M, N, NZ;      // Matrix dimensions
00545   int me = comm.MyPID();
00546 
00547   Epetra_Time timer(comm);
00548 
00549   // Make sure domain and range maps are either both defined or undefined
00550   if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
00551     EPETRA_CHK_ERR(-3);
00552   }
00553 
00554   // check maps to see if domain and range are 1-to-1
00555 
00556   if (domainMap!=0) {
00557     if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00558     if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00559   }
00560   else {
00561     // If domain and range are not specified,
00562     // rowMap becomes both and must be 1-to-1 if specified
00563     if (rowMap!=0) {
00564       if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
00565     }
00566   }
00567 
00568   FILE * handle = 0;
00569   if (me == 0) {
00570     if (verbose) std::cout << "Reading MatrixMarket file " << filename << std::endl;
00571     handle = fopen(filename,"r");  // Open file
00572     if (handle == 0)
00573       EPETRA_CHK_ERR(-1); // file not found
00574 
00575     // Check first line, which should be
00576     // %%MatrixMarket matrix coordinate real general
00577     if(fgets(line, headerlineLength, handle)==0) {
00578       if (handle!=0) fclose(handle);
00579       EPETRA_CHK_ERR(-1);
00580     }
00581     if(sscanf(line, "%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
00582       if (handle!=0) fclose(handle);
00583       EPETRA_CHK_ERR(-1);
00584     }
00585     if (strcmp(token1, "%%MatrixMarket") ||
00586         strcmp(token2, "matrix") ||
00587         strcmp(token3, "coordinate") ||
00588         strcmp(token4, "real") ||
00589         strcmp(token5, "general")) {
00590       if (handle!=0) fclose(handle);
00591       EPETRA_CHK_ERR(-1);
00592     }
00593 
00594     // Next, strip off header lines (which start with "%")
00595     do {
00596       if(fgets(line, headerlineLength, handle)==0) {
00597         if (handle!=0) fclose(handle);
00598         EPETRA_CHK_ERR(-1);
00599       }
00600     } while (line[0] == '%');
00601 
00602     // Next get problem dimensions: M, N, NZ
00603     if(sscanf(line, "%lld %lld %lld", &M, &N, &NZ)==0) {
00604       if (handle!=0) fclose(handle);
00605       EPETRA_CHK_ERR(-1);
00606     }
00607   }
00608   comm.Broadcast(&M, 1, 0);
00609   comm.Broadcast(&N, 1, 0);
00610   comm.Broadcast(&NZ, 1, 0);
00611 
00612   // Now create matrix using user maps if provided.
00613 
00614 
00615   // Now read in chunks of triplets and broadcast them to other processors.
00616   char *buffer = new char[chunk_read*lineLength];
00617   long long nchunk;
00618   int nmillion = 0;
00619   long long nread = 0;
00620   int rlen;
00621 
00622   // Storage for this processor's nonzeros.
00623   const int localblock = 100000;
00624   int localsize = (int) (NZ / comm.NumProc()) + localblock;
00625   long long *iv = (long long *) malloc(localsize * sizeof(long long));
00626   long long *jv = (long long *) malloc(localsize * sizeof(long long));
00627   double *vv = (double *) malloc(localsize * sizeof(double));
00628   int lnz = 0;   //  Number of non-zeros on this processor.
00629 
00630   if (!iv || !jv || !vv)
00631     EPETRA_CHK_ERR(-1);
00632 
00633   Epetra_Map *rowMap1;
00634   bool allocatedHere=false;
00635   if (rowMap != 0)
00636     rowMap1 = const_cast<Epetra_Map *>(rowMap);
00637   else {
00638     rowMap1 = new Epetra_Map(M, 0, comm);
00639     allocatedHere = true;
00640   }
00641   long long ioffset = rowMap1->IndexBase64()-1;
00642   long long joffset = (colMap != 0 ? colMap->IndexBase64()-1 : ioffset);
00643 
00644   int rowmajor = 1;  // Assume non-zeros are listed in row-major order,
00645   long long prevrow = -1;  // but test to detect otherwise.  If non-zeros
00646                      // are row-major, we can avoid the sort.
00647 
00648   while (nread < NZ) {
00649     if (NZ-nread > chunk_read) nchunk = chunk_read;
00650     else nchunk = NZ - nread;
00651 
00652     if (me == 0) {
00653       char *eof;
00654       rlen = 0;
00655       for (int i = 0; i < nchunk; i++) {
00656         eof = fgets(&buffer[rlen],lineLength,handle);
00657         if (eof == NULL) {
00658           fprintf(stderr, "%s", "Unexpected end of matrix file.");
00659           EPETRA_CHK_ERR(-1);
00660         }
00661         rlen += strlen(&buffer[rlen]);
00662       }
00663       buffer[rlen++]='\n';
00664     }
00665     comm.Broadcast(&rlen, 1, 0);
00666     comm.Broadcast(buffer, rlen, 0);
00667 
00668     buffer[rlen++] = '\0';
00669     nread += nchunk;
00670 
00671     // Process the received data, saving non-zeros belonging on this proc.
00672     char *lineptr = buffer;
00673     for (rlen = 0; rlen < nchunk; rlen++) {
00674       char *next = strchr(lineptr,'\n');
00675       char *endp;
00676       const int base = 10;
00677 #if defined(_MSC_VER)
00678       long long I = _strtoi64(strtok(lineptr," \t\n"), &endp, base) + ioffset;
00679       long long J = _strtoi64(strtok(NULL," \t\n"), &endp, base) + joffset;
00680 #else
00681 #if defined(__PGI)
00682       long long I, J;
00683       std::istringstream ssI(strtok(lineptr," \t\n"));
00684       ssI >> I; I += ioffset;
00685       std::istringstream ssJ(strtok(NULL," \t\n"));
00686       ssJ >> J; J += joffset;
00687 #else
00688       long long I = strtoll(strtok(lineptr," \t\n"), &endp, base) + ioffset;
00689       long long J = strtoll(strtok(NULL," \t\n"), &endp, base) + joffset;
00690 #endif
00691 #endif
00692       double V = atof(strtok(NULL," \t\n"));
00693       lineptr = next + 1;
00694       if (transpose) {
00695         // swap I and J indices.
00696         long long tmp = I;
00697         I = J;
00698         J = tmp;
00699       }
00700       if (rowMap1->MyGID(I)) {
00701         //  This processor keeps this non-zero.
00702         if (lnz >= localsize) {
00703           // Need more memory to store nonzeros.
00704           localsize += localblock;
00705           iv = (long long *) realloc(iv, localsize * sizeof(long long));
00706           jv = (long long *) realloc(jv, localsize * sizeof(long long));
00707           vv = (double *) realloc(vv, localsize * sizeof(double));
00708         }
00709         iv[lnz] = I;
00710         jv[lnz] = J;
00711         vv[lnz] = V;
00712         lnz++;
00713         if (I < prevrow) rowmajor = 0;
00714         prevrow = I;
00715       }
00716     }
00717 
00718     // Status check
00719     if (nread / 1000000 > nmillion) {
00720       nmillion++;
00721       if (verbose && me == 0) std::cout << nmillion << "M ";
00722     }
00723   }
00724 
00725   delete [] buffer;
00726 
00727   if (!rowmajor) {
00728     // Sort into row-major order (by iv) so can insert entire rows at once.
00729     // Reorder jv and vv to parallel iv.
00730     if (verbose && me == 0) std::cout << std::endl << "   Sorting local nonzeros" << std::endl;
00731     sort_three(iv, jv, vv, 0, lnz-1);
00732   }
00733 
00734   //  Compute number of entries per local row for use in constructor;
00735   //  saves reallocs in FillComplete.
00736 
00737   //  Now construct the matrix.
00738   //  Count number of entries in each row so can use StaticProfile=2.
00739   if (verbose && me == 0) std::cout << std::endl << "   Constructing the matrix" << std::endl;
00740   int numRows = rowMap1->NumMyElements();
00741   int *numNonzerosPerRow = new int[numRows];
00742   for (int i = 0; i < numRows; i++) numNonzerosPerRow[i] = 0;
00743   for (int i = 0; i < lnz; i++)
00744     numNonzerosPerRow[rowMap1->LID(iv[i])]++;
00745 
00746   if (rowMap!=0 && colMap !=0)
00747     A = new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
00748   else if (rowMap!=0) {
00749     // Construct with StaticProfile=true since we know numNonzerosPerRow.
00750     // Less memory will be needed in FillComplete.
00751     A = new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow, true);
00752   }
00753   else {
00754     // Construct with StaticProfile=true since we know numNonzerosPerRow.
00755     // Less memory will be needed in FillComplete.
00756     A = new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow, true);
00757   }
00758   A->SetTracebackMode(1);
00759 
00760   // Rows are inserted in ascending global number, and the insertion uses numNonzerosPerRow.
00761   // Therefore numNonzerosPerRow must be permuted, as it was created in ascending local order.
00762   long long *gidList = new long long[numRows];
00763   for (int i = 0; i < numRows; i++) gidList[i] = rowMap1->GID64(i);
00764   Epetra_Util::Sort(true,numRows,gidList,0,0,1,&numNonzerosPerRow,0,0);
00765   delete [] gidList;
00766 
00767   //  Insert the global values into the matrix row-by-row.
00768   if (verbose && me == 0) std::cout << "   Inserting global values" << std::endl;
00769   {
00770     int i = 0;
00771     for (int sum = 0; i < numRows; i++) {
00772       if (numNonzerosPerRow[i]) {
00773         int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
00774                                          &vv[sum], &jv[sum]);
00775         if (ierr<0) EPETRA_CHK_ERR(ierr);
00776         sum += numNonzerosPerRow[i];
00777       }
00778     }
00779   }
00780 
00781   delete [] numNonzerosPerRow;
00782   free(iv);
00783   free(jv);
00784   free(vv);
00785 
00786   if (verbose && me == 0) std::cout << "   Completing matrix fill" << std::endl;
00787   if (rangeMap != 0 && domainMap != 0) {
00788     EPETRA_CHK_ERR(A->FillComplete(*domainMap, *rangeMap));
00789   }
00790   else if (M!=N) {
00791     Epetra_Map newDomainMap(N,rowMap1->IndexBase64(), comm);
00792     EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
00793   }
00794   else {
00795     EPETRA_CHK_ERR(A->FillComplete());
00796   }
00797 
00798   if (allocatedHere) delete rowMap1;
00799 
00800   if (handle!=0) fclose(handle);
00801   double dt = timer.ElapsedTime();
00802   if (verbose && me == 0) std::cout << "File Read time (secs):  " << dt << std::endl;
00803   return(0);
00804 }
00805 #endif
00806 
00808 // Sorting values in array list in increasing order. Criteria is int.
00809 // Also rearrange values in arrays parlista and partlistb
00810 // to match the new order of list.
00811 
00812 template<typename int_type>
00813 static void quickpart_list_inc_int (
00814   int_type *list, int_type *parlista, double *parlistb,
00815   int start, int end, int *equal, int *larger)
00816 {
00817 int i;
00818 int_type key, itmp;
00819 double dtmp;
00820 
00821   key = list ? list[(end+start)/2] : 1;
00822 
00823   *equal = *larger = start;
00824   for (i = start; i <= end; i++)
00825     if (list[i] < key) {
00826       itmp                = parlista[i];
00827       parlista[i]         = parlista[*larger];
00828       parlista[(*larger)] = parlista[*equal];
00829       parlista[(*equal)]  = itmp;
00830       dtmp                = parlistb[i];
00831       parlistb[i]         = parlistb[*larger];
00832       parlistb[(*larger)] = parlistb[*equal];
00833       parlistb[(*equal)]  = dtmp;
00834       itmp                = list[i];
00835       list[i]             = list[*larger];
00836       list[(*larger)++]   = list[*equal];
00837       list[(*equal)++]    = itmp;
00838     }
00839     else if (list[i] == key) {
00840       itmp                = parlista[i];
00841       parlista[i]         = parlista[*larger];
00842       parlista[(*larger)] = itmp;
00843       dtmp                = parlistb[i];
00844       parlistb[i]         = parlistb[*larger];
00845       parlistb[(*larger)] = dtmp;
00846       list[i]             = list[*larger];
00847       list[(*larger)++]   = key;
00848     }
00849 }
00850 
00852 template<typename int_type>
00853 static void sort_three(int_type* list, int_type *parlista, double *parlistb,
00854                        int start, int end)
00855 {
00856 int  equal, larger;
00857 
00858   if (start < end) {
00859     quickpart_list_inc_int(list, parlista, parlistb, start, end,
00860                            &equal, &larger);
00861     sort_three(list, parlista, parlistb, start,  equal-1);
00862     sort_three(list, parlista, parlistb, larger, end);
00863   }
00864 }
00865 
00867 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00868 int MatlabFileToCrsMatrix(const char *filename,
00869                                 const Epetra_Comm & comm,
00870                                 Epetra_CrsMatrix * & A)
00871 {
00872   const int lineLength = 1025;
00873   char line[lineLength];
00874   int I, J;
00875   double V;
00876 
00877   FILE * handle = 0;
00878 
00879   handle = fopen(filename,"r");  // Open file
00880   if (handle == 0)
00881     EPETRA_CHK_ERR(-1); // file not found
00882 
00883   int numGlobalRows = 0;
00884   int numGlobalCols = 0;
00885   while(fgets(line, lineLength, handle)!=0) {
00886     if(sscanf(line, "%d %d %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
00887     if (I>numGlobalRows) numGlobalRows = I;
00888     if (J>numGlobalCols) numGlobalCols = J;
00889   }
00890 
00891   if (handle!=0) fclose(handle);
00892   Epetra_Map rangeMap(numGlobalRows, 0, comm);
00893   Epetra_Map domainMap(numGlobalCols, 0, comm);
00894   A = new Epetra_CrsMatrix(Copy, rangeMap, 0);
00895 
00896   // Now read in each triplet and store to the local portion of the matrix if the row is owned.
00897   const Epetra_Map & rowMap1 = A->RowMap();
00898 
00899   handle = 0;
00900 
00901   handle = fopen(filename,"r");  // Open file
00902   if (handle == 0)
00903     EPETRA_CHK_ERR(-1); // file not found
00904 
00905   while (fgets(line, lineLength, handle)!=0) {
00906     if(sscanf(line, "%d %d %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
00907     I--; J--; // Convert to Zero based
00908     if (rowMap1.MyGID(I)) {
00909       int ierr = A->InsertGlobalValues(I, 1, &V, &J);
00910       if (ierr<0) EPETRA_CHK_ERR(ierr);
00911     }
00912   }
00913 
00914   EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
00915 
00916   if (handle!=0) fclose(handle);
00917   return(0);
00918 }
00919 #endif
00920 
00921 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00922 int MatlabFileToCrsMatrix64(const char *filename,
00923                                 const Epetra_Comm & comm,
00924                                 Epetra_CrsMatrix * & A)
00925 {
00926   const int lineLength = 1025;
00927   char line[lineLength];
00928   long long I, J;
00929   double V;
00930 
00931   FILE * handle = 0;
00932 
00933   handle = fopen(filename,"r");  // Open file
00934   if (handle == 0)
00935     EPETRA_CHK_ERR(-1); // file not found
00936 
00937   long long numGlobalRows = 0;
00938   long long numGlobalCols = 0;
00939   while(fgets(line, lineLength, handle)!=0) {
00940     if(sscanf(line, "%lld %lld %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
00941     if (I>numGlobalRows) numGlobalRows = I;
00942     if (J>numGlobalCols) numGlobalCols = J;
00943   }
00944 
00945   if (handle!=0) fclose(handle);
00946   Epetra_Map rangeMap(numGlobalRows, 0, comm);
00947   Epetra_Map domainMap(numGlobalCols, 0, comm);
00948   A = new Epetra_CrsMatrix(Copy, rangeMap, 0);
00949 
00950   // Now read in each triplet and store to the local portion of the matrix if the row is owned.
00951   const Epetra_Map & rowMap1 = A->RowMap();
00952 
00953   handle = 0;
00954 
00955   handle = fopen(filename,"r");  // Open file
00956   if (handle == 0)
00957     EPETRA_CHK_ERR(-1); // file not found
00958 
00959   while (fgets(line, lineLength, handle)!=0) {
00960     if(sscanf(line, "%lld %lld %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
00961     I--; J--; // Convert to Zero based
00962     if (rowMap1.MyGID(I)) {
00963       int ierr = A->InsertGlobalValues(I, 1, &V, &J);
00964       if (ierr<0) EPETRA_CHK_ERR(ierr);
00965     }
00966   }
00967 
00968   EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
00969 
00970   if (handle!=0) fclose(handle);
00971   return(0);
00972 }
00973 #endif
00974 
00975 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00976 int HypreFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix){
00977   int MyPID = comm.MyPID();
00978   // This double will be in the format we want for the extension besides the leading zero
00979   double filePID = (double)MyPID/(double)100000;
00980   std::ostringstream stream;
00981   // Using setprecision() puts it in the std::string
00982   stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
00983   // Then just ignore the first character
00984   std::string fileName(filename);
00985   fileName += stream.str().substr(1,7);
00986   // Open the file
00987   std::ifstream file(fileName.c_str());
00988   std::string line;
00989   if(file.is_open()){
00990     std::getline(file, line);
00991     int ilower, iupper;
00992     std::istringstream istream(line);
00993     // The first line of the file has the beginning and ending rows
00994     istream >> ilower;
00995     istream >> iupper;
00996     // Using those we can create a row map
00997     Epetra_Map RowMap(-1, iupper-ilower+1, 0, comm);
00998     Matrix = new Epetra_CrsMatrix(Copy, RowMap, 0);
00999     int currRow = -1;
01000     int counter = 0;
01001     std::vector<int> indices;
01002     std::vector<double> values;
01003     while(!file.eof()){
01004       std::getline(file, line);
01005       std::istringstream lineStr(line);
01006       int row, col;
01007       double val;
01008       lineStr >> row;
01009       lineStr >> col;
01010       lineStr >> val;
01011       if(currRow == -1) currRow = row; // First line
01012       if(row == currRow){
01013         // add to the vector
01014         counter = counter + 1;
01015         indices.push_back(col);
01016         values.push_back(val);
01017       } else {
01018         Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
01019         indices.clear();
01020         values.clear();
01021         counter = 0;
01022         currRow = row;
01023         // make a new vector
01024         indices.push_back(col);
01025         values.push_back(val);
01026         counter = counter + 1;
01027       }
01028     }
01029     Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
01030     Matrix->Comm().Barrier();
01031     Matrix->FillComplete();
01032     file.close();
01033     return 0;
01034   } else {
01035     std::cout << "\nERROR:\nCouldn't open " << fileName << ".\n";
01036     return -1;
01037   }
01038 }
01039 #endif
01040 
01041 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
01042 int HypreFileToCrsMatrix64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix){
01043   int MyPID = comm.MyPID();
01044   // This double will be in the format we want for the extension besides the leading zero
01045   double filePID = (double)MyPID/(double)100000;
01046   std::ostringstream stream;
01047   // Using setprecision() puts it in the std::string
01048   stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
01049   // Then just ignore the first character
01050   std::string fileName(filename);
01051   fileName += stream.str().substr(1,7);
01052   // Open the file
01053   std::ifstream file(fileName.c_str());
01054   std::string line;
01055   if(file.is_open()){
01056     std::getline(file, line);
01057     int ilower, iupper;
01058     std::istringstream istream(line);
01059     // The first line of the file has the beginning and ending rows
01060     istream >> ilower;
01061     istream >> iupper;
01062     // Using those we can create a row map
01063     Epetra_Map RowMap(-1LL, iupper-ilower+1, 0LL, comm);
01064     Matrix = new Epetra_CrsMatrix(Copy, RowMap, 0);
01065     long long currRow = -1;
01066     int counter = 0;
01067     std::vector<long long> indices;
01068     std::vector<double> values;
01069     while(!file.eof()){
01070       std::getline(file, line);
01071       std::istringstream lineStr(line);
01072       long long row, col;
01073       double val;
01074       lineStr >> row;
01075       lineStr >> col;
01076       lineStr >> val;
01077       if(currRow == -1) currRow = row; // First line
01078       if(row == currRow){
01079         // add to the vector
01080         counter = counter + 1;
01081         indices.push_back(col);
01082         values.push_back(val);
01083       } else {
01084         Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
01085         indices.clear();
01086         values.clear();
01087         counter = 0;
01088         currRow = row;
01089         // make a new vector
01090         indices.push_back(col);
01091         values.push_back(val);
01092         counter = counter + 1;
01093       }
01094     }
01095     Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
01096     Matrix->Comm().Barrier();
01097     Matrix->FillComplete();
01098     file.close();
01099     return 0;
01100   } else {
01101     std::cout << "\nERROR:\nCouldn't open " << fileName << ".\n";
01102     return -1;
01103   }
01104 }
01105 #endif
01106 
01107 } // namespace EpetraExt
01108 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines