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\ 6*ba4e87b7SRey 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]; 33*ba4e87b7SRey 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; 37*ba4e87b7SRey Koki PetscBool sametype,flag,symmetric = PETSC_FALSE,skew = PETSC_FALSE,real = PETSC_FALSE,pattern = PETSC_FALSE,aijonly = PETSC_FALSE, permute = PETSC_FALSE; 38*ba4e87b7SRey Koki IS rowperm = NULL,colperm = NULL; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscInitialize(&argc,&argv,(char *)0,help); 41ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 42*ba4e87b7SRey Koki PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 43c4762a1bSJed Brown 44*ba4e87b7SRey Koki ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Matrix Market example options","");CHKERRQ(ierr); 45*ba4e87b7SRey Koki { 46*ba4e87b7SRey Koki ierr = PetscOptionsString("-fin","Input Matrix Market file","",filein,filein,sizeof(filein),&flag);CHKERRQ(ierr); 47*ba4e87b7SRey Koki PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Please use -fin <filename> to specify the input file name!"); 48*ba4e87b7SRey Koki ierr = PetscOptionsString("-fout","Output file in petsc sparse binary format","",fileout,fileout,sizeof(fileout),&flag);CHKERRQ(ierr); 49*ba4e87b7SRey Koki PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Please use -fout <filename> to specify the output file name!"); 50*ba4e87b7SRey Koki ierr = PetscOptionsBool("-aij_only","Use MATAIJ for all cases","",aijonly,&aijonly,NULL);CHKERRQ(ierr); 51*ba4e87b7SRey Koki ierr = PetscOptionsFList("-permute","Permute matrix and vector to solving in new ordering","",MatOrderingList,ordering,ordering,sizeof(ordering),&permute);CHKERRQ(ierr); 52*ba4e87b7SRey Koki } 53*ba4e87b7SRey Koki ierr = PetscOptionsEnd();CHKERRQ(ierr); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Read in matrix */ 56c4762a1bSJed Brown ierr = PetscFOpen(PETSC_COMM_SELF,filein,"r",&file);CHKERRQ(ierr); 57c4762a1bSJed Brown 58*ba4e87b7SRey 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. */ 62*ba4e87b7SRey 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 .... */ 70*ba4e87b7SRey 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 72c4762a1bSJed Brown ierr = mm_write_banner(stdout, matcode);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"M: %d, N: %d, nnz: %d\n",M,N,nz);CHKERRQ(ierr); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Reseve memory for matrices */ 76c4762a1bSJed Brown ierr = PetscMalloc4(nz,&ia,nz,&ja,nz,&val,M,&rownz);CHKERRQ(ierr); 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]); 862c71b3e2SJacob Faibussowitsch PetscCheckFalse(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]); 902c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 } 100c4762a1bSJed Brown ierr = PetscFClose(PETSC_COMM_SELF,file);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Reading matrix completes.\n");CHKERRQ(ierr); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Create, preallocate, and then assemble the matrix */ 104c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 106c4762a1bSJed Brown 107c4762a1bSJed Brown if (symmetric && !aijonly) { 108c4762a1bSJed Brown ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = MatSeqSBAIJSetPreallocation(A,1,0,rownz);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sametype);CHKERRQ(ierr); 1132c71b3e2SJacob Faibussowitsch PetscCheckFalse(!sametype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only AIJ and SBAIJ are supported. Your mattype is not supported"); 114c4762a1bSJed Brown } else { 115c4762a1bSJed Brown ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,0,rownz);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&sametype);CHKERRQ(ierr); 1202c71b3e2SJacob Faibussowitsch PetscCheckFalse(!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 */ 124c4762a1bSJed Brown for (j=0; j<M; j++) { 125c4762a1bSJed Brown ierr = MatSetValues(A,1,&j,1,&j,&zero,INSERT_VALUES);CHKERRQ(ierr); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown /* Add values to the matrix, these correspond to lower triangular part for symmetric or skew matrices */ 128c4762a1bSJed Brown for (j=0; j<nz; j++) { 129c4762a1bSJed Brown ierr = MatSetValues(A,1,&ia[j],1,&ja[j],&val[j],INSERT_VALUES);CHKERRQ(ierr); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Add values to upper triangular part for some cases */ 133c4762a1bSJed Brown if (symmetric && aijonly) { 134c4762a1bSJed Brown /* MatrixMarket matrix stores symm matrix in lower triangular part. Take its transpose */ 135c4762a1bSJed Brown for (j=0; j<nz; j++) { 136c4762a1bSJed Brown ierr = MatSetValues(A,1,&ja[j],1,&ia[j],&val[j],INSERT_VALUES);CHKERRQ(ierr); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 139c4762a1bSJed Brown if (skew) { 140c4762a1bSJed Brown for (j=0; j<nz; j++) { 141c4762a1bSJed Brown val[j] = -val[j]; 142c4762a1bSJed Brown ierr = MatSetValues(A,1,&ja[j],1,&ia[j],&val[j],INSERT_VALUES);CHKERRQ(ierr); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 148c4762a1bSJed Brown 149*ba4e87b7SRey Koki if (permute) { 150*ba4e87b7SRey Koki Mat Aperm; 151*ba4e87b7SRey Koki ierr = MatGetOrdering(A,ordering,&rowperm,&colperm);CHKERRQ(ierr); 152*ba4e87b7SRey Koki ierr = MatPermute(A,rowperm,colperm,&Aperm);CHKERRQ(ierr); 153*ba4e87b7SRey Koki ierr = MatDestroy(&A);CHKERRQ(ierr); 154*ba4e87b7SRey Koki A = Aperm; /* Replace original operator with permuted version */ 155*ba4e87b7SRey Koki } 156*ba4e87b7SRey Koki 157c4762a1bSJed Brown /* Write out matrix */ 158c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Writing matrix to binary file %s using PETSc %s format ...\n",fileout,(symmetric && !aijonly)?"SBAIJ":"AIJ");CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,fileout,FILE_MODE_WRITE,&view);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = MatView(A,view);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = PetscViewerDestroy(&view);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Writing matrix completes.\n");CHKERRQ(ierr); 163c4762a1bSJed Brown 164c4762a1bSJed Brown ierr = PetscFree4(ia,ja,val,rownz);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 166*ba4e87b7SRey Koki ierr = ISDestroy(&rowperm);CHKERRQ(ierr); 167*ba4e87b7SRey Koki ierr = ISDestroy(&colperm);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = PetscFinalize();CHKERRQ(ierr); 169c4762a1bSJed Brown return 0; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /*TEST 173c4762a1bSJed Brown 174c4762a1bSJed Brown build: 175dfd57a17SPierre Jolivet requires: !complex double !defined(PETSC_USE_64BIT_INDICES) 176c4762a1bSJed Brown depends: ex72mmio.c 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 179c4762a1bSJed Brown suffix: 1 180c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij 181c4762a1bSJed Brown output_file: output/ex72_1.out 182c4762a1bSJed Brown 183c4762a1bSJed Brown test: 184c4762a1bSJed Brown suffix: 2 185c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/LFAT5.mtx -fout petscmat.sbaij 186c4762a1bSJed Brown output_file: output/ex72_2.out 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown suffix: 3 190c4762a1bSJed Brown args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/m_05_05_crk.mtx -fout petscmat2.aij 191c4762a1bSJed Brown output_file: output/ex72_3.out 192*ba4e87b7SRey Koki 193*ba4e87b7SRey Koki test: 194*ba4e87b7SRey Koki suffix: 4 195*ba4e87b7SRey Koki args: -fin ${wPETSC_DIR}/share/petsc/datafiles/matrices/amesos2_test_mat0.mtx -fout petscmat.aij -permute rcm 196*ba4e87b7SRey Koki output_file: output/ex72_4.out 197c4762a1bSJed Brown TEST*/ 198