EpetraExt  Development
EpetraExt_CrsMatrixIn.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 #include "Epetra_ConfigDefs.h"
42 #include "EpetraExt_CrsMatrixIn.h"
43 #include "Epetra_Comm.h"
44 #include "Epetra_CrsMatrix.h"
45 #include "Epetra_Map.h"
46 #include "Epetra_IntVector.h"
47 #include "Epetra_IntSerialDenseVector.h"
48 #include "Epetra_Import.h"
49 #include "Epetra_Time.h"
50 #include "Epetra_Util.h"
51 
52 #if defined(__PGI)
53 #include <sstream>
54 #endif
55 
58 // Functions to read a MatrixMarket file and load it into a matrix.
59 // Adapted from
60 // Trilinos/packages/epetraext/src/inout/EpetraExt_CrsMatrixIn.cpp
61 // Modified by Jon Berry and Karen Devine to make matrix reallocations
62 // more efficient; rather than insert each non-zero separately, we
63 // collect rows of non-zeros for insertion.
64 // Modified by Karen Devine and Steve Plimpton to prevent all processors
65 // from reading the same file at the same time (ouch!); chunks of the file
66 // are read and broadcast by processor zero; each processor then extracts
67 // its nonzeros from the chunk, sorts them by row, and inserts them into
68 // the matrix.
69 // The variable "chunk_read" can be changed to modify the size of the chunks
70 // read from the file. Larger values of chunk_read lead to faster execution
71 // and greater memory use.
74 
75 using namespace EpetraExt;
76 namespace EpetraExt {
77 
78 template<typename int_type>
79 static void sort_three(int_type* list, int_type *parlista, double *parlistb,
80  int start, int end);
81 
83 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
84 int MatrixMarketFileToCrsMatrix(const char *filename,
85  const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
86 {
87  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A));
88  return(0);
89 }
90 
92 int MatrixMarketFileToCrsMatrix(const char *filename,
93  const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose)
94 {
95  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A,
96  0, 0, 0, 0, transpose));
97  return(0);
98 }
100 
101 int MatrixMarketFileToCrsMatrix(const char *filename,
102  const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose,
103  const bool verbose)
104 {
105  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename, comm, A,
106  0, 0, 0, 0,
107  transpose, verbose));
108  return(0);
109 }
110 
112 int MatrixMarketFileToCrsMatrix(const char *filename,
113  const Epetra_Map & rowMap, const Epetra_Map& rangeMap,
114  const Epetra_Map& domainMap, Epetra_CrsMatrix * & A, const bool transpose,
115  const bool verbose)
116 {
117  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
118  rowMap.Comm(), A,
119  &rowMap, 0,
120  &rangeMap, &domainMap,
121  transpose, verbose));
122  return(0);
123 }
124 
126 int MatrixMarketFileToCrsMatrix(const char *filename,
127  const Epetra_Map & rowMap, Epetra_CrsMatrix * & A, const bool transpose,
128  const bool verbose)
129 {
130  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
131  rowMap.Comm(), A,
132  &rowMap, 0, 0, 0,
133  transpose, verbose));
134  return(0);
135 }
136 
138 int MatrixMarketFileToCrsMatrix(const char *filename,
139  const Epetra_Map & rowMap, const Epetra_Map & colMap,
140  Epetra_CrsMatrix * & A, const bool transpose,
141  const bool verbose)
142 {
143  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
144  rowMap.Comm(), A,
145  &rowMap, &colMap, 0, 0,
146  transpose, verbose));
147  return(0);
148 }
149 
151 int MatrixMarketFileToCrsMatrix(const char *filename,
152  const Epetra_Map & rowMap, const Epetra_Map & colMap,
153  const Epetra_Map& rangeMap, const Epetra_Map& domainMap,
154  Epetra_CrsMatrix * & A, const bool transpose,
155  const bool verbose)
156 {
157  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle(filename,
158  rowMap.Comm(), A,
159  &rowMap, &colMap,
160  &rangeMap, &domainMap,
161  transpose, verbose));
162  return(0);
163 }
164 #endif
165 
166 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
167 int MatrixMarketFileToCrsMatrix64(const char *filename,
168  const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
169 {
170  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A));
171  return(0);
172 }
173 
175 int MatrixMarketFileToCrsMatrix64(const char *filename,
176  const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose)
177 {
178  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A,
179  0, 0, 0, 0, transpose));
180  return(0);
181 }
183 
184 int MatrixMarketFileToCrsMatrix64(const char *filename,
185  const Epetra_Comm & comm, Epetra_CrsMatrix * & A, const bool transpose,
186  const bool verbose)
187 {
188  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename, comm, A,
189  0, 0, 0, 0,
190  transpose, verbose));
191  return(0);
192 }
193 
195 int MatrixMarketFileToCrsMatrix64(const char *filename,
196  const Epetra_Map & rowMap, const Epetra_Map& rangeMap,
197  const Epetra_Map& domainMap, Epetra_CrsMatrix * & A, const bool transpose,
198  const bool verbose)
199 {
200  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
201  rowMap.Comm(), A,
202  &rowMap, 0,
203  &rangeMap, &domainMap,
204  transpose, verbose));
205  return(0);
206 }
207 
209 int MatrixMarketFileToCrsMatrix64(const char *filename,
210  const Epetra_Map & rowMap, Epetra_CrsMatrix * & A, const bool transpose,
211  const bool verbose)
212 {
213  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
214  rowMap.Comm(), A,
215  &rowMap, 0, 0, 0,
216  transpose, verbose));
217  return(0);
218 }
219 
221 int MatrixMarketFileToCrsMatrix64(const char *filename,
222  const Epetra_Map & rowMap, const Epetra_Map & colMap,
223  Epetra_CrsMatrix * & A, const bool transpose,
224  const bool verbose)
225 {
226  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
227  rowMap.Comm(), A,
228  &rowMap, &colMap, 0, 0,
229  transpose, verbose));
230  return(0);
231 }
232 
234 int MatrixMarketFileToCrsMatrix64(const char *filename,
235  const Epetra_Map & rowMap, const Epetra_Map & colMap,
236  const Epetra_Map& rangeMap, const Epetra_Map& domainMap,
237  Epetra_CrsMatrix * & A, const bool transpose,
238  const bool verbose)
239 {
240  EPETRA_CHK_ERR(MatrixMarketFileToCrsMatrixHandle64(filename,
241  rowMap.Comm(), A,
242  &rowMap, &colMap,
243  &rangeMap, &domainMap,
244  transpose, verbose));
245  return(0);
246 }
247 #endif
248 
250 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
251 int MatrixMarketFileToCrsMatrixHandle(const char *filename,
252  const Epetra_Comm & comm,
253  Epetra_CrsMatrix * & A,
254  const Epetra_Map * rowMap,
255  const Epetra_Map * colMap,
256  const Epetra_Map * rangeMap,
257  const Epetra_Map * domainMap,
258  const bool transpose,
259  const bool verbose
260 )
261 {
262  const int chunk_read = 500000; // Modify this variable to change the
263  // size of the chunks read from the file.
264  const int headerlineLength = 257;
265  const int lineLength = 81;
266  const int tokenLength = 35;
267  char line[headerlineLength];
268  char token1[tokenLength];
269  char token2[tokenLength];
270  char token3[tokenLength];
271  char token4[tokenLength];
272  char token5[tokenLength];
273  int M, N, NZ; // Matrix dimensions
274  int me = comm.MyPID();
275 
276  Epetra_Time timer(comm);
277 
278  // Make sure domain and range maps are either both defined or undefined
279  if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
280  EPETRA_CHK_ERR(-3);
281  }
282 
283  // check maps to see if domain and range are 1-to-1
284 
285  if (domainMap!=0) {
286  if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
287  if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
288  }
289  else {
290  // If domain and range are not specified,
291  // rowMap becomes both and must be 1-to-1 if specified
292  if (rowMap!=0) {
293  if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
294  }
295  }
296 
297  FILE * handle = 0;
298  int rv=0;
299  if (me == 0) {
300  if (verbose) std::cout << "Reading MatrixMarket file " << filename << std::endl;
301  handle = fopen(filename,"r"); // Open file
302  if (handle == 0) {
303  rv = -1; // file not found
304  }
305 
306  // Check first line, which should be
307  // %%MatrixMarket matrix coordinate real general
308  if( (rv != -1) && fgets(line, headerlineLength, handle)==0 ) {
309  if (handle!=0) fclose(handle);
310  rv = -1;
311  }
312  if( (rv != -1) && sscanf(line, "%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
313  if (handle!=0) fclose(handle);
314  rv = -1;
315  }
316  if ( (rv != -1) && (strcmp(token1, "%%MatrixMarket") ||
317  strcmp(token2, "matrix") ||
318  strcmp(token3, "coordinate") ||
319  strcmp(token4, "real") ||
320  strcmp(token5, "general")) ) {
321  if (handle!=0) fclose(handle);
322  rv = -1;
323  }
324 
325  // Next, strip off header lines (which start with "%")
326  if (rv != -1) {
327  do {
328  if(fgets(line, headerlineLength, handle)==0) {
329  if (handle!=0) fclose(handle);
330  rv = -1;
331  break;
332  }
333  } while (line[0] == '%');
334  }
335 
336  // Next get problem dimensions: M, N, NZ
337  if((rv != -1) && sscanf(line, "%d %d %d", &M, &N, &NZ)==0) {
338  if (handle!=0) fclose(handle);
339  rv = -1;
340  }
341  }
342  // If there's an error reading the file, broadcast to all processes in communicator.
343  comm.Broadcast(&rv, 1, 0);
344  if (rv==-1)
345  EPETRA_CHK_ERR(-1);
346  comm.Broadcast(&M, 1, 0);
347  comm.Broadcast(&N, 1, 0);
348  comm.Broadcast(&NZ, 1, 0);
349 
350  // Now create matrix using user maps if provided.
351 
352 
353  // Now read in chunks of triplets and broadcast them to other processors.
354  char *buffer = new char[chunk_read*lineLength];
355  int nchunk;
356  int nmillion = 0;
357  int nread = 0;
358  int rlen;
359 
360  // Storage for this processor's nonzeros.
361  const int localblock = 100000;
362  int localsize = NZ / comm.NumProc() + localblock;
363  int *iv = (int *) malloc(localsize * sizeof(int));
364  int *jv = (int *) malloc(localsize * sizeof(int));
365  double *vv = (double *) malloc(localsize * sizeof(double));
366  int lnz = 0; // Number of non-zeros on this processor.
367 
368  if (!iv || !jv || !vv)
369  EPETRA_CHK_ERR(-1);
370 
371  Epetra_Map *rowMap1;
372  bool allocatedHere=false;
373  if (rowMap != 0)
374  rowMap1 = const_cast<Epetra_Map *>(rowMap);
375  else {
376  rowMap1 = new Epetra_Map(M, 0, comm);
377  allocatedHere = true;
378  }
379  int ioffset = rowMap1->IndexBase()-1;
380  int joffset = (colMap != 0 ? colMap->IndexBase()-1 : ioffset);
381 
382  int rowmajor = 1; // Assume non-zeros are listed in row-major order,
383  int prevrow = -1; // but test to detect otherwise. If non-zeros
384  // are row-major, we can avoid the sort.
385 
386  while (nread < NZ) {
387  if (NZ-nread > chunk_read) nchunk = chunk_read;
388  else nchunk = NZ - nread;
389 
390  if (me == 0) {
391  char *eof;
392  rlen = 0;
393  for (int i = 0; i < nchunk; i++) {
394  eof = fgets(&buffer[rlen],lineLength,handle);
395  if (eof == NULL) {
396  fprintf(stderr, "%s", "Unexpected end of matrix file.");
397  EPETRA_CHK_ERR(-1);
398  }
399  rlen += strlen(&buffer[rlen]);
400  }
401  buffer[rlen++]='\n';
402  }
403  comm.Broadcast(&rlen, 1, 0);
404  comm.Broadcast(buffer, rlen, 0);
405 
406  buffer[rlen++] = '\0';
407  nread += nchunk;
408 
409  // Process the received data, saving non-zeros belonging on this proc.
410  char *lineptr = buffer;
411  for (rlen = 0; rlen < nchunk; rlen++) {
412  char *next = strchr(lineptr,'\n');
413  int I = atoi(strtok(lineptr," \t\n")) + ioffset;
414  int J = atoi(strtok(NULL," \t\n")) + joffset;
415  double V = atof(strtok(NULL," \t\n"));
416  lineptr = next + 1;
417  if (transpose) {
418  // swap I and J indices.
419  int tmp = I;
420  I = J;
421  J = tmp;
422  }
423  if (rowMap1->MyGID(I)) {
424  // This processor keeps this non-zero.
425  if (lnz >= localsize) {
426  // Need more memory to store nonzeros.
427  localsize += localblock;
428  iv = (int *) realloc(iv, localsize * sizeof(int));
429  jv = (int *) realloc(jv, localsize * sizeof(int));
430  vv = (double *) realloc(vv, localsize * sizeof(double));
431  }
432  iv[lnz] = I;
433  jv[lnz] = J;
434  vv[lnz] = V;
435  lnz++;
436  if (I < prevrow) rowmajor = 0;
437  prevrow = I;
438  }
439  }
440 
441  // Status check
442  if (nread / 1000000 > nmillion) {
443  nmillion++;
444  if (verbose && me == 0) std::cout << nmillion << "M ";
445  }
446  }
447 
448  delete [] buffer;
449 
450  if (!rowmajor) {
451  // Sort into row-major order (by iv) so can insert entire rows at once.
452  // Reorder jv and vv to parallel iv.
453  if (verbose && me == 0) std::cout << std::endl << " Sorting local nonzeros" << std::endl;
454  sort_three(iv, jv, vv, 0, lnz-1);
455  }
456 
457  // Compute number of entries per local row for use in constructor;
458  // saves reallocs in FillComplete.
459 
460  // Now construct the matrix.
461  // Count number of entries in each row so can use StaticProfile=2.
462  if (verbose && me == 0) std::cout << std::endl << " Constructing the matrix" << std::endl;
463  int numRows = rowMap1->NumMyElements();
464  int *numNonzerosPerRow = new int[numRows];
465  for (int i = 0; i < numRows; i++) numNonzerosPerRow[i] = 0;
466  for (int i = 0; i < lnz; i++)
467  numNonzerosPerRow[rowMap1->LID(iv[i])]++;
468 
469  if (rowMap!=0 && colMap !=0)
470  A = new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
471  else if (rowMap!=0) {
472  // Construct with StaticProfile=true since we know numNonzerosPerRow.
473  // Less memory will be needed in FillComplete.
474  A = new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow, true);
475  }
476  else {
477  // Construct with StaticProfile=true since we know numNonzerosPerRow.
478  // Less memory will be needed in FillComplete.
479  A = new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow, true);
480  }
481  A->SetTracebackMode(1);
482 
483  // Rows are inserted in ascending global number, and the insertion uses numNonzerosPerRow.
484  // Therefore numNonzerosPerRow must be permuted, as it was created in ascending local order.
485  int *gidList = new int[numRows];
486  for (int i = 0; i < numRows; i++) gidList[i] = rowMap1->GID(i);
487  Epetra_Util::Sort(true,numRows,gidList,0,0,1,&numNonzerosPerRow);
488  delete [] gidList;
489 
490  // Insert the global values into the matrix row-by-row.
491  if (verbose && me == 0) std::cout << " Inserting global values" << std::endl;
492  {
493  int i = 0;
494  for (int sum = 0; i < numRows; i++) {
495  if (numNonzerosPerRow[i]) {
496  int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
497  &vv[sum], &jv[sum]);
498  if (ierr<0) EPETRA_CHK_ERR(ierr);
499  sum += numNonzerosPerRow[i];
500  }
501  }
502  }
503 
504  delete [] numNonzerosPerRow;
505  free(iv);
506  free(jv);
507  free(vv);
508 
509  if (verbose && me == 0) std::cout << " Completing matrix fill" << std::endl;
510  if (rangeMap != 0 && domainMap != 0) {
511  EPETRA_CHK_ERR(A->FillComplete(*domainMap, *rangeMap));
512  }
513  else if (M!=N) {
514  Epetra_Map newDomainMap(N,rowMap1->IndexBase(), comm);
515  EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
516  }
517  else {
518  EPETRA_CHK_ERR(A->FillComplete());
519  }
520 
521  if (allocatedHere) delete rowMap1;
522 
523  if (handle!=0) fclose(handle);
524  double dt = timer.ElapsedTime();
525  if (verbose && me == 0) std::cout << "File Read time (secs): " << dt << std::endl;
526  return(0);
527 }
528 #endif
529 
530 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
531 int MatrixMarketFileToCrsMatrixHandle64(const char *filename,
532  const Epetra_Comm & comm,
533  Epetra_CrsMatrix * & A,
534  const Epetra_Map * rowMap,
535  const Epetra_Map * colMap,
536  const Epetra_Map * rangeMap,
537  const Epetra_Map * domainMap,
538  const bool transpose,
539  const bool verbose
540 )
541 {
542  const int chunk_read = 500000; // Modify this variable to change the
543  // size of the chunks read from the file.
544  const int headerlineLength = 257;
545  const int lineLength = 81;
546  const int tokenLength = 35;
547  char line[headerlineLength];
548  char token1[tokenLength];
549  char token2[tokenLength];
550  char token3[tokenLength];
551  char token4[tokenLength];
552  char token5[tokenLength];
553  long long M, N, NZ; // Matrix dimensions
554  int me = comm.MyPID();
555 
556  Epetra_Time timer(comm);
557 
558  // Make sure domain and range maps are either both defined or undefined
559  if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
560  EPETRA_CHK_ERR(-3);
561  }
562 
563  // check maps to see if domain and range are 1-to-1
564 
565  if (domainMap!=0) {
566  if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
567  if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
568  }
569  else {
570  // If domain and range are not specified,
571  // rowMap becomes both and must be 1-to-1 if specified
572  if (rowMap!=0) {
573  if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
574  }
575  }
576 
577  FILE * handle = 0;
578  if (me == 0) {
579  if (verbose) std::cout << "Reading MatrixMarket file " << filename << std::endl;
580  handle = fopen(filename,"r"); // Open file
581  if (handle == 0)
582  EPETRA_CHK_ERR(-1); // file not found
583 
584  // Check first line, which should be
585  // %%MatrixMarket matrix coordinate real general
586  if(fgets(line, headerlineLength, handle)==0) {
587  if (handle!=0) fclose(handle);
588  EPETRA_CHK_ERR(-1);
589  }
590  if(sscanf(line, "%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
591  if (handle!=0) fclose(handle);
592  EPETRA_CHK_ERR(-1);
593  }
594  if (strcmp(token1, "%%MatrixMarket") ||
595  strcmp(token2, "matrix") ||
596  strcmp(token3, "coordinate") ||
597  strcmp(token4, "real") ||
598  strcmp(token5, "general")) {
599  if (handle!=0) fclose(handle);
600  EPETRA_CHK_ERR(-1);
601  }
602 
603  // Next, strip off header lines (which start with "%")
604  do {
605  if(fgets(line, headerlineLength, handle)==0) {
606  if (handle!=0) fclose(handle);
607  EPETRA_CHK_ERR(-1);
608  }
609  } while (line[0] == '%');
610 
611  // Next get problem dimensions: M, N, NZ
612  if(sscanf(line, "%lld %lld %lld", &M, &N, &NZ)==0) {
613  if (handle!=0) fclose(handle);
614  EPETRA_CHK_ERR(-1);
615  }
616  }
617  comm.Broadcast(&M, 1, 0);
618  comm.Broadcast(&N, 1, 0);
619  comm.Broadcast(&NZ, 1, 0);
620 
621  // Now create matrix using user maps if provided.
622 
623 
624  // Now read in chunks of triplets and broadcast them to other processors.
625  char *buffer = new char[chunk_read*lineLength];
626  long long nchunk;
627  int nmillion = 0;
628  long long nread = 0;
629  int rlen;
630 
631  // Storage for this processor's nonzeros.
632  const int localblock = 100000;
633  int localsize = (int) (NZ / comm.NumProc()) + localblock;
634  long long *iv = (long long *) malloc(localsize * sizeof(long long));
635  long long *jv = (long long *) malloc(localsize * sizeof(long long));
636  double *vv = (double *) malloc(localsize * sizeof(double));
637  int lnz = 0; // Number of non-zeros on this processor.
638 
639  if (!iv || !jv || !vv)
640  EPETRA_CHK_ERR(-1);
641 
642  Epetra_Map *rowMap1;
643  bool allocatedHere=false;
644  if (rowMap != 0)
645  rowMap1 = const_cast<Epetra_Map *>(rowMap);
646  else {
647  rowMap1 = new Epetra_Map(M, 0, comm);
648  allocatedHere = true;
649  }
650  long long ioffset = rowMap1->IndexBase64()-1;
651  long long joffset = (colMap != 0 ? colMap->IndexBase64()-1 : ioffset);
652 
653  int rowmajor = 1; // Assume non-zeros are listed in row-major order,
654  long long prevrow = -1; // but test to detect otherwise. If non-zeros
655  // are row-major, we can avoid the sort.
656 
657  while (nread < NZ) {
658  if (NZ-nread > chunk_read) nchunk = chunk_read;
659  else nchunk = NZ - nread;
660 
661  if (me == 0) {
662  char *eof;
663  rlen = 0;
664  for (int i = 0; i < nchunk; i++) {
665  eof = fgets(&buffer[rlen],lineLength,handle);
666  if (eof == NULL) {
667  fprintf(stderr, "%s", "Unexpected end of matrix file.");
668  EPETRA_CHK_ERR(-1);
669  }
670  rlen += strlen(&buffer[rlen]);
671  }
672  buffer[rlen++]='\n';
673  }
674  comm.Broadcast(&rlen, 1, 0);
675  comm.Broadcast(buffer, rlen, 0);
676 
677  buffer[rlen++] = '\0';
678  nread += nchunk;
679 
680  // Process the received data, saving non-zeros belonging on this proc.
681  char *lineptr = buffer;
682  for (rlen = 0; rlen < nchunk; rlen++) {
683  char *next = strchr(lineptr,'\n');
684  char *endp;
685  const int base = 10;
686 #if defined(_MSC_VER)
687  long long I = _strtoi64(strtok(lineptr," \t\n"), &endp, base) + ioffset;
688  long long J = _strtoi64(strtok(NULL," \t\n"), &endp, base) + joffset;
689 #else
690 #if defined(__PGI)
691  long long I, J;
692  std::istringstream ssI(strtok(lineptr," \t\n"));
693  ssI >> I; I += ioffset;
694  std::istringstream ssJ(strtok(NULL," \t\n"));
695  ssJ >> J; J += joffset;
696 #else
697  long long I = strtoll(strtok(lineptr," \t\n"), &endp, base) + ioffset;
698  long long J = strtoll(strtok(NULL," \t\n"), &endp, base) + joffset;
699 #endif
700 #endif
701  double V = atof(strtok(NULL," \t\n"));
702  lineptr = next + 1;
703  if (transpose) {
704  // swap I and J indices.
705  long long tmp = I;
706  I = J;
707  J = tmp;
708  }
709  if (rowMap1->MyGID(I)) {
710  // This processor keeps this non-zero.
711  if (lnz >= localsize) {
712  // Need more memory to store nonzeros.
713  localsize += localblock;
714  iv = (long long *) realloc(iv, localsize * sizeof(long long));
715  jv = (long long *) realloc(jv, localsize * sizeof(long long));
716  vv = (double *) realloc(vv, localsize * sizeof(double));
717  }
718  iv[lnz] = I;
719  jv[lnz] = J;
720  vv[lnz] = V;
721  lnz++;
722  if (I < prevrow) rowmajor = 0;
723  prevrow = I;
724  }
725  }
726 
727  // Status check
728  if (nread / 1000000 > nmillion) {
729  nmillion++;
730  if (verbose && me == 0) std::cout << nmillion << "M ";
731  }
732  }
733 
734  delete [] buffer;
735 
736  if (!rowmajor) {
737  // Sort into row-major order (by iv) so can insert entire rows at once.
738  // Reorder jv and vv to parallel iv.
739  if (verbose && me == 0) std::cout << std::endl << " Sorting local nonzeros" << std::endl;
740  sort_three(iv, jv, vv, 0, lnz-1);
741  }
742 
743  // Compute number of entries per local row for use in constructor;
744  // saves reallocs in FillComplete.
745 
746  // Now construct the matrix.
747  // Count number of entries in each row so can use StaticProfile=2.
748  if (verbose && me == 0) std::cout << std::endl << " Constructing the matrix" << std::endl;
749  int numRows = rowMap1->NumMyElements();
750  int *numNonzerosPerRow = new int[numRows];
751  for (int i = 0; i < numRows; i++) numNonzerosPerRow[i] = 0;
752  for (int i = 0; i < lnz; i++)
753  numNonzerosPerRow[rowMap1->LID(iv[i])]++;
754 
755  if (rowMap!=0 && colMap !=0)
756  A = new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
757  else if (rowMap!=0) {
758  // Construct with StaticProfile=true since we know numNonzerosPerRow.
759  // Less memory will be needed in FillComplete.
760  A = new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow, true);
761  }
762  else {
763  // Construct with StaticProfile=true since we know numNonzerosPerRow.
764  // Less memory will be needed in FillComplete.
765  A = new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow, true);
766  }
767  A->SetTracebackMode(1);
768 
769  // Rows are inserted in ascending global number, and the insertion uses numNonzerosPerRow.
770  // Therefore numNonzerosPerRow must be permuted, as it was created in ascending local order.
771  long long *gidList = new long long[numRows];
772  for (int i = 0; i < numRows; i++) gidList[i] = rowMap1->GID64(i);
773  Epetra_Util::Sort(true,numRows,gidList,0,0,1,&numNonzerosPerRow,0,0);
774  delete [] gidList;
775 
776  // Insert the global values into the matrix row-by-row.
777  if (verbose && me == 0) std::cout << " Inserting global values" << std::endl;
778  {
779  int i = 0;
780  for (int sum = 0; i < numRows; i++) {
781  if (numNonzerosPerRow[i]) {
782  int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
783  &vv[sum], &jv[sum]);
784  if (ierr<0) EPETRA_CHK_ERR(ierr);
785  sum += numNonzerosPerRow[i];
786  }
787  }
788  }
789 
790  delete [] numNonzerosPerRow;
791  free(iv);
792  free(jv);
793  free(vv);
794 
795  if (verbose && me == 0) std::cout << " Completing matrix fill" << std::endl;
796  if (rangeMap != 0 && domainMap != 0) {
797  EPETRA_CHK_ERR(A->FillComplete(*domainMap, *rangeMap));
798  }
799  else if (M!=N) {
800  Epetra_Map newDomainMap(N,rowMap1->IndexBase64(), comm);
801  EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
802  }
803  else {
804  EPETRA_CHK_ERR(A->FillComplete());
805  }
806 
807  if (allocatedHere) delete rowMap1;
808 
809  if (handle!=0) fclose(handle);
810  double dt = timer.ElapsedTime();
811  if (verbose && me == 0) std::cout << "File Read time (secs): " << dt << std::endl;
812  return(0);
813 }
814 #endif
815 
817 // Sorting values in array list in increasing order. Criteria is int.
818 // Also rearrange values in arrays parlista and partlistb
819 // to match the new order of list.
820 
821 template<typename int_type>
823  int_type *list, int_type *parlista, double *parlistb,
824  int start, int end, int *equal, int *larger)
825 {
826 int i;
827 int_type key, itmp;
828 double dtmp;
829 
830  key = list ? list[(end+start)/2] : 1;
831 
832  *equal = *larger = start;
833  for (i = start; i <= end; i++)
834  if (list[i] < key) {
835  itmp = parlista[i];
836  parlista[i] = parlista[*larger];
837  parlista[(*larger)] = parlista[*equal];
838  parlista[(*equal)] = itmp;
839  dtmp = parlistb[i];
840  parlistb[i] = parlistb[*larger];
841  parlistb[(*larger)] = parlistb[*equal];
842  parlistb[(*equal)] = dtmp;
843  itmp = list[i];
844  list[i] = list[*larger];
845  list[(*larger)++] = list[*equal];
846  list[(*equal)++] = itmp;
847  }
848  else if (list[i] == key) {
849  itmp = parlista[i];
850  parlista[i] = parlista[*larger];
851  parlista[(*larger)] = itmp;
852  dtmp = parlistb[i];
853  parlistb[i] = parlistb[*larger];
854  parlistb[(*larger)] = dtmp;
855  list[i] = list[*larger];
856  list[(*larger)++] = key;
857  }
858 }
859 
861 template<typename int_type>
862 static void sort_three(int_type* list, int_type *parlista, double *parlistb,
863  int start, int end)
864 {
865 int equal, larger;
866 
867  if (start < end) {
868  quickpart_list_inc_int(list, parlista, parlistb, start, end,
869  &equal, &larger);
870  sort_three(list, parlista, parlistb, start, equal-1);
871  sort_three(list, parlista, parlistb, larger, end);
872  }
873 }
874 
876 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
877 int MatlabFileToCrsMatrix(const char *filename,
878  const Epetra_Comm & comm,
879  Epetra_CrsMatrix * & A)
880 {
881  const int lineLength = 1025;
882  char line[lineLength];
883  int I, J;
884  double V;
885 
886  FILE * handle = 0;
887 
888  handle = fopen(filename,"r"); // Open file
889  if (handle == 0)
890  EPETRA_CHK_ERR(-1); // file not found
891 
892  int numGlobalRows = 0;
893  int numGlobalCols = 0;
894  while(fgets(line, lineLength, handle)!=0) {
895  if(sscanf(line, "%d %d %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
896  if (I>numGlobalRows) numGlobalRows = I;
897  if (J>numGlobalCols) numGlobalCols = J;
898  }
899 
900  if (handle!=0) fclose(handle);
901  Epetra_Map rangeMap(numGlobalRows, 0, comm);
902  Epetra_Map domainMap(numGlobalCols, 0, comm);
903  A = new Epetra_CrsMatrix(Copy, rangeMap, 0);
904 
905  // Now read in each triplet and store to the local portion of the matrix if the row is owned.
906  const Epetra_Map & rowMap1 = A->RowMap();
907 
908  handle = 0;
909 
910  handle = fopen(filename,"r"); // Open file
911  if (handle == 0)
912  EPETRA_CHK_ERR(-1); // file not found
913 
914  while (fgets(line, lineLength, handle)!=0) {
915  if(sscanf(line, "%d %d %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
916  I--; J--; // Convert to Zero based
917  if (rowMap1.MyGID(I)) {
918  int ierr = A->InsertGlobalValues(I, 1, &V, &J);
919  if (ierr<0) EPETRA_CHK_ERR(ierr);
920  }
921  }
922 
923  EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
924 
925  if (handle!=0) fclose(handle);
926  return(0);
927 }
928 #endif
929 
930 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
931 int MatlabFileToCrsMatrix64(const char *filename,
932  const Epetra_Comm & comm,
933  Epetra_CrsMatrix * & A)
934 {
935  const int lineLength = 1025;
936  char line[lineLength];
937  long long I, J;
938  double V;
939 
940  FILE * handle = 0;
941 
942  handle = fopen(filename,"r"); // Open file
943  if (handle == 0)
944  EPETRA_CHK_ERR(-1); // file not found
945 
946  long long numGlobalRows = 0;
947  long long numGlobalCols = 0;
948  while(fgets(line, lineLength, handle)!=0) {
949  if(sscanf(line, "%lld %lld %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
950  if (I>numGlobalRows) numGlobalRows = I;
951  if (J>numGlobalCols) numGlobalCols = J;
952  }
953 
954  if (handle!=0) fclose(handle);
955  Epetra_Map rangeMap(numGlobalRows, 0, comm);
956  Epetra_Map domainMap(numGlobalCols, 0, comm);
957  A = new Epetra_CrsMatrix(Copy, rangeMap, 0);
958 
959  // Now read in each triplet and store to the local portion of the matrix if the row is owned.
960  const Epetra_Map & rowMap1 = A->RowMap();
961 
962  handle = 0;
963 
964  handle = fopen(filename,"r"); // Open file
965  if (handle == 0)
966  EPETRA_CHK_ERR(-1); // file not found
967 
968  while (fgets(line, lineLength, handle)!=0) {
969  if(sscanf(line, "%lld %lld %lg\n", &I, &J, &V)==0) {if (handle!=0) fclose(handle); EPETRA_CHK_ERR(-1);}
970  I--; J--; // Convert to Zero based
971  if (rowMap1.MyGID(I)) {
972  int ierr = A->InsertGlobalValues(I, 1, &V, &J);
973  if (ierr<0) EPETRA_CHK_ERR(ierr);
974  }
975  }
976 
977  EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
978 
979  if (handle!=0) fclose(handle);
980  return(0);
981 }
982 #endif
983 
984 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
985 int HypreFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix){
986  int MyPID = comm.MyPID();
987  // This double will be in the format we want for the extension besides the leading zero
988  double filePID = (double)MyPID/(double)100000;
989  std::ostringstream stream;
990  // Using setprecision() puts it in the std::string
991  stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
992  // Then just ignore the first character
993  std::string fileName(filename);
994  fileName += stream.str().substr(1,7);
995  // Open the file
996  std::ifstream file(fileName.c_str());
997  std::string line;
998  if(file.is_open()){
999  std::getline(file, line);
1000  int ilower, iupper;
1001  std::istringstream istream(line);
1002  // The first line of the file has the beginning and ending rows
1003  istream >> ilower;
1004  istream >> iupper;
1005  // Using those we can create a row map
1006  Epetra_Map RowMap(-1, iupper-ilower+1, 0, comm);
1007  Matrix = new Epetra_CrsMatrix(Copy, RowMap, 0);
1008  int currRow = -1;
1009  int counter = 0;
1010  std::vector<int> indices;
1011  std::vector<double> values;
1012  while(!file.eof()){
1013  std::getline(file, line);
1014  std::istringstream lineStr(line);
1015  int row, col;
1016  double val;
1017  lineStr >> row;
1018  lineStr >> col;
1019  lineStr >> val;
1020  if(currRow == -1) currRow = row; // First line
1021  if(row == currRow){
1022  // add to the vector
1023  counter = counter + 1;
1024  indices.push_back(col);
1025  values.push_back(val);
1026  } else {
1027  Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1028  indices.clear();
1029  values.clear();
1030  counter = 0;
1031  currRow = row;
1032  // make a new vector
1033  indices.push_back(col);
1034  values.push_back(val);
1035  counter = counter + 1;
1036  }
1037  }
1038  Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1039  Matrix->Comm().Barrier();
1040  Matrix->FillComplete();
1041  file.close();
1042  return 0;
1043  } else {
1044  std::cout << "\nERROR:\nCouldn't open " << fileName << ".\n";
1045  return -1;
1046  }
1047 }
1048 #endif
1049 
1050 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
1051 int HypreFileToCrsMatrix64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix){
1052  int MyPID = comm.MyPID();
1053  // This double will be in the format we want for the extension besides the leading zero
1054  double filePID = (double)MyPID/(double)100000;
1055  std::ostringstream stream;
1056  // Using setprecision() puts it in the std::string
1057  stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
1058  // Then just ignore the first character
1059  std::string fileName(filename);
1060  fileName += stream.str().substr(1,7);
1061  // Open the file
1062  std::ifstream file(fileName.c_str());
1063  std::string line;
1064  if(file.is_open()){
1065  std::getline(file, line);
1066  int ilower, iupper;
1067  std::istringstream istream(line);
1068  // The first line of the file has the beginning and ending rows
1069  istream >> ilower;
1070  istream >> iupper;
1071  // Using those we can create a row map
1072  Epetra_Map RowMap(-1LL, iupper-ilower+1, 0LL, comm);
1073  Matrix = new Epetra_CrsMatrix(Copy, RowMap, 0);
1074  long long currRow = -1;
1075  int counter = 0;
1076  std::vector<long long> indices;
1077  std::vector<double> values;
1078  while(!file.eof()){
1079  std::getline(file, line);
1080  std::istringstream lineStr(line);
1081  long long row, col;
1082  double val;
1083  lineStr >> row;
1084  lineStr >> col;
1085  lineStr >> val;
1086  if(currRow == -1) currRow = row; // First line
1087  if(row == currRow){
1088  // add to the vector
1089  counter = counter + 1;
1090  indices.push_back(col);
1091  values.push_back(val);
1092  } else {
1093  Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1094  indices.clear();
1095  values.clear();
1096  counter = 0;
1097  currRow = row;
1098  // make a new vector
1099  indices.push_back(col);
1100  values.push_back(val);
1101  counter = counter + 1;
1102  }
1103  }
1104  Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1105  Matrix->Comm().Barrier();
1106  Matrix->FillComplete();
1107  file.close();
1108  return 0;
1109  } else {
1110  std::cout << "\nERROR:\nCouldn't open " << fileName << ".\n";
1111  return -1;
1112  }
1113 }
1114 #endif
1115 
1116 } // namespace EpetraExt
1117 
int MatrixMarketFileToCrsMatrix64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
int HypreFileToCrsMatrix64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix)
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
static void quickpart_list_inc_int(int_type *list, int_type *parlista, double *parlistb, int start, int end, int *equal, int *larger)
int MatrixMarketFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
int MatlabFileToCrsMatrix64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
int MatrixMarketFileToCrsMatrixHandle(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A, const Epetra_Map *rowMap, const Epetra_Map *colMap, const Epetra_Map *rangeMap, const Epetra_Map *domainMap, const bool transpose, const bool verbose)
int MatlabFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
Constructs an Epetra_CrsMatrix object from a Matlab format file, distributes rows evenly across proce...
int HypreFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix)
Constructs an Epetra_CrsMatrix object from a Hypre Matrix Print command, the row map is specified...
int MatrixMarketFileToCrsMatrixHandle64(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A, const Epetra_Map *rowMap, const Epetra_Map *colMap, const Epetra_Map *rangeMap, const Epetra_Map *domainMap, const bool transpose, const bool verbose)
static void sort_three(int_type *list, int_type *parlista, double *parlistb, int start, int end)