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