41 #include "Epetra_ConfigDefs.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" 78 template<
typename int_type>
79 static void sort_three(int_type* list, int_type *parlista,
double *parlistb,
83 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 85 const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
93 const Epetra_Comm & comm, Epetra_CrsMatrix * & A,
const bool transpose)
96 0, 0, 0, 0, transpose));
102 const Epetra_Comm & comm, Epetra_CrsMatrix * & A,
const bool transpose,
107 transpose, verbose));
113 const Epetra_Map & rowMap,
const Epetra_Map& rangeMap,
114 const Epetra_Map& domainMap, Epetra_CrsMatrix * & A,
const bool transpose,
120 &rangeMap, &domainMap,
121 transpose, verbose));
127 const Epetra_Map & rowMap, Epetra_CrsMatrix * & A,
const bool transpose,
133 transpose, verbose));
139 const Epetra_Map & rowMap,
const Epetra_Map & colMap,
140 Epetra_CrsMatrix * & A,
const bool transpose,
145 &rowMap, &colMap, 0, 0,
146 transpose, verbose));
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,
160 &rangeMap, &domainMap,
161 transpose, verbose));
166 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES 168 const Epetra_Comm & comm, Epetra_CrsMatrix * & A)
176 const Epetra_Comm & comm, Epetra_CrsMatrix * & A,
const bool transpose)
179 0, 0, 0, 0, transpose));
185 const Epetra_Comm & comm, Epetra_CrsMatrix * & A,
const bool transpose,
190 transpose, verbose));
196 const Epetra_Map & rowMap,
const Epetra_Map& rangeMap,
197 const Epetra_Map& domainMap, Epetra_CrsMatrix * & A,
const bool transpose,
203 &rangeMap, &domainMap,
204 transpose, verbose));
210 const Epetra_Map & rowMap, Epetra_CrsMatrix * & A,
const bool transpose,
216 transpose, verbose));
222 const Epetra_Map & rowMap,
const Epetra_Map & colMap,
223 Epetra_CrsMatrix * & A,
const bool transpose,
228 &rowMap, &colMap, 0, 0,
229 transpose, verbose));
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,
243 &rangeMap, &domainMap,
244 transpose, verbose));
250 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 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,
262 const int chunk_read = 500000;
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];
274 int me = comm.MyPID();
276 Epetra_Time timer(comm);
279 if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
286 if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
287 if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
293 if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
300 if (verbose) std::cout <<
"Reading MatrixMarket file " << filename << std::endl;
301 handle = fopen(filename,
"r");
308 if( (rv != -1) && fgets(line, headerlineLength, handle)==0 ) {
309 if (handle!=0) fclose(handle);
312 if( (rv != -1) && sscanf(line,
"%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
313 if (handle!=0) fclose(handle);
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);
328 if(fgets(line, headerlineLength, handle)==0) {
329 if (handle!=0) fclose(handle);
333 }
while (line[0] ==
'%');
337 if((rv != -1) && sscanf(line,
"%d %d %d", &M, &N, &NZ)==0) {
338 if (handle!=0) fclose(handle);
343 comm.Broadcast(&rv, 1, 0);
346 comm.Broadcast(&M, 1, 0);
347 comm.Broadcast(&N, 1, 0);
348 comm.Broadcast(&NZ, 1, 0);
354 char *buffer =
new char[chunk_read*lineLength];
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));
368 if (!iv || !jv || !vv)
372 bool allocatedHere=
false;
374 rowMap1 =
const_cast<Epetra_Map *
>(rowMap);
376 rowMap1 =
new Epetra_Map(M, 0, comm);
377 allocatedHere =
true;
379 int ioffset = rowMap1->IndexBase()-1;
380 int joffset = (colMap != 0 ? colMap->IndexBase()-1 : ioffset);
387 if (NZ-nread > chunk_read) nchunk = chunk_read;
388 else nchunk = NZ - nread;
393 for (
int i = 0; i < nchunk; i++) {
394 eof = fgets(&buffer[rlen],lineLength,handle);
396 fprintf(stderr,
"%s",
"Unexpected end of matrix file.");
399 rlen += strlen(&buffer[rlen]);
403 comm.Broadcast(&rlen, 1, 0);
404 comm.Broadcast(buffer, rlen, 0);
406 buffer[rlen++] =
'\0';
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"));
423 if (rowMap1->MyGID(I)) {
425 if (lnz >= localsize) {
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));
436 if (I < prevrow) rowmajor = 0;
442 if (nread / 1000000 > nmillion) {
444 if (verbose && me == 0) std::cout << nmillion <<
"M ";
453 if (verbose && me == 0) std::cout << std::endl <<
" Sorting local nonzeros" << std::endl;
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])]++;
469 if (rowMap!=0 && colMap !=0)
470 A =
new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
471 else if (rowMap!=0) {
474 A =
new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow,
true);
479 A =
new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow,
true);
481 A->SetTracebackMode(1);
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);
491 if (verbose && me == 0) std::cout <<
" Inserting global values" << std::endl;
494 for (
int sum = 0; i < numRows; i++) {
495 if (numNonzerosPerRow[i]) {
496 int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
498 if (ierr<0) EPETRA_CHK_ERR(ierr);
499 sum += numNonzerosPerRow[i];
504 delete [] numNonzerosPerRow;
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));
514 Epetra_Map newDomainMap(N,rowMap1->IndexBase(), comm);
515 EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
518 EPETRA_CHK_ERR(A->FillComplete());
521 if (allocatedHere)
delete rowMap1;
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;
530 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES 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,
542 const int chunk_read = 500000;
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];
554 int me = comm.MyPID();
556 Epetra_Time timer(comm);
559 if ((domainMap!=0 && rangeMap==0) || (domainMap==0 && rangeMap!=0)) {
566 if (!domainMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
567 if (!rangeMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
573 if (!rowMap->UniqueGIDs()) {EPETRA_CHK_ERR(-2);}
579 if (verbose) std::cout <<
"Reading MatrixMarket file " << filename << std::endl;
580 handle = fopen(filename,
"r");
586 if(fgets(line, headerlineLength, handle)==0) {
587 if (handle!=0) fclose(handle);
590 if(sscanf(line,
"%s %s %s %s %s", token1,token2,token3,token4,token5)==0) {
591 if (handle!=0) fclose(handle);
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);
605 if(fgets(line, headerlineLength, handle)==0) {
606 if (handle!=0) fclose(handle);
609 }
while (line[0] ==
'%');
612 if(sscanf(line,
"%lld %lld %lld", &M, &N, &NZ)==0) {
613 if (handle!=0) fclose(handle);
617 comm.Broadcast(&M, 1, 0);
618 comm.Broadcast(&N, 1, 0);
619 comm.Broadcast(&NZ, 1, 0);
625 char *buffer =
new char[chunk_read*lineLength];
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));
639 if (!iv || !jv || !vv)
643 bool allocatedHere=
false;
645 rowMap1 =
const_cast<Epetra_Map *
>(rowMap);
647 rowMap1 =
new Epetra_Map(M, 0, comm);
648 allocatedHere =
true;
650 long long ioffset = rowMap1->IndexBase64()-1;
651 long long joffset = (colMap != 0 ? colMap->IndexBase64()-1 : ioffset);
654 long long prevrow = -1;
658 if (NZ-nread > chunk_read) nchunk = chunk_read;
659 else nchunk = NZ - nread;
664 for (
int i = 0; i < nchunk; i++) {
665 eof = fgets(&buffer[rlen],lineLength,handle);
667 fprintf(stderr,
"%s",
"Unexpected end of matrix file.");
670 rlen += strlen(&buffer[rlen]);
674 comm.Broadcast(&rlen, 1, 0);
675 comm.Broadcast(buffer, rlen, 0);
677 buffer[rlen++] =
'\0';
681 char *lineptr = buffer;
682 for (rlen = 0; rlen < nchunk; rlen++) {
683 char *next = strchr(lineptr,
'\n');
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;
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;
697 long long I = strtoll(strtok(lineptr,
" \t\n"), &endp, base) + ioffset;
698 long long J = strtoll(strtok(NULL,
" \t\n"), &endp, base) + joffset;
701 double V = atof(strtok(NULL,
" \t\n"));
709 if (rowMap1->MyGID(I)) {
711 if (lnz >= localsize) {
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));
722 if (I < prevrow) rowmajor = 0;
728 if (nread / 1000000 > nmillion) {
730 if (verbose && me == 0) std::cout << nmillion <<
"M ";
739 if (verbose && me == 0) std::cout << std::endl <<
" Sorting local nonzeros" << std::endl;
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])]++;
755 if (rowMap!=0 && colMap !=0)
756 A =
new Epetra_CrsMatrix(Copy, *rowMap, *colMap, numNonzerosPerRow);
757 else if (rowMap!=0) {
760 A =
new Epetra_CrsMatrix(Copy, *rowMap, numNonzerosPerRow,
true);
765 A =
new Epetra_CrsMatrix(Copy, *rowMap1, numNonzerosPerRow,
true);
767 A->SetTracebackMode(1);
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);
777 if (verbose && me == 0) std::cout <<
" Inserting global values" << std::endl;
780 for (
int sum = 0; i < numRows; i++) {
781 if (numNonzerosPerRow[i]) {
782 int ierr = A->InsertGlobalValues(iv[sum], numNonzerosPerRow[i],
784 if (ierr<0) EPETRA_CHK_ERR(ierr);
785 sum += numNonzerosPerRow[i];
790 delete [] numNonzerosPerRow;
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));
800 Epetra_Map newDomainMap(N,rowMap1->IndexBase64(), comm);
801 EPETRA_CHK_ERR(A->FillComplete(newDomainMap, *rowMap1));
804 EPETRA_CHK_ERR(A->FillComplete());
807 if (allocatedHere)
delete rowMap1;
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;
821 template<
typename int_type>
823 int_type *list, int_type *parlista,
double *parlistb,
824 int start,
int end,
int *equal,
int *larger)
830 key = list ? list[(end+start)/2] : 1;
832 *equal = *larger = start;
833 for (i = start; i <= end; i++)
836 parlista[i] = parlista[*larger];
837 parlista[(*larger)] = parlista[*equal];
838 parlista[(*equal)] = itmp;
840 parlistb[i] = parlistb[*larger];
841 parlistb[(*larger)] = parlistb[*equal];
842 parlistb[(*equal)] = dtmp;
844 list[i] = list[*larger];
845 list[(*larger)++] = list[*equal];
846 list[(*equal)++] = itmp;
848 else if (list[i] == key) {
850 parlista[i] = parlista[*larger];
851 parlista[(*larger)] = itmp;
853 parlistb[i] = parlistb[*larger];
854 parlistb[(*larger)] = dtmp;
855 list[i] = list[*larger];
856 list[(*larger)++] = key;
861 template<
typename int_type>
862 static void sort_three(int_type* list, int_type *parlista,
double *parlistb,
870 sort_three(list, parlista, parlistb, start, equal-1);
871 sort_three(list, parlista, parlistb, larger, end);
876 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 878 const Epetra_Comm & comm,
879 Epetra_CrsMatrix * & A)
881 const int lineLength = 1025;
882 char line[lineLength];
888 handle = fopen(filename,
"r");
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;
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);
906 const Epetra_Map & rowMap1 = A->RowMap();
910 handle = fopen(filename,
"r");
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);}
917 if (rowMap1.MyGID(I)) {
918 int ierr = A->InsertGlobalValues(I, 1, &V, &J);
919 if (ierr<0) EPETRA_CHK_ERR(ierr);
923 EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
925 if (handle!=0) fclose(handle);
930 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES 932 const Epetra_Comm & comm,
933 Epetra_CrsMatrix * & A)
935 const int lineLength = 1025;
936 char line[lineLength];
942 handle = fopen(filename,
"r");
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;
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);
960 const Epetra_Map & rowMap1 = A->RowMap();
964 handle = fopen(filename,
"r");
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);}
971 if (rowMap1.MyGID(I)) {
972 int ierr = A->InsertGlobalValues(I, 1, &V, &J);
973 if (ierr<0) EPETRA_CHK_ERR(ierr);
977 EPETRA_CHK_ERR(A->FillComplete(domainMap, rangeMap));
979 if (handle!=0) fclose(handle);
984 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 986 int MyPID = comm.MyPID();
988 double filePID = (double)MyPID/(
double)100000;
989 std::ostringstream stream;
991 stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
993 std::string fileName(filename);
994 fileName += stream.str().substr(1,7);
996 std::ifstream file(fileName.c_str());
999 std::getline(file, line);
1001 std::istringstream istream(line);
1006 Epetra_Map RowMap(-1, iupper-ilower+1, 0, comm);
1007 Matrix =
new Epetra_CrsMatrix(Copy, RowMap, 0);
1010 std::vector<int> indices;
1011 std::vector<double> values;
1013 std::getline(file, line);
1014 std::istringstream lineStr(line);
1020 if(currRow == -1) currRow = row;
1023 counter = counter + 1;
1024 indices.push_back(col);
1025 values.push_back(val);
1027 Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1033 indices.push_back(col);
1034 values.push_back(val);
1035 counter = counter + 1;
1038 Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1039 Matrix->Comm().Barrier();
1040 Matrix->FillComplete();
1044 std::cout <<
"\nERROR:\nCouldn't open " << fileName <<
".\n";
1050 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES 1052 int MyPID = comm.MyPID();
1054 double filePID = (double)MyPID/(
double)100000;
1055 std::ostringstream stream;
1057 stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
1059 std::string fileName(filename);
1060 fileName += stream.str().substr(1,7);
1062 std::ifstream file(fileName.c_str());
1065 std::getline(file, line);
1067 std::istringstream istream(line);
1072 Epetra_Map RowMap(-1LL, iupper-ilower+1, 0LL, comm);
1073 Matrix =
new Epetra_CrsMatrix(Copy, RowMap, 0);
1074 long long currRow = -1;
1076 std::vector<long long> indices;
1077 std::vector<double> values;
1079 std::getline(file, line);
1080 std::istringstream lineStr(line);
1086 if(currRow == -1) currRow = row;
1089 counter = counter + 1;
1090 indices.push_back(col);
1091 values.push_back(val);
1093 Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1099 indices.push_back(col);
1100 values.push_back(val);
1101 counter = counter + 1;
1104 Matrix->InsertGlobalValues(currRow, counter, &values[0], &indices[0]);
1105 Matrix->Comm().Barrier();
1106 Matrix->FillComplete();
1110 std::cout <<
"\nERROR:\nCouldn't open " << fileName <<
".\n";
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)