1 #include "mmloader.h" 2 3 PetscErrorCode MatCreateFromMTX(Mat *A, const char *filein, PetscBool aijonly) 4 { 5 MM_typecode matcode; 6 FILE *file; 7 PetscInt M, N, ninput; 8 PetscInt *ia, *ja; 9 PetscInt i, j, nz, *rownz; 10 PetscScalar *val; 11 PetscBool sametype, symmetric = PETSC_FALSE, skew = PETSC_FALSE; 12 13 /* Read in matrix */ 14 PetscCall(PetscFOpen(PETSC_COMM_SELF, filein, "r", &file)); 15 PetscCheck(mm_read_banner(file, &matcode) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not process Matrix Market banner."); 16 /* This is how one can screen matrix types if their application */ 17 /* only supports a subset of the Matrix Market data types. */ 18 PetscCheck(mm_is_matrix(matcode) && mm_is_sparse(matcode) && (mm_is_real(matcode) || mm_is_integer(matcode)), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input must be a sparse real or integer matrix. Market Market type: [%s]", mm_typecode_to_str(matcode)); 19 20 if (mm_is_symmetric(matcode)) symmetric = PETSC_TRUE; 21 if (mm_is_skew(matcode)) skew = PETSC_TRUE; 22 23 /* Find out size of sparse matrix .... */ 24 PetscCheck(mm_read_mtx_crd_size(file, &M, &N, &nz) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Size of sparse matrix is wrong."); 25 26 /* Reserve memory for matrices */ 27 PetscCall(PetscMalloc4(nz, &ia, nz, &ja, nz, &val, M, &rownz)); 28 for (i = 0; i < M; i++) rownz[i] = 0; 29 30 /* NOTE: when reading in doubles, ANSI C requires the use of the "l" */ 31 /* specifier as in "%lg", "%lf", "%le", otherwise errors will occur */ 32 /* (ANSI C X3.159-1989, Sec. 4.9.6.2, p. 136 lines 13-15) */ 33 for (i = 0; i < nz; i++) { 34 ninput = fscanf(file, "%d %d %lg\n", &ia[i], &ja[i], &val[i]); 35 PetscCheck(ninput >= 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Badly formatted input file"); 36 ia[i]--; 37 ja[i]--; /* adjust from 1-based to 0-based */ 38 if ((symmetric && aijonly) || skew) { /* transpose */ 39 rownz[ia[i]]++; 40 rownz[ja[i]]++; 41 } else rownz[ia[i]]++; 42 } 43 PetscCall(PetscFClose(PETSC_COMM_SELF, file)); 44 45 /* Create, preallocate, and then assemble the matrix */ 46 PetscCall(MatCreate(PETSC_COMM_SELF, A)); 47 PetscCall(MatSetSizes(*A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 48 49 if (symmetric && !aijonly) { 50 PetscCall(MatSetType(*A, MATSEQSBAIJ)); 51 PetscCall(MatSetFromOptions(*A)); 52 PetscCall(MatSetUp(*A)); 53 PetscCall(MatSeqSBAIJSetPreallocation(*A, 1, 0, rownz)); 54 PetscCall(PetscObjectTypeCompare((PetscObject)(*A), MATSEQSBAIJ, &sametype)); 55 PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported"); 56 } else { 57 PetscCall(MatSetType(*A, MATSEQAIJ)); 58 PetscCall(MatSetFromOptions(*A)); 59 PetscCall(MatSetUp(*A)); 60 PetscCall(MatSeqAIJSetPreallocation(*A, 0, rownz)); 61 PetscCall(PetscObjectTypeCompare((PetscObject)(*A), MATSEQAIJ, &sametype)); 62 PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported"); 63 } 64 /* Add values to the matrix, these correspond to lower triangular part for symmetric or skew matrices */ 65 for (j = 0; j < nz; j++) PetscCall(MatSetValues(*A, 1, &ia[j], 1, &ja[j], &val[j], INSERT_VALUES)); 66 67 /* Add values to upper triangular part for some cases */ 68 if (symmetric && aijonly) { 69 /* MatrixMarket matrix stores symm matrix in lower triangular part. Take its transpose */ 70 for (j = 0; j < nz; j++) PetscCall(MatSetValues(*A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES)); 71 } 72 if (skew) { 73 for (j = 0; j < nz; j++) { 74 val[j] = -val[j]; 75 PetscCall(MatSetValues(*A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES)); 76 } 77 } 78 PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY)); 79 PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY)); 80 PetscCall(PetscFree4(ia, ja, val, rownz)); 81 return 0; 82 } 83