1 static char help[] = "Read a non-complex sparse matrix from a Matrix Market (v. 2.0) file\n\ 2 and write it to a file in petsc sparse binary format. If the matrix is symmetric, the binary file is in \n\ 3 PETSc MATSBAIJ format, otherwise it is in MATAIJ format \n\ 4 Usage: ./ex72 -fin <infile> -fout <outfile> \n\ 5 (See https://math.nist.gov/MatrixMarket/ for details.)\n\ 6 The option -permute <natural,rcm,nd,...> permutes the matrix using the ordering type.\n\ 7 The option -aij_only allows to use MATAIJ for all cases.\n\\n"; 8 9 /* 10 NOTES: 11 12 1) Matrix Market files are always 1-based, i.e. the index of the first 13 element of a matrix is (1,1), not (0,0) as in C. ADJUST THESE 14 OFFSETS ACCORDINGLY offsets accordingly when reading and writing 15 to files. 16 17 2) ANSI C requires one to use the "l" format modifier when reading 18 double precision floating point numbers in scanf() and 19 its variants. For example, use "%lf", "%lg", or "%le" 20 when reading doubles, otherwise errors will occur. 21 */ 22 #include <petscmat.h> 23 #include "ex72mmio.h" 24 25 int main(int argc, char **argv) { 26 MM_typecode matcode; 27 FILE *file; 28 PetscInt M, N, ninput; 29 PetscInt *ia, *ja; 30 Mat A; 31 char filein[PETSC_MAX_PATH_LEN], fileout[PETSC_MAX_PATH_LEN]; 32 char ordering[256] = MATORDERINGRCM; 33 PetscInt i, j, nz, *rownz; 34 PetscScalar *val, zero = 0.0; 35 PetscViewer view; 36 PetscBool sametype, flag, symmetric = PETSC_FALSE, skew = PETSC_FALSE, real = PETSC_FALSE, pattern = PETSC_FALSE, aijonly = PETSC_FALSE, permute = PETSC_FALSE; 37 IS rowperm = NULL, colperm = NULL; 38 PetscMPIInt size; 39 40 PetscInitialize(&argc, &argv, (char *)0, help); 41 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 42 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 43 44 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Matrix Market example options", ""); 45 { 46 PetscCall(PetscOptionsString("-fin", "Input Matrix Market file", "", filein, filein, sizeof(filein), &flag)); 47 PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Please use -fin <filename> to specify the input file name!"); 48 PetscCall(PetscOptionsString("-fout", "Output file in petsc sparse binary format", "", fileout, fileout, sizeof(fileout), &flag)); 49 PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Please use -fout <filename> to specify the output file name!"); 50 PetscCall(PetscOptionsBool("-aij_only", "Use MATAIJ for all cases", "", aijonly, &aijonly, NULL)); 51 PetscCall(PetscOptionsFList("-permute", "Permute matrix and vector to solving in new ordering", "", MatOrderingList, ordering, ordering, sizeof(ordering), &permute)); 52 } 53 PetscOptionsEnd(); 54 55 /* Read in matrix */ 56 PetscCall(PetscFOpen(PETSC_COMM_SELF, filein, "r", &file)); 57 58 PetscCheck(mm_read_banner(file, &matcode) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not process Matrix Market banner."); 59 60 /* This is how one can screen matrix types if their application */ 61 /* only supports a subset of the Matrix Market data types. */ 62 PetscCheck(mm_is_matrix(matcode) && mm_is_sparse(matcode), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input must be a sparse matrix. Market Market type: [%s]", mm_typecode_to_str(matcode)); 63 64 if (mm_is_symmetric(matcode)) symmetric = PETSC_TRUE; 65 if (mm_is_skew(matcode)) skew = PETSC_TRUE; 66 if (mm_is_real(matcode)) real = PETSC_TRUE; 67 if (mm_is_pattern(matcode)) pattern = PETSC_TRUE; 68 69 /* Find out size of sparse matrix .... */ 70 PetscCheck(mm_read_mtx_crd_size(file, &M, &N, &nz) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Size of sparse matrix is wrong."); 71 72 PetscCall(mm_write_banner(stdout, matcode)); 73 PetscCall(PetscPrintf(PETSC_COMM_SELF, "M: %d, N: %d, nnz: %d\n", M, N, nz)); 74 75 /* Reseve memory for matrices */ 76 PetscCall(PetscMalloc4(nz, &ia, nz, &ja, nz, &val, M, &rownz)); 77 for (i = 0; i < M; i++) rownz[i] = 1; /* Since we will add 0.0 to diagonal entries */ 78 79 /* NOTE: when reading in doubles, ANSI C requires the use of the "l" */ 80 /* specifier as in "%lg", "%lf", "%le", otherwise errors will occur */ 81 /* (ANSI C X3.159-1989, Sec. 4.9.6.2, p. 136 lines 13-15) */ 82 83 for (i = 0; i < nz; i++) { 84 if (pattern) { 85 ninput = fscanf(file, "%d %d\n", &ia[i], &ja[i]); 86 PetscCheck(ninput >= 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Badly formatted input file"); 87 val[i] = 1.0; 88 } else if (real) { 89 ninput = fscanf(file, "%d %d %lg\n", &ia[i], &ja[i], &val[i]); 90 PetscCheck(ninput >= 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Badly formatted input file"); 91 } 92 ia[i]--; 93 ja[i]--; /* adjust from 1-based to 0-based */ 94 if (ia[i] != ja[i]) { /* already counted the diagonals above */ 95 if ((symmetric && aijonly) || skew) { /* transpose */ 96 rownz[ia[i]]++; 97 rownz[ja[i]]++; 98 } else rownz[ia[i]]++; 99 } 100 } 101 PetscCall(PetscFClose(PETSC_COMM_SELF, file)); 102 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Reading matrix completes.\n")); 103 104 /* Create, preallocate, and then assemble the matrix */ 105 PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 106 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 107 108 if (symmetric && !aijonly) { 109 PetscCall(MatSetType(A, MATSEQSBAIJ)); 110 PetscCall(MatSetFromOptions(A)); 111 PetscCall(MatSetUp(A)); 112 PetscCall(MatSeqSBAIJSetPreallocation(A, 1, 0, rownz)); 113 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &sametype)); 114 PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported"); 115 } else { 116 PetscCall(MatSetType(A, MATSEQAIJ)); 117 PetscCall(MatSetFromOptions(A)); 118 PetscCall(MatSetUp(A)); 119 PetscCall(MatSeqAIJSetPreallocation(A, 0, rownz)); 120 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &sametype)); 121 PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported"); 122 } 123 124 /* Add zero to diagonals, in case the matrix missing diagonals */ 125 for (j = 0; j < M; j++) PetscCall(MatSetValues(A, 1, &j, 1, &j, &zero, INSERT_VALUES)); 126 /* Add values to the matrix, these correspond to lower triangular part for symmetric or skew matrices */ 127 for (j = 0; j < nz; j++) PetscCall(MatSetValues(A, 1, &ia[j], 1, &ja[j], &val[j], INSERT_VALUES)); 128 129 /* Add values to upper triangular part for some cases */ 130 if (symmetric && aijonly) { 131 /* MatrixMarket matrix stores symm matrix in lower triangular part. Take its transpose */ 132 for (j = 0; j < nz; j++) PetscCall(MatSetValues(A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES)); 133 } 134 if (skew) { 135 for (j = 0; j < nz; j++) { 136 val[j] = -val[j]; 137 PetscCall(MatSetValues(A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES)); 138 } 139 } 140 141 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 142 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 143 144 if (permute) { 145 Mat Aperm; 146 PetscCall(MatGetOrdering(A, ordering, &rowperm, &colperm)); 147 PetscCall(MatPermute(A, rowperm, colperm, &Aperm)); 148 PetscCall(MatDestroy(&A)); 149 A = Aperm; /* Replace original operator with permuted version */ 150 } 151 152 /* Write out matrix */ 153 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Writing matrix to binary file %s using PETSc %s format ...\n", fileout, (symmetric && !aijonly) ? "SBAIJ" : "AIJ")); 154 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, fileout, FILE_MODE_WRITE, &view)); 155 PetscCall(MatView(A, view)); 156 PetscCall(PetscViewerDestroy(&view)); 157 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Writing matrix completes.\n")); 158 159 PetscCall(PetscFree4(ia, ja, val, rownz)); 160 PetscCall(MatDestroy(&A)); 161 PetscCall(ISDestroy(&rowperm)); 162 PetscCall(ISDestroy(&colperm)); 163 PetscCall(PetscFinalize()); 164 return 0; 165 } 166 167 /*TEST 168 169 build: 170 requires: !complex double !defined(PETSC_USE_64BIT_INDICES) 171 depends: ex72mmio.c 172 173 test: 174 suffix: 1 175 args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij 176 output_file: output/ex72_1.out 177 178 test: 179 suffix: 2 180 args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/LFAT5.mtx -fout petscmat.sbaij 181 output_file: output/ex72_2.out 182 183 test: 184 suffix: 3 185 args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/m_05_05_crk.mtx -fout petscmat2.aij 186 output_file: output/ex72_3.out 187 188 test: 189 suffix: 4 190 args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij -permute rcm 191 output_file: output/ex72_4.out 192 TEST*/ 193