xref: /petsc/src/mat/tests/ex55.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg_loadmat));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-testseqaij",&flg));
25c4762a1bSJed Brown   if (flg) {
26c4762a1bSJed Brown     if (size == 1) {
27c4762a1bSJed Brown       type[0] = MATSEQAIJ;
28c4762a1bSJed Brown     } else {
29c4762a1bSJed Brown       type[0] = MATMPIAIJ;
30c4762a1bSJed Brown     }
31c4762a1bSJed Brown   } else {
32c4762a1bSJed Brown     type[0] = MATAIJ;
33c4762a1bSJed Brown   }
34c4762a1bSJed Brown   if (size == 1) {
35c4762a1bSJed Brown     ntypes  = 3;
36c4762a1bSJed Brown     type[1] = MATSEQBAIJ;
37c4762a1bSJed Brown     type[2] = MATSEQSBAIJ;
38c4762a1bSJed Brown   } else {
39c4762a1bSJed Brown     ntypes  = 3;
40c4762a1bSJed Brown     type[1] = MATMPIBAIJ;
41c4762a1bSJed Brown     type[2] = MATMPISBAIJ;
42c4762a1bSJed Brown   }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* input matrix C */
45c4762a1bSJed Brown   if (flg_loadmat) {
46c4762a1bSJed Brown     /* Open binary file. */
475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown     /* Load a baij matrix, then destroy the viewer. */
505f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
51c4762a1bSJed Brown     if (size == 1) {
525f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetType(C,MATSEQBAIJ));
53c4762a1bSJed Brown     } else {
545f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetType(C,MATMPIBAIJ));
55c4762a1bSJed Brown     }
565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(C));
575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLoad(C,fd));
585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&fd));
59c4762a1bSJed Brown   } else { /* Create a baij mat with bs>1  */
60c4762a1bSJed Brown     bs   = 2; mbs=8;
615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
632c71b3e2SJacob Faibussowitsch     PetscCheckFalse(bs <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," bs must be >1 in this case");
64c4762a1bSJed Brown     m    = mbs*bs;
655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,m,m,d_nz,NULL,o_nz,NULL,&C));
66c4762a1bSJed Brown     for (block=0; block<mbs; block++) {
67c4762a1bSJed Brown       /* diagonal blocks */
68c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
69c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
70c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
715f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES));
72c4762a1bSJed Brown       }
73c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
74c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
755f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
78c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
795f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES));
80c4762a1bSJed Brown     }
81c4762a1bSJed Brown     /* off-diagonal blocks */
82c4762a1bSJed Brown     value[0]=-1.0;
83c4762a1bSJed Brown     for (i=0; i<(mbs-1)*bs; i++) {
84c4762a1bSJed Brown       col[0]=i+bs;
855f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,1,col,value,INSERT_VALUES));
86c4762a1bSJed Brown       col[0]=i; row=i+bs;
875f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&row,1,col,value,INSERT_VALUES));
88c4762a1bSJed Brown     }
895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   {
93c4762a1bSJed Brown     /* Check the symmetry of C because it will be converted to a sbaij matrix */
94c4762a1bSJed Brown     Mat Ctrans;
955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&Ctrans));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(MatEqual(C,Ctrans,&flg));
97c4762a1bSJed Brown /*
98c4762a1bSJed Brown     {
995f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAXPY(C,1.,Ctrans,DIFFERENT_NONZERO_PATTERN));
100c4762a1bSJed Brown       flg  = PETSC_TRUE;
101c4762a1bSJed Brown     }
102c4762a1bSJed Brown */
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(C,MAT_SYMMETRIC,flg));
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Ctrans));
105c4762a1bSJed Brown   }
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIsSymmetric(C,0.0,&issymmetric));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(C,NULL,"-view_mat"));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* convert C to other formats */
110c4762a1bSJed Brown   for (i=0; i<ntypes; i++) {
111c4762a1bSJed Brown     PetscBool ismpisbaij,isseqsbaij;
112c4762a1bSJed Brown 
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(type[i],MATMPISBAIJ,&ismpisbaij));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(type[i],MATMPISBAIJ,&isseqsbaij));
115c4762a1bSJed Brown     if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(C,type[i],MAT_INITIAL_MATRIX,&A));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,C,10,&equal));
11828b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from BAIJ to %s",type[i]);
119c4762a1bSJed Brown     for (j=i+1; j<ntypes; j++) {
1205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij));
1215f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij));
122c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
123c4762a1bSJed Brown       if (verbose>0) {
1245f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," \n[%d] test conversion between %s and %s\n",rank,type[i],type[j]));
125c4762a1bSJed Brown       }
126c4762a1bSJed Brown 
127dd400576SPatrick Sanan       if (rank == 0 && verbose) printf("Convert %s A to %s B\n",type[i],type[j]);
1285f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,type[j],MAT_INITIAL_MATRIX,&B));
129c4762a1bSJed Brown       /*
130c4762a1bSJed Brown       if (j == 2) {
1315f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF," A: %s\n",type[i]));
1325f80ce2aSJacob Faibussowitsch         CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
1335f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF," B: %s\n",type[j]));
1345f80ce2aSJacob Faibussowitsch         CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
135c4762a1bSJed Brown       }
136c4762a1bSJed Brown        */
1375f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultEqual(A,B,10,&equal));
13828b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[i],type[j]);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown       if (size == 1 || j != 2) { /* Matconvert from mpisbaij mat to other formats are not supported */
141dd400576SPatrick Sanan         if (rank == 0 && verbose) printf("Convert %s B to %s D\n",type[j],type[i]);
1425f80ce2aSJacob Faibussowitsch         CHKERRQ(MatConvert(B,type[i],MAT_INITIAL_MATRIX,&D));
1435f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultEqual(B,D,10,&equal));
14428b400f6SJacob Faibussowitsch         PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[j],type[i]);
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch         CHKERRQ(MatDestroy(&D));
147c4762a1bSJed Brown       }
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B));
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&D));
150c4762a1bSJed Brown     }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     /* Test in-place convert */
153c4762a1bSJed Brown     if (size == 1) { /* size > 1 is not working yet! */
154c4762a1bSJed Brown       j = (i+1)%ntypes;
1555f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij));
1565f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij));
157c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
158c4762a1bSJed Brown       /* printf("[%d] i: %d, j: %d\n",rank,i,j); */
1595f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,type[j],MAT_INPLACE_MATRIX,&A));
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown 
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /* test BAIJ to MATIS */
166c4762a1bSJed Brown   if (size > 1) {
167c4762a1bSJed Brown     MatType ctype;
168c4762a1bSJed Brown 
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetType(C,&ctype));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(C,MATIS,MAT_INITIAL_MATRIX,&A));
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,C,10,&equal));
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(A,NULL,"-view_conv"));
173c4762a1bSJed Brown     if (!equal) {
174c4762a1bSJed Brown       Mat C2;
175c4762a1bSJed Brown 
1765f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
1775f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
1785f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
1795f80ce2aSJacob Faibussowitsch       CHKERRQ(MatChop(C2,PETSC_SMALL));
1805f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_err"));
1815f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&C2));
182c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
183c4762a1bSJed Brown     }
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(C,MATIS,MAT_REUSE_MATRIX,&A));
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,C,10,&equal));
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(A,NULL,"-view_conv"));
187c4762a1bSJed Brown     if (!equal) {
188c4762a1bSJed Brown       Mat C2;
189c4762a1bSJed Brown 
1905f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
1915f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
1925f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
1935f80ce2aSJacob Faibussowitsch       CHKERRQ(MatChop(C2,PETSC_SMALL));
1945f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_err"));
1955f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&C2));
196c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
197c4762a1bSJed Brown     }
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
1995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&A));
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATIS,MAT_INPLACE_MATRIX,&A));
2015f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(A,NULL,"-view_conv"));
2025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,C,10,&equal));
203c4762a1bSJed Brown     if (!equal) {
204c4762a1bSJed Brown       Mat C2;
205c4762a1bSJed Brown 
2065f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(A,NULL,"-view_conv"));
2075f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2));
2085f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_conv_assembled"));
2095f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN));
2105f80ce2aSJacob Faibussowitsch       CHKERRQ(MatChop(C2,PETSC_SMALL));
2115f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(C2,NULL,"-view_err"));
2125f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&C2));
213c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
214c4762a1bSJed Brown     }
2155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
216c4762a1bSJed Brown   }
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
218c4762a1bSJed Brown 
219*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
220*b122ec5aSJacob Faibussowitsch   return 0;
221c4762a1bSJed Brown }
222c4762a1bSJed Brown 
223c4762a1bSJed Brown /*TEST
224c4762a1bSJed Brown 
225c4762a1bSJed Brown    test:
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       suffix: 2
229c4762a1bSJed Brown       nsize: 3
230c4762a1bSJed Brown 
231c4762a1bSJed Brown    testset:
232c4762a1bSJed Brown       requires: parmetis
233c4762a1bSJed Brown       output_file: output/ex55_1.out
234c4762a1bSJed Brown       nsize: 3
235c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis
236c4762a1bSJed Brown       test:
237c4762a1bSJed Brown         suffix: matis_baij_parmetis_nd
238c4762a1bSJed Brown       test:
239c4762a1bSJed Brown         suffix: matis_aij_parmetis_nd
240c4762a1bSJed Brown         args: -testseqaij
241c4762a1bSJed Brown       test:
242dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
243c4762a1bSJed Brown         suffix: matis_poisson1_parmetis_nd
244c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
245c4762a1bSJed Brown 
246c4762a1bSJed Brown    testset:
247dfd57a17SPierre Jolivet       requires: ptscotch defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
248c4762a1bSJed Brown       output_file: output/ex55_1.out
249c4762a1bSJed Brown       nsize: 4
250c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch
251c4762a1bSJed Brown       test:
252c4762a1bSJed Brown         suffix: matis_baij_ptscotch_nd
253c4762a1bSJed Brown       test:
254c4762a1bSJed Brown         suffix: matis_aij_ptscotch_nd
255c4762a1bSJed Brown         args: -testseqaij
256c4762a1bSJed Brown       test:
257dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
258c4762a1bSJed Brown         suffix: matis_poisson1_ptscotch_nd
259c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
260c4762a1bSJed Brown 
261c4762a1bSJed Brown TEST*/
262