xref: /petsc/src/mat/tests/ex55.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown /* Usage: mpiexec -n <np> ex55 -verbose <0 or 1> */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            C,A,B,D;
10c4762a1bSJed Brown   PetscInt       i,j,ntypes,bs,mbs,m,block,d_nz=6, o_nz=3,col[3],row,verbose=0;
11c4762a1bSJed Brown   PetscMPIInt    size,rank;
12c4762a1bSJed Brown   MatType        type[9];
13c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
14c4762a1bSJed Brown   PetscViewer    fd;
15c4762a1bSJed Brown   PetscBool      equal,flg_loadmat,flg,issymmetric;
16c4762a1bSJed Brown   PetscScalar    value[3];
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg_loadmat));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-testseqaij",&flg));
26c4762a1bSJed Brown   if (flg) {
27c4762a1bSJed Brown     if (size == 1) {
28c4762a1bSJed Brown       type[0] = MATSEQAIJ;
29c4762a1bSJed Brown     } else {
30c4762a1bSJed Brown       type[0] = MATMPIAIJ;
31c4762a1bSJed Brown     }
32c4762a1bSJed Brown   } else {
33c4762a1bSJed Brown     type[0] = MATAIJ;
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   if (size == 1) {
36c4762a1bSJed Brown     ntypes  = 3;
37c4762a1bSJed Brown     type[1] = MATSEQBAIJ;
38c4762a1bSJed Brown     type[2] = MATSEQSBAIJ;
39c4762a1bSJed Brown   } else {
40c4762a1bSJed Brown     ntypes  = 3;
41c4762a1bSJed Brown     type[1] = MATMPIBAIJ;
42c4762a1bSJed Brown     type[2] = MATMPISBAIJ;
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* input matrix C */
46c4762a1bSJed Brown   if (flg_loadmat) {
47c4762a1bSJed Brown     /* Open binary file. */
489566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown     /* Load a baij matrix, then destroy the viewer. */
519566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
52c4762a1bSJed Brown     if (size == 1) {
539566063dSJacob Faibussowitsch       PetscCall(MatSetType(C,MATSEQBAIJ));
54c4762a1bSJed Brown     } else {
559566063dSJacob Faibussowitsch       PetscCall(MatSetType(C,MATMPIBAIJ));
56c4762a1bSJed Brown     }
579566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
589566063dSJacob Faibussowitsch     PetscCall(MatLoad(C,fd));
599566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
60c4762a1bSJed Brown   } else { /* Create a baij mat with bs>1  */
61c4762a1bSJed Brown     bs   = 2; mbs=8;
629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
6408401ef6SPierre Jolivet     PetscCheck(bs > 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," bs must be >1 in this case");
65c4762a1bSJed Brown     m    = mbs*bs;
669566063dSJacob Faibussowitsch     PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,m,m,d_nz,NULL,o_nz,NULL,&C));
67c4762a1bSJed Brown     for (block=0; block<mbs; block++) {
68c4762a1bSJed Brown       /* diagonal blocks */
69c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
70c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
71c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
729566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES));
73c4762a1bSJed Brown       }
74c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
75c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
769566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
79c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
809566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
81c4762a1bSJed Brown     }
82c4762a1bSJed Brown     /* off-diagonal blocks */
83c4762a1bSJed Brown     value[0]=-1.0;
84c4762a1bSJed Brown     for (i=0; i<(mbs-1)*bs; i++) {
85c4762a1bSJed Brown       col[0]=i+bs;
869566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&i,1,col,value,INSERT_VALUES));
87c4762a1bSJed Brown       col[0]=i; row=i+bs;
889566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&row,1,col,value,INSERT_VALUES));
89c4762a1bSJed Brown     }
909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown   }
93c4762a1bSJed Brown   {
94c4762a1bSJed Brown     /* Check the symmetry of C because it will be converted to a sbaij matrix */
95c4762a1bSJed Brown     Mat Ctrans;
969566063dSJacob Faibussowitsch     PetscCall(MatTranspose(C,MAT_INITIAL_MATRIX,&Ctrans));
979566063dSJacob Faibussowitsch     PetscCall(MatEqual(C,Ctrans,&flg));
98c4762a1bSJed Brown /*
99c4762a1bSJed Brown     {
1009566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C,1.,Ctrans,DIFFERENT_NONZERO_PATTERN));
101c4762a1bSJed Brown       flg  = PETSC_TRUE;
102c4762a1bSJed Brown     }
103c4762a1bSJed Brown */
1049566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_SYMMETRIC,flg));
1059566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Ctrans));
106c4762a1bSJed Brown   }
1079566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(C,0.0,&issymmetric));
1089566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(C,NULL,"-view_mat"));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* convert C to other formats */
111c4762a1bSJed Brown   for (i=0; i<ntypes; i++) {
112c4762a1bSJed Brown     PetscBool ismpisbaij,isseqsbaij;
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type[i],MATMPISBAIJ,&ismpisbaij));
1159566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type[i],MATMPISBAIJ,&isseqsbaij));
116c4762a1bSJed Brown     if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
1179566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,type[i],MAT_INITIAL_MATRIX,&A));
1189566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,C,10,&equal));
11928b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from BAIJ to %s",type[i]);
120c4762a1bSJed Brown     for (j=i+1; j<ntypes; j++) {
1219566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij));
1229566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij));
123c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
124c4762a1bSJed Brown       if (verbose>0) {
1259566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD," \n[%d] test conversion between %s and %s\n",rank,type[i],type[j]));
126c4762a1bSJed Brown       }
127c4762a1bSJed Brown 
128dd400576SPatrick Sanan       if (rank == 0 && verbose) printf("Convert %s A to %s B\n",type[i],type[j]);
1299566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,type[j],MAT_INITIAL_MATRIX,&B));
130c4762a1bSJed Brown       /*
131c4762a1bSJed Brown       if (j == 2) {
1329566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF," A: %s\n",type[i]));
1339566063dSJacob Faibussowitsch         PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
1349566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF," B: %s\n",type[j]));
1359566063dSJacob Faibussowitsch         PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
136c4762a1bSJed Brown       }
137c4762a1bSJed Brown        */
1389566063dSJacob Faibussowitsch       PetscCall(MatMultEqual(A,B,10,&equal));
13928b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[i],type[j]);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown       if (size == 1 || j != 2) { /* Matconvert from mpisbaij mat to other formats are not supported */
142dd400576SPatrick Sanan         if (rank == 0 && verbose) printf("Convert %s B to %s D\n",type[j],type[i]);
1439566063dSJacob Faibussowitsch         PetscCall(MatConvert(B,type[i],MAT_INITIAL_MATRIX,&D));
1449566063dSJacob Faibussowitsch         PetscCall(MatMultEqual(B,D,10,&equal));
14528b400f6SJacob Faibussowitsch         PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[j],type[i]);
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&D));
148c4762a1bSJed Brown       }
1499566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
1509566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&D));
151c4762a1bSJed Brown     }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown     /* Test in-place convert */
154c4762a1bSJed Brown     if (size == 1) { /* size > 1 is not working yet! */
155c4762a1bSJed Brown       j = (i+1)%ntypes;
1569566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij));
1579566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij));
158c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
159c4762a1bSJed Brown       /* printf("[%d] i: %d, j: %d\n",rank,i,j); */
1609566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,type[j],MAT_INPLACE_MATRIX,&A));
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown 
1639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* test BAIJ to MATIS */
167c4762a1bSJed Brown   if (size > 1) {
168c4762a1bSJed Brown     MatType ctype;
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch     PetscCall(MatGetType(C,&ctype));
1719566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,MATIS,MAT_INITIAL_MATRIX,&A));
1729566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,C,10,&equal));
1739566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A,NULL,"-view_conv"));
174c4762a1bSJed Brown     if (!equal) {
175c4762a1bSJed Brown       Mat C2;
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
1789566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
1799566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
1809566063dSJacob Faibussowitsch       PetscCall(MatChop(C2,PETSC_SMALL));
1819566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_err"));
1829566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
183c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
184c4762a1bSJed Brown     }
1859566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,MATIS,MAT_REUSE_MATRIX,&A));
1869566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,C,10,&equal));
1879566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A,NULL,"-view_conv"));
188c4762a1bSJed Brown     if (!equal) {
189c4762a1bSJed Brown       Mat C2;
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
1929566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
1939566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
1949566063dSJacob Faibussowitsch       PetscCall(MatChop(C2,PETSC_SMALL));
1959566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_err"));
1969566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
197c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
198c4762a1bSJed Brown     }
1999566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
2009566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&A));
2019566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATIS,MAT_INPLACE_MATRIX,&A));
2029566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A,NULL,"-view_conv"));
2039566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,C,10,&equal));
204c4762a1bSJed Brown     if (!equal) {
205c4762a1bSJed Brown       Mat C2;
206c4762a1bSJed Brown 
2079566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(A,NULL,"-view_conv"));
2089566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
2099566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
2109566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
2119566063dSJacob Faibussowitsch       PetscCall(MatChop(C2,PETSC_SMALL));
2129566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2,NULL,"-view_err"));
2139566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
214c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
215c4762a1bSJed Brown     }
2169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
217c4762a1bSJed Brown   }
2189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
219c4762a1bSJed Brown 
2209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
221b122ec5aSJacob Faibussowitsch   return 0;
222c4762a1bSJed Brown }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown /*TEST
225c4762a1bSJed Brown 
226c4762a1bSJed Brown    test:
227c4762a1bSJed Brown 
228c4762a1bSJed Brown    test:
229c4762a1bSJed Brown       suffix: 2
230c4762a1bSJed Brown       nsize: 3
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    testset:
233c4762a1bSJed Brown       requires: parmetis
234c4762a1bSJed Brown       output_file: output/ex55_1.out
235c4762a1bSJed Brown       nsize: 3
236c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis
237c4762a1bSJed Brown       test:
238c4762a1bSJed Brown         suffix: matis_baij_parmetis_nd
239c4762a1bSJed Brown       test:
240c4762a1bSJed Brown         suffix: matis_aij_parmetis_nd
241c4762a1bSJed Brown         args: -testseqaij
242c4762a1bSJed Brown       test:
243dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
244c4762a1bSJed Brown         suffix: matis_poisson1_parmetis_nd
245c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    testset:
248dfd57a17SPierre Jolivet       requires: ptscotch defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
249c4762a1bSJed Brown       output_file: output/ex55_1.out
250c4762a1bSJed Brown       nsize: 4
251c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch
252c4762a1bSJed Brown       test:
253c4762a1bSJed Brown         suffix: matis_baij_ptscotch_nd
254c4762a1bSJed Brown       test:
255c4762a1bSJed Brown         suffix: matis_aij_ptscotch_nd
256c4762a1bSJed Brown         args: -testseqaij
257c4762a1bSJed Brown       test:
258dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
259c4762a1bSJed Brown         suffix: matis_poisson1_ptscotch_nd
260c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
261c4762a1bSJed Brown 
262c4762a1bSJed Brown TEST*/
263