xref: /petsc/src/mat/tests/ex55.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown /* Usage: mpiexec -n <np> ex55 -verbose <0 or 1> */
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc,char **args)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   Mat            C,A,B,D;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscInt       i,j,ntypes,bs,mbs,m,block,d_nz=6, o_nz=3,col[3],row,verbose=0;
12*c4762a1bSJed Brown   PetscMPIInt    size,rank;
13*c4762a1bSJed Brown   MatType        type[9];
14*c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
15*c4762a1bSJed Brown   PetscViewer    fd;
16*c4762a1bSJed Brown   PetscBool      equal,flg_loadmat,flg,issymmetric;
17*c4762a1bSJed Brown   PetscScalar    value[3];
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg_loadmat);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-testseqaij",&flg);CHKERRQ(ierr);
26*c4762a1bSJed Brown   if (flg) {
27*c4762a1bSJed Brown     if (size == 1) {
28*c4762a1bSJed Brown       type[0] = MATSEQAIJ;
29*c4762a1bSJed Brown     } else {
30*c4762a1bSJed Brown       type[0] = MATMPIAIJ;
31*c4762a1bSJed Brown     }
32*c4762a1bSJed Brown   } else {
33*c4762a1bSJed Brown     type[0] = MATAIJ;
34*c4762a1bSJed Brown   }
35*c4762a1bSJed Brown   if (size == 1) {
36*c4762a1bSJed Brown     ntypes  = 3;
37*c4762a1bSJed Brown     type[1] = MATSEQBAIJ;
38*c4762a1bSJed Brown     type[2] = MATSEQSBAIJ;
39*c4762a1bSJed Brown   } else {
40*c4762a1bSJed Brown     ntypes  = 3;
41*c4762a1bSJed Brown     type[1] = MATMPIBAIJ;
42*c4762a1bSJed Brown     type[2] = MATMPISBAIJ;
43*c4762a1bSJed Brown   }
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown   /* input matrix C */
46*c4762a1bSJed Brown   if (flg_loadmat) {
47*c4762a1bSJed Brown     /* Open binary file. */
48*c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown     /* Load a baij matrix, then destroy the viewer. */
51*c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
52*c4762a1bSJed Brown     if (size == 1) {
53*c4762a1bSJed Brown       ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
54*c4762a1bSJed Brown     } else {
55*c4762a1bSJed Brown       ierr = MatSetType(C,MATMPIBAIJ);CHKERRQ(ierr);
56*c4762a1bSJed Brown     }
57*c4762a1bSJed Brown     ierr = MatSetFromOptions(C);CHKERRQ(ierr);
58*c4762a1bSJed Brown     ierr = MatLoad(C,fd);CHKERRQ(ierr);
59*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
60*c4762a1bSJed Brown   } else { /* Create a baij mat with bs>1  */
61*c4762a1bSJed Brown     bs   = 2; mbs=8;
62*c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr);
63*c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
64*c4762a1bSJed Brown     if (bs <= 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," bs must be >1 in this case");
65*c4762a1bSJed Brown     m    = mbs*bs;
66*c4762a1bSJed Brown     ierr = MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,m,m,d_nz,NULL,o_nz,NULL,&C);CHKERRQ(ierr);
67*c4762a1bSJed Brown     for (block=0; block<mbs; block++) {
68*c4762a1bSJed Brown       /* diagonal blocks */
69*c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
70*c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
71*c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
72*c4762a1bSJed Brown         ierr   = MatSetValues(C,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
73*c4762a1bSJed Brown       }
74*c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
75*c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
76*c4762a1bSJed Brown       ierr    = MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
79*c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
80*c4762a1bSJed Brown       ierr    = MatSetValues(C,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
81*c4762a1bSJed Brown     }
82*c4762a1bSJed Brown     /* off-diagonal blocks */
83*c4762a1bSJed Brown     value[0]=-1.0;
84*c4762a1bSJed Brown     for (i=0; i<(mbs-1)*bs; i++) {
85*c4762a1bSJed Brown       col[0]=i+bs;
86*c4762a1bSJed Brown       ierr  = MatSetValues(C,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
87*c4762a1bSJed Brown       col[0]=i; row=i+bs;
88*c4762a1bSJed Brown       ierr  = MatSetValues(C,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
89*c4762a1bSJed Brown     }
90*c4762a1bSJed Brown     ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91*c4762a1bSJed Brown     ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92*c4762a1bSJed Brown   }
93*c4762a1bSJed Brown   {
94*c4762a1bSJed Brown     /* Check the symmetry of C because it will be converted to a sbaij matrix */
95*c4762a1bSJed Brown     Mat Ctrans;
96*c4762a1bSJed Brown     ierr = MatTranspose(C,MAT_INITIAL_MATRIX,&Ctrans);CHKERRQ(ierr);
97*c4762a1bSJed Brown     ierr = MatEqual(C,Ctrans,&flg);CHKERRQ(ierr);
98*c4762a1bSJed Brown /*
99*c4762a1bSJed Brown     {
100*c4762a1bSJed Brown       ierr = MatAXPY(C,1.,Ctrans,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
101*c4762a1bSJed Brown       flg  = PETSC_TRUE;
102*c4762a1bSJed Brown     }
103*c4762a1bSJed Brown */
104*c4762a1bSJed Brown     ierr = MatSetOption(C,MAT_SYMMETRIC,flg);CHKERRQ(ierr);
105*c4762a1bSJed Brown     ierr = MatDestroy(&Ctrans);CHKERRQ(ierr);
106*c4762a1bSJed Brown   }
107*c4762a1bSJed Brown   ierr = MatIsSymmetric(C,0.0,&issymmetric);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatViewFromOptions(C,NULL,"-view_mat");CHKERRQ(ierr);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   /* convert C to other formats */
111*c4762a1bSJed Brown   for (i=0; i<ntypes; i++) {
112*c4762a1bSJed Brown     PetscBool ismpisbaij,isseqsbaij;
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown     ierr = PetscStrcmp(type[i],MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
115*c4762a1bSJed Brown     ierr = PetscStrcmp(type[i],MATMPISBAIJ,&isseqsbaij);CHKERRQ(ierr);
116*c4762a1bSJed Brown     if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
117*c4762a1bSJed Brown     ierr = MatConvert(C,type[i],MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
118*c4762a1bSJed Brown     ierr = MatMultEqual(A,C,10,&equal);CHKERRQ(ierr);
119*c4762a1bSJed Brown     if (!equal) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from BAIJ to %s",type[i]);
120*c4762a1bSJed Brown     for (j=i+1; j<ntypes; j++) {
121*c4762a1bSJed Brown       ierr = PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
122*c4762a1bSJed Brown       ierr = PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij);CHKERRQ(ierr);
123*c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
124*c4762a1bSJed Brown       if (verbose>0) {
125*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD," \n[%d] test conversion between %s and %s\n",rank,type[i],type[j]);CHKERRQ(ierr);
126*c4762a1bSJed Brown       }
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown       if (!rank && verbose) printf("Convert %s A to %s B\n",type[i],type[j]);
129*c4762a1bSJed Brown       ierr = MatConvert(A,type[j],MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
130*c4762a1bSJed Brown       /*
131*c4762a1bSJed Brown       if (j == 2) {
132*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," A: %s\n",type[i]);CHKERRQ(ierr);
133*c4762a1bSJed Brown         ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
134*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," B: %s\n",type[j]);CHKERRQ(ierr);
135*c4762a1bSJed Brown         ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
136*c4762a1bSJed Brown       }
137*c4762a1bSJed Brown        */
138*c4762a1bSJed Brown       ierr = MatMultEqual(A,B,10,&equal);CHKERRQ(ierr);
139*c4762a1bSJed Brown       if (!equal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[i],type[j]);
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown       if (size == 1 || j != 2) { /* Matconvert from mpisbaij mat to other formats are not supported */
142*c4762a1bSJed Brown         if (!rank && verbose) printf("Convert %s B to %s D\n",type[j],type[i]);
143*c4762a1bSJed Brown         ierr = MatConvert(B,type[i],MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr);
144*c4762a1bSJed Brown         ierr = MatMultEqual(B,D,10,&equal);CHKERRQ(ierr);
145*c4762a1bSJed Brown         if (!equal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[j],type[i]);
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown         ierr = MatDestroy(&D);CHKERRQ(ierr);
148*c4762a1bSJed Brown       }
149*c4762a1bSJed Brown       ierr = MatDestroy(&B);CHKERRQ(ierr);
150*c4762a1bSJed Brown       ierr = MatDestroy(&D);CHKERRQ(ierr);
151*c4762a1bSJed Brown     }
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown     /* Test in-place convert */
154*c4762a1bSJed Brown     if (size == 1) { /* size > 1 is not working yet! */
155*c4762a1bSJed Brown       j = (i+1)%ntypes;
156*c4762a1bSJed Brown       ierr = PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij);CHKERRQ(ierr);
157*c4762a1bSJed Brown       ierr = PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij);CHKERRQ(ierr);
158*c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
159*c4762a1bSJed Brown       /* printf("[%d] i: %d, j: %d\n",rank,i,j); */
160*c4762a1bSJed Brown       ierr = MatConvert(A,type[j],MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
161*c4762a1bSJed Brown     }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
164*c4762a1bSJed Brown   }
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   /* test BAIJ to MATIS */
167*c4762a1bSJed Brown   if (size > 1) {
168*c4762a1bSJed Brown     MatType ctype;
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown     ierr = MatGetType(C,&ctype);CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr = MatConvert(C,MATIS,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
172*c4762a1bSJed Brown     ierr = MatMultEqual(A,C,10,&equal);CHKERRQ(ierr);
173*c4762a1bSJed Brown     ierr = MatViewFromOptions(A,NULL,"-view_conv");CHKERRQ(ierr);
174*c4762a1bSJed Brown     if (!equal) {
175*c4762a1bSJed Brown       Mat C2;
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown       ierr = MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2);CHKERRQ(ierr);
178*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_conv_assembled");CHKERRQ(ierr);
179*c4762a1bSJed Brown       ierr = MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
180*c4762a1bSJed Brown       ierr = MatChop(C2,PETSC_SMALL);CHKERRQ(ierr);
181*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_err");CHKERRQ(ierr);
182*c4762a1bSJed Brown       ierr = MatDestroy(&C2);CHKERRQ(ierr);
183*c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
184*c4762a1bSJed Brown     }
185*c4762a1bSJed Brown     ierr = MatConvert(C,MATIS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
186*c4762a1bSJed Brown     ierr = MatMultEqual(A,C,10,&equal);CHKERRQ(ierr);
187*c4762a1bSJed Brown     ierr = MatViewFromOptions(A,NULL,"-view_conv");CHKERRQ(ierr);
188*c4762a1bSJed Brown     if (!equal) {
189*c4762a1bSJed Brown       Mat C2;
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown       ierr = MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2);CHKERRQ(ierr);
192*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_conv_assembled");CHKERRQ(ierr);
193*c4762a1bSJed Brown       ierr = MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
194*c4762a1bSJed Brown       ierr = MatChop(C2,PETSC_SMALL);CHKERRQ(ierr);
195*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_err");CHKERRQ(ierr);
196*c4762a1bSJed Brown       ierr = MatDestroy(&C2);CHKERRQ(ierr);
197*c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
198*c4762a1bSJed Brown     }
199*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
200*c4762a1bSJed Brown     ierr = MatDuplicate(C,MAT_COPY_VALUES,&A);CHKERRQ(ierr);
201*c4762a1bSJed Brown     ierr = MatConvert(A,MATIS,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
202*c4762a1bSJed Brown     ierr = MatViewFromOptions(A,NULL,"-view_conv");CHKERRQ(ierr);
203*c4762a1bSJed Brown     ierr = MatMultEqual(A,C,10,&equal);CHKERRQ(ierr);
204*c4762a1bSJed Brown     if (!equal) {
205*c4762a1bSJed Brown       Mat C2;
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown       ierr = MatViewFromOptions(A,NULL,"-view_conv");CHKERRQ(ierr);
208*c4762a1bSJed Brown       ierr = MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2);CHKERRQ(ierr);
209*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_conv_assembled");CHKERRQ(ierr);
210*c4762a1bSJed Brown       ierr = MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
211*c4762a1bSJed Brown       ierr = MatChop(C2,PETSC_SMALL);CHKERRQ(ierr);
212*c4762a1bSJed Brown       ierr = MatViewFromOptions(C2,NULL,"-view_err");CHKERRQ(ierr);
213*c4762a1bSJed Brown       ierr = MatDestroy(&C2);CHKERRQ(ierr);
214*c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS");
215*c4762a1bSJed Brown     }
216*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
217*c4762a1bSJed Brown   }
218*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown   ierr = PetscFinalize();
221*c4762a1bSJed Brown   return ierr;
222*c4762a1bSJed Brown }
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown /*TEST
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown    test:
227*c4762a1bSJed Brown 
228*c4762a1bSJed Brown    test:
229*c4762a1bSJed Brown       suffix: 2
230*c4762a1bSJed Brown       nsize: 3
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown    testset:
233*c4762a1bSJed Brown       requires: parmetis
234*c4762a1bSJed Brown       output_file: output/ex55_1.out
235*c4762a1bSJed Brown       nsize: 3
236*c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis
237*c4762a1bSJed Brown       test:
238*c4762a1bSJed Brown         suffix: matis_baij_parmetis_nd
239*c4762a1bSJed Brown       test:
240*c4762a1bSJed Brown         suffix: matis_aij_parmetis_nd
241*c4762a1bSJed Brown         args: -testseqaij
242*c4762a1bSJed Brown       test:
243*c4762a1bSJed Brown         requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
244*c4762a1bSJed Brown         suffix: matis_poisson1_parmetis_nd
245*c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown    testset:
248*c4762a1bSJed Brown       requires: ptscotch define(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
249*c4762a1bSJed Brown       output_file: output/ex55_1.out
250*c4762a1bSJed Brown       nsize: 4
251*c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch
252*c4762a1bSJed Brown       test:
253*c4762a1bSJed Brown         suffix: matis_baij_ptscotch_nd
254*c4762a1bSJed Brown       test:
255*c4762a1bSJed Brown         suffix: matis_aij_ptscotch_nd
256*c4762a1bSJed Brown         args: -testseqaij
257*c4762a1bSJed Brown       test:
258*c4762a1bSJed Brown         requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
259*c4762a1bSJed Brown         suffix: matis_poisson1_ptscotch_nd
260*c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown TEST*/
263