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