189928cc5SHong Zhang #include "mmloader.h"
289928cc5SHong Zhang
MatCreateFromMTX(Mat * A,const char * filein,PetscBool aijonly)389928cc5SHong Zhang PetscErrorCode MatCreateFromMTX(Mat *A, const char *filein, PetscBool aijonly)
489928cc5SHong Zhang {
589928cc5SHong Zhang MM_typecode matcode;
689928cc5SHong Zhang FILE *file;
789928cc5SHong Zhang PetscInt M, N, ninput;
889928cc5SHong Zhang PetscInt *ia, *ja;
989928cc5SHong Zhang PetscInt i, j, nz, *rownz;
1089928cc5SHong Zhang PetscScalar *val;
1189928cc5SHong Zhang PetscBool sametype, symmetric = PETSC_FALSE, skew = PETSC_FALSE;
1289928cc5SHong Zhang
1389928cc5SHong Zhang /* Read in matrix */
143ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
1589928cc5SHong Zhang PetscCall(PetscFOpen(PETSC_COMM_SELF, filein, "r", &file));
1689928cc5SHong Zhang PetscCheck(mm_read_banner(file, &matcode) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not process Matrix Market banner.");
1789928cc5SHong Zhang /* This is how one can screen matrix types if their application */
1889928cc5SHong Zhang /* only supports a subset of the Matrix Market data types. */
1989928cc5SHong 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));
2089928cc5SHong Zhang
2189928cc5SHong Zhang if (mm_is_symmetric(matcode)) symmetric = PETSC_TRUE;
2289928cc5SHong Zhang if (mm_is_skew(matcode)) skew = PETSC_TRUE;
2389928cc5SHong Zhang
2489928cc5SHong Zhang /* Find out size of sparse matrix .... */
2589928cc5SHong 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.");
2689928cc5SHong Zhang
2789928cc5SHong Zhang /* Reserve memory for matrices */
2889928cc5SHong Zhang PetscCall(PetscMalloc4(nz, &ia, nz, &ja, nz, &val, M, &rownz));
2989928cc5SHong Zhang for (i = 0; i < M; i++) rownz[i] = 0;
3089928cc5SHong Zhang
3189928cc5SHong Zhang /* NOTE: when reading in doubles, ANSI C requires the use of the "l" */
3289928cc5SHong Zhang /* specifier as in "%lg", "%lf", "%le", otherwise errors will occur */
3389928cc5SHong Zhang /* (ANSI C X3.159-1989, Sec. 4.9.6.2, p. 136 lines 13-15) */
3489928cc5SHong Zhang for (i = 0; i < nz; i++) {
3589928cc5SHong Zhang ninput = fscanf(file, "%d %d %lg\n", &ia[i], &ja[i], &val[i]);
3689928cc5SHong Zhang PetscCheck(ninput >= 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Badly formatted input file");
3789928cc5SHong Zhang ia[i]--;
3889928cc5SHong Zhang ja[i]--; /* adjust from 1-based to 0-based */
3989928cc5SHong Zhang if ((symmetric && aijonly) || skew) { /* transpose */
4089928cc5SHong Zhang rownz[ia[i]]++;
419f114c04SHong Zhang if (ja[i] != ia[i]) rownz[ja[i]]++;
42c33121ecSPierre Jolivet } else {
43c33121ecSPierre Jolivet if (symmetric) rownz[ja[i]]++;
44c33121ecSPierre Jolivet else rownz[ia[i]]++;
45c33121ecSPierre Jolivet }
4689928cc5SHong Zhang }
4789928cc5SHong Zhang PetscCall(PetscFClose(PETSC_COMM_SELF, file));
4889928cc5SHong Zhang
4989928cc5SHong Zhang /* Create, preallocate, and then assemble the matrix */
5089928cc5SHong Zhang PetscCall(MatCreate(PETSC_COMM_SELF, A));
5189928cc5SHong Zhang PetscCall(MatSetSizes(*A, PETSC_DECIDE, PETSC_DECIDE, M, N));
5289928cc5SHong Zhang
5389928cc5SHong Zhang if (symmetric && !aijonly) {
5489928cc5SHong Zhang PetscCall(MatSetType(*A, MATSEQSBAIJ));
5589928cc5SHong Zhang PetscCall(MatSetFromOptions(*A));
5689928cc5SHong Zhang PetscCall(MatSeqSBAIJSetPreallocation(*A, 1, 0, rownz));
57*f4f49eeaSPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)*A, MATSEQSBAIJ, &sametype));
5889928cc5SHong Zhang PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported");
5989928cc5SHong Zhang } else {
6089928cc5SHong Zhang PetscCall(MatSetType(*A, MATSEQAIJ));
6189928cc5SHong Zhang PetscCall(MatSetFromOptions(*A));
6289928cc5SHong Zhang PetscCall(MatSeqAIJSetPreallocation(*A, 0, rownz));
63*f4f49eeaSPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)*A, MATSEQAIJ, &sametype));
6489928cc5SHong Zhang PetscCheck(sametype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Only AIJ and SBAIJ are supported. Your mattype is not supported");
6589928cc5SHong Zhang }
6689928cc5SHong Zhang /* Add values to the matrix, these correspond to lower triangular part for symmetric or skew matrices */
67c33121ecSPierre Jolivet if (!(symmetric && !aijonly))
6889928cc5SHong Zhang for (j = 0; j < nz; j++) PetscCall(MatSetValues(*A, 1, &ia[j], 1, &ja[j], &val[j], INSERT_VALUES));
6989928cc5SHong Zhang
7089928cc5SHong Zhang /* Add values to upper triangular part for some cases */
71c33121ecSPierre Jolivet if (symmetric && !aijonly) {
7289928cc5SHong Zhang /* MatrixMarket matrix stores symm matrix in lower triangular part. Take its transpose */
7389928cc5SHong Zhang for (j = 0; j < nz; j++) PetscCall(MatSetValues(*A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES));
7489928cc5SHong Zhang }
7589928cc5SHong Zhang if (skew) {
7689928cc5SHong Zhang for (j = 0; j < nz; j++) {
7789928cc5SHong Zhang val[j] = -val[j];
7889928cc5SHong Zhang PetscCall(MatSetValues(*A, 1, &ja[j], 1, &ia[j], &val[j], INSERT_VALUES));
7989928cc5SHong Zhang }
8089928cc5SHong Zhang }
8189928cc5SHong Zhang PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY));
8289928cc5SHong Zhang PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY));
8389928cc5SHong Zhang PetscCall(PetscFree4(ia, ja, val, rownz));
843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8589928cc5SHong Zhang }
86