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