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