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; 34*d0609cedSBarry Smith PetscInt i,j,nz,*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; 39*d0609cedSBarry Smith PetscMPIInt size; 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscInitialize(&argc,&argv,(char *)0,help); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 43ba4e87b7SRey Koki PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 44c4762a1bSJed Brown 45*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Matrix Market example options",""); 46ba4e87b7SRey Koki { 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-fin","Input Matrix Market file","",filein,filein,sizeof(filein),&flag)); 48ba4e87b7SRey Koki PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Please use -fin <filename> to specify the input file name!"); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-fout","Output file in petsc sparse binary format","",fileout,fileout,sizeof(fileout),&flag)); 50ba4e87b7SRey Koki PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Please use -fout <filename> to specify the output file name!"); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-aij_only","Use MATAIJ for all cases","",aijonly,&aijonly,NULL)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-permute","Permute matrix and vector to solving in new ordering","",MatOrderingList,ordering,ordering,sizeof(ordering),&permute)); 53ba4e87b7SRey Koki } 54*d0609cedSBarry Smith PetscOptionsEnd(); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Read in matrix */ 579566063dSJacob Faibussowitsch PetscCall(PetscFOpen(PETSC_COMM_SELF,filein,"r",&file)); 58c4762a1bSJed Brown 59ba4e87b7SRey Koki PetscCheck(mm_read_banner(file, &matcode) == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Could not process Matrix Market banner."); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* This is how one can screen matrix types if their application */ 62c4762a1bSJed Brown /* only supports a subset of the Matrix Market data types. */ 63ba4e87b7SRey 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)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown if (mm_is_symmetric(matcode)) symmetric = PETSC_TRUE; 66c4762a1bSJed Brown if (mm_is_skew(matcode)) skew = PETSC_TRUE; 67c4762a1bSJed Brown if (mm_is_real(matcode)) real = PETSC_TRUE; 68c4762a1bSJed Brown if (mm_is_pattern(matcode)) pattern = PETSC_TRUE; 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Find out size of sparse matrix .... */ 71ba4e87b7SRey 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."); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(mm_write_banner(stdout, matcode)); 749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"M: %d, N: %d, nnz: %d\n",M,N,nz)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Reseve memory for matrices */ 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(nz,&ia,nz,&ja,nz,&val,M,&rownz)); 78c4762a1bSJed Brown for (i=0; i<M; i++) rownz[i] = 1; /* Since we will add 0.0 to diagonal entries */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* NOTE: when reading in doubles, ANSI C requires the use of the "l" */ 81c4762a1bSJed Brown /* specifier as in "%lg", "%lf", "%le", otherwise errors will occur */ 82c4762a1bSJed Brown /* (ANSI C X3.159-1989, Sec. 4.9.6.2, p. 136 lines 13-15) */ 83c4762a1bSJed Brown 84c4762a1bSJed Brown for (i=0; i<nz; i++) { 85c4762a1bSJed Brown if (pattern) { 86c4762a1bSJed Brown ninput = fscanf(file, "%d %d\n", &ia[i], &ja[i]); 875f80ce2aSJacob Faibussowitsch PetscCheck(ninput >= 2,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 88c4762a1bSJed Brown val[i] = 1.0; 89c4762a1bSJed Brown } else if (real) { 90c4762a1bSJed Brown ninput = fscanf(file, "%d %d %lg\n", &ia[i], &ja[i], &val[i]); 915f80ce2aSJacob Faibussowitsch PetscCheck(ninput >= 3,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown ia[i]--; ja[i]--; /* adjust from 1-based to 0-based */ 94c4762a1bSJed Brown if (ia[i] != ja[i]) { /* already counted the diagonals above */ 95c4762a1bSJed Brown if ((symmetric && aijonly) || skew) { /* transpose */ 96c4762a1bSJed Brown rownz[ia[i]]++; 97c4762a1bSJed Brown rownz[ja[i]]++; 98c4762a1bSJed Brown } else rownz[ia[i]]++; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 1019566063dSJacob Faibussowitsch PetscCall(PetscFClose(PETSC_COMM_SELF,file)); 1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Reading matrix completes.\n")); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Create, preallocate, and then assemble the matrix */ 1059566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 1069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown if (symmetric && !aijonly) { 1099566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQSBAIJ)); 1109566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1119566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1129566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A,1,0,rownz)); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sametype)); 1145f80ce2aSJacob Faibussowitsch PetscCheck(sametype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only AIJ and SBAIJ are supported. Your mattype is not supported"); 115c4762a1bSJed Brown } else { 1169566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 1179566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1189566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1199566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,0,rownz)); 1209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&sametype)); 1215f80ce2aSJacob Faibussowitsch PetscCheck(sametype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only AIJ and SBAIJ are supported. Your mattype is not supported"); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* Add zero to diagonals, in case the matrix missing diagonals */ 1259566063dSJacob Faibussowitsch for (j=0; j<M; j++) PetscCall(MatSetValues(A,1,&j,1,&j,&zero,INSERT_VALUES)); 126c4762a1bSJed Brown /* Add values to the matrix, these correspond to lower triangular part for symmetric or skew matrices */ 1279566063dSJacob Faibussowitsch for (j=0; j<nz; j++) PetscCall(MatSetValues(A,1,&ia[j],1,&ja[j],&val[j],INSERT_VALUES)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Add values to upper triangular part for some cases */ 130c4762a1bSJed Brown if (symmetric && aijonly) { 131c4762a1bSJed Brown /* MatrixMarket matrix stores symm matrix in lower triangular part. Take its transpose */ 1329566063dSJacob Faibussowitsch for (j=0; j<nz; j++) PetscCall(MatSetValues(A,1,&ja[j],1,&ia[j],&val[j],INSERT_VALUES)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown if (skew) { 135c4762a1bSJed Brown for (j=0; j<nz; j++) { 136c4762a1bSJed Brown val[j] = -val[j]; 1379566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&ja[j],1,&ia[j],&val[j],INSERT_VALUES)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 143c4762a1bSJed Brown 144ba4e87b7SRey Koki if (permute) { 145ba4e87b7SRey Koki Mat Aperm; 1469566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,ordering,&rowperm,&colperm)); 1479566063dSJacob Faibussowitsch PetscCall(MatPermute(A,rowperm,colperm,&Aperm)); 1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 149ba4e87b7SRey Koki A = Aperm; /* Replace original operator with permuted version */ 150ba4e87b7SRey Koki } 151ba4e87b7SRey Koki 152c4762a1bSJed Brown /* Write out matrix */ 1539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Writing matrix to binary file %s using PETSc %s format ...\n",fileout,(symmetric && !aijonly)?"SBAIJ":"AIJ")); 1549566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,fileout,FILE_MODE_WRITE,&view)); 1559566063dSJacob Faibussowitsch PetscCall(MatView(A,view)); 1569566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Writing matrix completes.\n")); 158c4762a1bSJed Brown 1599566063dSJacob Faibussowitsch PetscCall(PetscFree4(ia,ja,val,rownz)); 1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rowperm)); 1629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&colperm)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 164c4762a1bSJed Brown return 0; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /*TEST 168c4762a1bSJed Brown 169c4762a1bSJed Brown build: 170dfd57a17SPierre Jolivet requires: !complex double !defined(PETSC_USE_64BIT_INDICES) 171c4762a1bSJed Brown depends: ex72mmio.c 172c4762a1bSJed Brown 173c4762a1bSJed Brown test: 174c4762a1bSJed Brown suffix: 1 175c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij 176c4762a1bSJed Brown output_file: output/ex72_1.out 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 179c4762a1bSJed Brown suffix: 2 180c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/LFAT5.mtx -fout petscmat.sbaij 181c4762a1bSJed Brown output_file: output/ex72_2.out 182c4762a1bSJed Brown 183c4762a1bSJed Brown test: 184c4762a1bSJed Brown suffix: 3 185c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/m_05_05_crk.mtx -fout petscmat2.aij 186c4762a1bSJed Brown output_file: output/ex72_3.out 187ba4e87b7SRey Koki 188ba4e87b7SRey Koki test: 189ba4e87b7SRey Koki suffix: 4 190ba4e87b7SRey Koki args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij -permute rcm 191ba4e87b7SRey Koki output_file: output/ex72_4.out 192c4762a1bSJed Brown TEST*/ 193