Actual source code: ex72.c

  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";

  9: /*
 10:    NOTES:

 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.

 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 "mmloader.h"

 25: int main(int argc, char **argv)
 26: {
 27:   MM_typecode matcode;
 28:   FILE       *file;
 29:   PetscInt    M, N, nz;
 30:   Mat         A;
 31:   char        filein[PETSC_MAX_PATH_LEN], fileout[PETSC_MAX_PATH_LEN];
 32:   char        ordering[256] = MATORDERINGRCM;
 33:   PetscViewer view;
 34:   PetscBool   flag, symmetric = PETSC_FALSE, skew = PETSC_FALSE, aijonly = PETSC_FALSE, permute = PETSC_FALSE;
 35:   IS          rowperm = NULL, colperm = NULL;
 36:   PetscMPIInt size;
 37:   MatInfo     info;

 39:   PetscFunctionBeginUser;
 40:   PetscCall(PetscInitialize(&argc, &argv, NULL, 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!");

 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();

 55:   PetscCall(MatCreateFromMTX(&A, filein, aijonly));
 56:   PetscCall(PetscFOpen(PETSC_COMM_SELF, filein, "r", &file));
 57:   PetscCallExternal(mm_read_banner, file, &matcode);
 58:   if (mm_is_symmetric(matcode)) symmetric = PETSC_TRUE;
 59:   if (mm_is_skew(matcode)) skew = PETSC_TRUE;
 60:   PetscCallExternal(mm_write_banner, stdout, matcode);
 61:   PetscCallExternal(mm_read_mtx_crd_size, file, &M, &N, &nz);
 62:   PetscCall(PetscFClose(PETSC_COMM_SELF, file));
 63:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "M: %d, N: %d, nnz: %d\n", M, N, nz));
 64:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Reading matrix completes.\n"));
 65:   if (permute) {
 66:     Mat Aperm;
 67:     PetscCall(MatGetOrdering(A, ordering, &rowperm, &colperm));
 68:     PetscCall(MatPermute(A, rowperm, colperm, &Aperm));
 69:     PetscCall(MatDestroy(&A));
 70:     A = Aperm; /* Replace original operator with permuted version */
 71:   }

 73:   /* Write out matrix */
 74:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Writing matrix to binary file %s using PETSc %s format ...\n", fileout, (symmetric && !aijonly) ? "SBAIJ" : "AIJ"));
 75:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, fileout, FILE_MODE_WRITE, &view));
 76:   PetscCall(MatView(A, view));
 77:   PetscCall(PetscViewerDestroy(&view));
 78:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Writing matrix completes.\n"));
 79:   PetscCall(MatGetInfo(A, MAT_GLOBAL_SUM, &info));
 80:   PetscCheck(info.nz_used == (!skew ? nz : 2 * nz), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Different number of nonzero in MM matrix and PETSc Mat, %d != %g", nz, (double)info.nz_used);
 81:   PetscCall(MatDestroy(&A));
 82:   PetscCall(ISDestroy(&rowperm));
 83:   PetscCall(ISDestroy(&colperm));
 84:   PetscCall(PetscFinalize());
 85:   return 0;
 86: }

 88: /*TEST

 90:    build:
 91:       requires: !complex double !defined(PETSC_USE_64BIT_INDICES)
 92:       depends: mmloader.c mmio.c

 94:    test:
 95:       suffix: 1
 96:       args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij
 97:       output_file: output/ex72_1.out

 99:    test:
100:       suffix: 2
101:       args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/LFAT5.mtx -fout petscmat.sbaij
102:       output_file: output/ex72_2.out

104:    test:
105:       suffix: 2_permute
106:       args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/LFAT5.mtx -fout petscmat.sbaij -permute rcm
107:       output_file: output/ex72_2.out

109:    test:
110:       suffix: 3
111:       args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/m_05_05_crk.mtx -fout petscmat2.aij
112:       output_file: output/ex72_3.out

114:    test:
115:       suffix: 4
116:       args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij -permute rcm
117:       output_file: output/ex72_4.out
118: TEST*/