xref: /petsc/src/mat/tests/ex81.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
18c78258cSHong Zhang static char help[] = "Tests MatOption MAT_FORCE_DIAGONAL_ENTRIES.\n\n";
28c78258cSHong Zhang 
38c78258cSHong Zhang #include <petscmat.h>
48c78258cSHong Zhang 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
78c78258cSHong Zhang   Mat         A, B;
88c78258cSHong Zhang   Vec         diag;
98c78258cSHong Zhang   PetscInt    i, n = 10, col[3], test;
108c78258cSHong Zhang   PetscScalar v[3];
118c78258cSHong Zhang 
12327415f7SBarry Smith   PetscFunctionBeginUser;
13*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
158c78258cSHong Zhang 
168c78258cSHong Zhang   /* Create A which has empty 0-th row and column */
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n));
199566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
228c78258cSHong Zhang 
239371c9d4SSatish Balay   v[0] = -1.;
249371c9d4SSatish Balay   v[1] = 2.;
259371c9d4SSatish Balay   v[2] = -1.;
268c78258cSHong Zhang   for (i = 2; i < n - 1; i++) {
279371c9d4SSatish Balay     col[0] = i - 1;
289371c9d4SSatish Balay     col[1] = i;
299371c9d4SSatish Balay     col[2] = i + 1;
309566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, col, v, INSERT_VALUES));
318c78258cSHong Zhang   }
329371c9d4SSatish Balay   i      = 1;
339371c9d4SSatish Balay   col[0] = 1;
349371c9d4SSatish Balay   col[1] = 2;
359566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, v + 1, INSERT_VALUES));
369371c9d4SSatish Balay   i      = n - 1;
379371c9d4SSatish Balay   col[0] = n - 2;
389371c9d4SSatish Balay   col[1] = n - 1;
399566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, v, INSERT_VALUES));
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
428c78258cSHong Zhang 
438c78258cSHong Zhang   for (test = 0; test < 2; test++) {
449566063dSJacob Faibussowitsch     PetscCall(MatProductCreate(A, A, NULL, &B));
458c78258cSHong Zhang 
468c78258cSHong Zhang     if (test == 0) {
478c78258cSHong Zhang       /* Compute B = A*A; B misses 0-th diagonal */
489566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(B, MATPRODUCT_AB));
499566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(B, "AB_"));
508c78258cSHong Zhang     } else {
518c78258cSHong Zhang       /* Compute B = A^t*A; B misses 0-th diagonal */
529566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(B, MATPRODUCT_AtB));
539566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(B, "AtB_"));
548c78258cSHong Zhang     }
558c78258cSHong Zhang 
568c78258cSHong Zhang     /* Force allocate missing diagonal entries of B */
579566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_FORCE_DIAGONAL_ENTRIES, PETSC_TRUE));
589566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
598c78258cSHong Zhang 
609566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
619566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(B));
628c78258cSHong Zhang 
639566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
648c78258cSHong Zhang 
658c78258cSHong Zhang     /* Insert entries to diagonal of B */
669566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, NULL, &diag));
679566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(B, diag));
689566063dSJacob Faibussowitsch     PetscCall(VecSetValue(diag, 0, 100.0, INSERT_VALUES));
699566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(diag));
709566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(diag));
718c78258cSHong Zhang 
729566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(B, diag, INSERT_VALUES));
7348a46eb9SPierre Jolivet     if (test == 1) PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
749566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
768c78258cSHong Zhang   }
778c78258cSHong Zhang 
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
80b122ec5aSJacob Faibussowitsch   return 0;
818c78258cSHong Zhang }
828c78258cSHong Zhang 
838c78258cSHong Zhang /*TEST
848c78258cSHong Zhang 
858c78258cSHong Zhang    test:
868c78258cSHong Zhang      output_file: output/ex81_1.out
878c78258cSHong Zhang 
888c78258cSHong Zhang    test:
898c78258cSHong Zhang      suffix: 2
903e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
918c78258cSHong Zhang      output_file: output/ex81_1.out
928c78258cSHong Zhang 
938c78258cSHong Zhang    test:
948c78258cSHong Zhang      suffix: 3
953e662e0bSHong Zhang      args: -AtB_mat_product_algorithm outerproduct
968c78258cSHong Zhang      output_file: output/ex81_1.out
978c78258cSHong Zhang 
988c78258cSHong Zhang    test:
998c78258cSHong Zhang      suffix: 4
1008c78258cSHong Zhang      nsize: 3
1013e662e0bSHong Zhang      args: -AtB_mat_product_algorithm nonscalable
1028c78258cSHong Zhang      output_file: output/ex81_3.out
1038c78258cSHong Zhang 
1048c78258cSHong Zhang    test:
1058c78258cSHong Zhang      suffix: 5
1068c78258cSHong Zhang      nsize: 3
1073e662e0bSHong Zhang      args: -AtB_mat_product_algorithm scalable
1088c78258cSHong Zhang      output_file: output/ex81_3.out
1098c78258cSHong Zhang 
1108c78258cSHong Zhang    test:
1118c78258cSHong Zhang      suffix: 6
1128c78258cSHong Zhang      nsize: 3
1133e662e0bSHong Zhang      args: -AtB_mat_product_algorithm at*b
1148c78258cSHong Zhang      output_file: output/ex81_3.out
1158c78258cSHong Zhang 
1168c78258cSHong Zhang TEST*/
117