xref: /petsc/src/mat/tests/ex109.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c20d7725SJed Brown   Mat         A, B, C, D, AT;
8c4762a1bSJed Brown   PetscInt    i, M, N, Istart, Iend, n = 7, j, J, Ii, m = 8, am, an;
9c4762a1bSJed Brown   PetscScalar v;
10c4762a1bSJed Brown   PetscRandom r;
113fbca975SStefano Zampini   PetscBool   equal = PETSC_FALSE, flg;
12c20d7725SJed Brown   PetscReal   fill  = 1.0, norm;
13c4762a1bSJed Brown   PetscMPIInt size;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
16c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL));
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r));
229566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create a aij matrix A */
25c4762a1bSJed Brown   M = N = m * n;
269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
289566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
309566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL));
319566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
34c4762a1bSJed Brown   am = Iend - Istart;
35c4762a1bSJed Brown   for (Ii = Istart; Ii < Iend; Ii++) {
369371c9d4SSatish Balay     v = -1.0;
379371c9d4SSatish Balay     i = Ii / n;
389371c9d4SSatish Balay     j = Ii - i * n;
399371c9d4SSatish Balay     if (i > 0) {
409371c9d4SSatish Balay       J = Ii - n;
419371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
429371c9d4SSatish Balay     }
439371c9d4SSatish Balay     if (i < m - 1) {
449371c9d4SSatish Balay       J = Ii + n;
459371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
469371c9d4SSatish Balay     }
479371c9d4SSatish Balay     if (j > 0) {
489371c9d4SSatish Balay       J = Ii - 1;
499371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
509371c9d4SSatish Balay     }
519371c9d4SSatish Balay     if (j < n - 1) {
529371c9d4SSatish Balay       J = Ii + 1;
539371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
549371c9d4SSatish Balay     }
559371c9d4SSatish Balay     v = 4.0;
569371c9d4SSatish Balay     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
57c4762a1bSJed Brown   }
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Create a dense matrix B */
629566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
639566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
649566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, an, PETSC_DECIDE, PETSC_DECIDE, M));
659566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATDENSE));
669566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(B, NULL));
679566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(B, NULL));
689566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
699566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B, r));
709566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
71c4762a1bSJed Brown 
72c20d7725SJed Brown   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
749566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
759566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, am, PETSC_DECIDE, PETSC_DECIDE, N));
769566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
779566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
789566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
799566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C));
809566063dSJacob Faibussowitsch   PetscCall(MatNorm(C, NORM_INFINITY, &norm));
819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
82c20d7725SJed Brown 
83c4762a1bSJed Brown   /* Test C = A*B (aij*dense) */
849566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C));
859566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C));
86c4762a1bSJed Brown 
87c20d7725SJed Brown   /* Test developer API */
889566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(A, B, NULL, &D));
899566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(D, MATPRODUCT_AB));
909566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(D, "default"));
919566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(D, fill));
929566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(D));
939566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(D));
9448a46eb9SPierre Jolivet   for (i = 0; i < 2; i++) PetscCall(MatProductNumeric(D));
959566063dSJacob Faibussowitsch   PetscCall(MatEqual(C, D, &equal));
9628b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "C != D");
979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
98c4762a1bSJed Brown 
99c20d7725SJed Brown   /* Test D = AT*B (transpose(aij)*dense) */
1009566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(A, &AT));
1019566063dSJacob Faibussowitsch   PetscCall(MatMatMult(AT, B, MAT_INITIAL_MATRIX, fill, &D));
1029566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(AT, B, D, 10, &equal));
10328b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != AT*B (transpose(aij)*dense)");
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AT));
106c20d7725SJed Brown 
107c4762a1bSJed Brown   /* Test D = C*A (dense*aij) */
1089566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &D));
1099566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C, A, MAT_REUSE_MATRIX, fill, &D));
1109566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(C, A, D, 10, &equal));
11128b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != C*A (dense*aij)");
1129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Test D = A*C (aij*dense) */
1159566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, fill, &D));
1169566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, fill, &D));
1179566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A, C, D, 10, &equal));
11828b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != A*C (aij*dense)");
1199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Test D = B*C (dense*dense) */
1229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
123c4762a1bSJed Brown   if (size == 1) {
1249566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B, C, MAT_INITIAL_MATRIX, fill, &D));
1259566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B, C, MAT_REUSE_MATRIX, fill, &D));
1269566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(B, C, D, 10, &equal));
12728b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != B*C (dense*dense)");
1289566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown 
131c20d7725SJed Brown   /* Test D = B*C^T (dense*dense) */
1329566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B, C, MAT_INITIAL_MATRIX, fill, &D));
1339566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B, C, MAT_REUSE_MATRIX, fill, &D));
1349566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMultEqual(B, C, D, 10, &equal));
13528b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != B*C^T (dense*dense)");
1369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
137c20d7725SJed Brown 
1383fbca975SStefano Zampini   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
1393fbca975SStefano Zampini   flg = PETSC_FALSE;
1409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_userAPI", &flg));
1413fbca975SStefano Zampini   if (flg) {
1423fbca975SStefano Zampini     /* user driver */
1439566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, fill, &B));
1443fbca975SStefano Zampini   } else {
1453fbca975SStefano Zampini     /* clear internal data structures related with previous products to avoid circular references */
1469566063dSJacob Faibussowitsch     PetscCall(MatProductClear(A));
1479566063dSJacob Faibussowitsch     PetscCall(MatProductClear(B));
1489566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
1499566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(A, C, NULL, B));
1509566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B, MATPRODUCT_AB));
1519566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
1529566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
1539566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(B));
1543fbca975SStefano Zampini   }
1553fbca975SStefano Zampini 
1569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1599566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
160b122ec5aSJacob Faibussowitsch   return 0;
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown /*TEST
164c4762a1bSJed Brown 
165c4762a1bSJed Brown    test:
166c4762a1bSJed Brown       args: -M 10 -N 10
167*3886731fSPierre Jolivet       output_file: output/empty.out
168c4762a1bSJed Brown 
169c4762a1bSJed Brown    test:
170c4762a1bSJed Brown       suffix: 2
171c4762a1bSJed Brown       nsize: 3
172*3886731fSPierre Jolivet       output_file: output/empty.out
173c4762a1bSJed Brown 
174c20d7725SJed Brown    test:
175c20d7725SJed Brown       suffix: 3
176c20d7725SJed Brown       nsize: 2
177c20d7725SJed Brown       args: -matmattransmult_mpidense_mpidense_via cyclic
178*3886731fSPierre Jolivet       output_file: output/empty.out
179c20d7725SJed Brown 
1803fbca975SStefano Zampini    test:
1813fbca975SStefano Zampini       suffix: 4
1823fbca975SStefano Zampini       args: -test_userAPI
183*3886731fSPierre Jolivet       output_file: output/empty.out
1843fbca975SStefano Zampini 
1853fbca975SStefano Zampini    test:
1863fbca975SStefano Zampini       suffix: 5
1873fbca975SStefano Zampini       nsize: 3
1883fbca975SStefano Zampini       args: -test_userAPI
189*3886731fSPierre Jolivet       output_file: output/empty.out
1903fbca975SStefano Zampini 
191c4762a1bSJed Brown TEST*/
192