xref: /petsc/src/mat/tests/ex55.c (revision a289a2f2dff3362f6bfc0c874680f4e4d851e010)
1c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown /* Usage: mpiexec -n <np> ex55 -verbose <0 or 1> */
5c4762a1bSJed Brown 
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C, A, B, D;
9c4762a1bSJed Brown   PetscInt    i, j, ntypes, bs, mbs, m, block, d_nz = 6, o_nz = 3, col[3], row, verbose = 0;
10c4762a1bSJed Brown   PetscMPIInt size, rank;
11c4762a1bSJed Brown   MatType     type[9];
12c4762a1bSJed Brown   char        file[PETSC_MAX_PATH_LEN];
13c4762a1bSJed Brown   PetscViewer fd;
14c4762a1bSJed Brown   PetscBool   equal, flg_loadmat, flg, issymmetric;
15c4762a1bSJed Brown   PetscScalar value[3];
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-verbose", &verbose, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg_loadmat));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCall(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. */
479566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown     /* Load a baij matrix, then destroy the viewer. */
509566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
51c4762a1bSJed Brown     if (size == 1) {
529566063dSJacob Faibussowitsch       PetscCall(MatSetType(C, MATSEQBAIJ));
53c4762a1bSJed Brown     } else {
549566063dSJacob Faibussowitsch       PetscCall(MatSetType(C, MATMPIBAIJ));
55c4762a1bSJed Brown     }
569566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(C));
579566063dSJacob Faibussowitsch     PetscCall(MatLoad(C, fd));
589566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
59c4762a1bSJed Brown   } else { /* Create a baij mat with bs>1  */
609371c9d4SSatish Balay     bs  = 2;
619371c9d4SSatish Balay     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 */
699371c9d4SSatish Balay       value[0] = -1.0;
709371c9d4SSatish Balay       value[1] = 4.0;
719371c9d4SSatish Balay       value[2] = -1.0;
72c4762a1bSJed Brown       for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) {
739371c9d4SSatish Balay         col[0] = i - 1;
749371c9d4SSatish Balay         col[1] = i;
759371c9d4SSatish Balay         col[2] = i + 1;
769566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C, 1, &i, 3, col, value, INSERT_VALUES));
77c4762a1bSJed Brown       }
789371c9d4SSatish Balay       i        = bs - 1 + block * bs;
799371c9d4SSatish Balay       col[0]   = bs - 2 + block * bs;
809371c9d4SSatish Balay       col[1]   = bs - 1 + block * bs;
819371c9d4SSatish Balay       value[0] = -1.0;
829371c9d4SSatish Balay       value[1] = 4.0;
839566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &i, 2, col, value, INSERT_VALUES));
84c4762a1bSJed Brown 
859371c9d4SSatish Balay       i        = 0 + block * bs;
869371c9d4SSatish Balay       col[0]   = 0 + block * bs;
879371c9d4SSatish Balay       col[1]   = 1 + block * bs;
889371c9d4SSatish Balay       value[0] = 4.0;
899371c9d4SSatish Balay       value[1] = -1.0;
909566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &i, 2, col, value, INSERT_VALUES));
91c4762a1bSJed Brown     }
92c4762a1bSJed Brown     /* off-diagonal blocks */
93c4762a1bSJed Brown     value[0] = -1.0;
94c4762a1bSJed Brown     for (i = 0; i < (mbs - 1) * bs; i++) {
95c4762a1bSJed Brown       col[0] = i + bs;
969566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &i, 1, col, value, INSERT_VALUES));
979371c9d4SSatish Balay       col[0] = i;
989371c9d4SSatish Balay       row    = i + bs;
999566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &row, 1, col, value, INSERT_VALUES));
100c4762a1bSJed Brown     }
1019566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   {
105c4762a1bSJed Brown     /* Check the symmetry of C because it will be converted to a sbaij matrix */
106c4762a1bSJed Brown     Mat Ctrans;
1079566063dSJacob Faibussowitsch     PetscCall(MatTranspose(C, MAT_INITIAL_MATRIX, &Ctrans));
1089566063dSJacob Faibussowitsch     PetscCall(MatEqual(C, Ctrans, &flg));
109c4762a1bSJed Brown     /*
110c4762a1bSJed Brown     {
1119566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C,1.,Ctrans,DIFFERENT_NONZERO_PATTERN));
112c4762a1bSJed Brown       flg  = PETSC_TRUE;
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown */
1159566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_SYMMETRIC, flg));
1169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Ctrans));
117c4762a1bSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(C, 0.0, &issymmetric));
1199566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(C, NULL, "-view_mat"));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* convert C to other formats */
122c4762a1bSJed Brown   for (i = 0; i < ntypes; i++) {
123c4762a1bSJed Brown     PetscBool ismpisbaij, isseqsbaij;
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(type[i], MATMPISBAIJ, &ismpisbaij));
126*2b2f8cc6SPierre Jolivet     PetscCall(PetscStrcmp(type[i], MATSEQSBAIJ, &isseqsbaij));
127c4762a1bSJed Brown     if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
1289566063dSJacob Faibussowitsch     PetscCall(MatConvert(C, type[i], MAT_INITIAL_MATRIX, &A));
1299566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, C, 10, &equal));
13028b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in conversion from BAIJ to %s", type[i]);
131c4762a1bSJed Brown     for (j = i + 1; j < ntypes; j++) {
1329566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j], MATMPISBAIJ, &ismpisbaij));
133*2b2f8cc6SPierre Jolivet       PetscCall(PetscStrcmp(type[j], MATSEQSBAIJ, &isseqsbaij));
134c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
13548a46eb9SPierre Jolivet       if (verbose > 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " \n[%d] test conversion between %s and %s\n", rank, type[i], type[j]));
136c4762a1bSJed Brown 
137dd400576SPatrick Sanan       if (rank == 0 && verbose) printf("Convert %s A to %s B\n", type[i], type[j]);
1389566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, type[j], MAT_INITIAL_MATRIX, &B));
139c4762a1bSJed Brown       /*
140c4762a1bSJed Brown       if (j == 2) {
1419566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF," A: %s\n",type[i]));
1429566063dSJacob Faibussowitsch         PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
1439566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF," B: %s\n",type[j]));
1449566063dSJacob Faibussowitsch         PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
145c4762a1bSJed Brown       }
146c4762a1bSJed Brown        */
1479566063dSJacob Faibussowitsch       PetscCall(MatMultEqual(A, B, 10, &equal));
14828b400f6SJacob Faibussowitsch       PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in conversion from %s to %s", type[i], type[j]);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown       if (size == 1 || j != 2) { /* Matconvert from mpisbaij mat to other formats are not supported */
151dd400576SPatrick Sanan         if (rank == 0 && verbose) printf("Convert %s B to %s D\n", type[j], type[i]);
1529566063dSJacob Faibussowitsch         PetscCall(MatConvert(B, type[i], MAT_INITIAL_MATRIX, &D));
1539566063dSJacob Faibussowitsch         PetscCall(MatMultEqual(B, D, 10, &equal));
15428b400f6SJacob Faibussowitsch         PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in conversion from %s to %s", type[j], type[i]);
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&D));
157c4762a1bSJed Brown       }
1589566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
1599566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&D));
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     /* Test in-place convert */
163c4762a1bSJed Brown     if (size == 1) { /* size > 1 is not working yet! */
164c4762a1bSJed Brown       j = (i + 1) % ntypes;
1659566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(type[j], MATMPISBAIJ, &ismpisbaij));
166*2b2f8cc6SPierre Jolivet       PetscCall(PetscStrcmp(type[j], MATSEQSBAIJ, &isseqsbaij));
167c4762a1bSJed Brown       if (!issymmetric && (ismpisbaij || isseqsbaij)) continue;
168c4762a1bSJed Brown       /* printf("[%d] i: %d, j: %d\n",rank,i,j); */
1699566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, type[j], MAT_INPLACE_MATRIX, &A));
170c4762a1bSJed Brown     }
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* test BAIJ to MATIS */
176c4762a1bSJed Brown   if (size > 1) {
177c4762a1bSJed Brown     MatType ctype;
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch     PetscCall(MatGetType(C, &ctype));
1809566063dSJacob Faibussowitsch     PetscCall(MatConvert(C, MATIS, MAT_INITIAL_MATRIX, &A));
1819566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, C, 10, &equal));
1829566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A, NULL, "-view_conv"));
183c4762a1bSJed Brown     if (!equal) {
184c4762a1bSJed Brown       Mat C2;
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, ctype, MAT_INITIAL_MATRIX, &C2));
1879566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_conv_assembled"));
1889566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2, -1., C, DIFFERENT_NONZERO_PATTERN));
1892ce66baaSPierre Jolivet       PetscCall(MatFilter(C2, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
1909566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_err"));
1919566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
192c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in conversion from BAIJ to MATIS");
193c4762a1bSJed Brown     }
1949566063dSJacob Faibussowitsch     PetscCall(MatConvert(C, MATIS, MAT_REUSE_MATRIX, &A));
1959566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, C, 10, &equal));
1969566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A, NULL, "-view_conv"));
197c4762a1bSJed Brown     if (!equal) {
198c4762a1bSJed Brown       Mat C2;
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, ctype, MAT_INITIAL_MATRIX, &C2));
2019566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_conv_assembled"));
2029566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2, -1., C, DIFFERENT_NONZERO_PATTERN));
2032ce66baaSPierre Jolivet       PetscCall(MatFilter(C2, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
2049566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_err"));
2059566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
206c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in conversion from BAIJ to MATIS");
207c4762a1bSJed Brown     }
2089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
2099566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &A));
2109566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATIS, MAT_INPLACE_MATRIX, &A));
2119566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(A, NULL, "-view_conv"));
2129566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, C, 10, &equal));
213c4762a1bSJed Brown     if (!equal) {
214c4762a1bSJed Brown       Mat C2;
215c4762a1bSJed Brown 
2169566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(A, NULL, "-view_conv"));
2179566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, ctype, MAT_INITIAL_MATRIX, &C2));
2189566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_conv_assembled"));
2199566063dSJacob Faibussowitsch       PetscCall(MatAXPY(C2, -1., C, DIFFERENT_NONZERO_PATTERN));
2202ce66baaSPierre Jolivet       PetscCall(MatFilter(C2, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
2219566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(C2, NULL, "-view_err"));
2229566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&C2));
223c4762a1bSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in conversion from BAIJ to MATIS");
224c4762a1bSJed Brown     }
2259566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
226c4762a1bSJed Brown   }
2279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
230b122ec5aSJacob Faibussowitsch   return 0;
231c4762a1bSJed Brown }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown /*TEST
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    test:
2363886731fSPierre Jolivet       output_file: output/empty.out
237c4762a1bSJed Brown 
238c4762a1bSJed Brown    test:
239c4762a1bSJed Brown       suffix: 2
240c4762a1bSJed Brown       nsize: 3
2413886731fSPierre Jolivet       output_file: output/empty.out
242c4762a1bSJed Brown 
243c4762a1bSJed Brown    testset:
244c4762a1bSJed Brown       requires: parmetis
2453886731fSPierre Jolivet       output_file: output/empty.out
246c4762a1bSJed Brown       nsize: 3
247c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis
248c4762a1bSJed Brown       test:
249c4762a1bSJed Brown         suffix: matis_baij_parmetis_nd
250c4762a1bSJed Brown       test:
251c4762a1bSJed Brown         suffix: matis_aij_parmetis_nd
252c4762a1bSJed Brown         args: -testseqaij
253c4762a1bSJed Brown       test:
254dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
255c4762a1bSJed Brown         suffix: matis_poisson1_parmetis_nd
256c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    testset:
259dfd57a17SPierre Jolivet       requires: ptscotch defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
2603886731fSPierre Jolivet       output_file: output/empty.out
261c4762a1bSJed Brown       nsize: 4
262c4762a1bSJed Brown       args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch
263c4762a1bSJed Brown       test:
264c4762a1bSJed Brown         suffix: matis_baij_ptscotch_nd
265c4762a1bSJed Brown       test:
266c4762a1bSJed Brown         suffix: matis_aij_ptscotch_nd
267c4762a1bSJed Brown         args: -testseqaij
268c4762a1bSJed Brown       test:
269dfd57a17SPierre Jolivet         requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
270c4762a1bSJed Brown         suffix: matis_poisson1_ptscotch_nd
271c4762a1bSJed Brown         args: -f ${DATAFILESPATH}/matrices/poisson1
272c4762a1bSJed Brown 
273c4762a1bSJed Brown TEST*/
274