xref: /petsc/src/mat/tests/ex101.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Testing PtAP for SeqMAIJ matrix, P, with SeqAIJ matrix, A.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
59371c9d4SSatish Balay int main(int argc, char **argv) {
6c4762a1bSJed Brown   Mat         pA, P, aijP;
7c4762a1bSJed Brown   PetscScalar pa[] = {1., -1., 0., 0., 1., -1., 0., 0., 1.};
8c4762a1bSJed Brown   PetscInt    i, pij[] = {0, 1, 2};
99371c9d4SSatish Balay   PetscInt    aij[3][3] = {
109371c9d4SSatish Balay        {0, 1, 2},
119371c9d4SSatish Balay        {3, 4, 5},
129371c9d4SSatish Balay        {6, 7, 8}
139371c9d4SSatish Balay   };
14c4762a1bSJed Brown   Mat         A, mC, C;
15c4762a1bSJed Brown   PetscScalar one = 1.;
16c4762a1bSJed Brown   PetscMPIInt size;
17c20d7725SJed Brown   PetscBool   flg;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create MAIJ matrix, P */
259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &pA));
269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(pA, 3, 3, 3, 3));
279566063dSJacob Faibussowitsch   PetscCall(MatSetType(pA, MATSEQAIJ));
289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(pA));
299566063dSJacob Faibussowitsch   PetscCall(MatSetOption(pA, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
309566063dSJacob Faibussowitsch   PetscCall(MatSetValues(pA, 3, pij, 3, pij, pa, ADD_VALUES));
319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(pA, MAT_FINAL_ASSEMBLY));
329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(pA, MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(pA, 3, &P));
349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pA));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* Create AIJ equivalent matrix, aijP, for comparison testing */
379566063dSJacob Faibussowitsch   PetscCall(MatConvert(P, MATSEQAIJ, MAT_INITIAL_MATRIX, &aijP));
38c4762a1bSJed Brown 
39c20d7725SJed Brown   /* Create AIJ matrix A */
409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, 9, 9, 9, 9));
429566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSEQAIJ));
439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
449566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
459566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, aij[0], 3, aij[0], pa, ADD_VALUES));
469566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, aij[1], 3, aij[1], pa, ADD_VALUES));
479566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, aij[2], 3, aij[2], pa, ADD_VALUES));
48*48a46eb9SPierre Jolivet   for (i = 0; i < 9; i++) PetscCall(MatSetValue(A, i, i, one, ADD_VALUES));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */
539566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, aijP, MAT_INITIAL_MATRIX, 1., &C));
54c20d7725SJed Brown 
55c20d7725SJed Brown   /* Perform PtAP_SeqAIJ_SeqMAIJ */
56c20d7725SJed Brown   /* Developer API */
579566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(A, P, NULL, &mC));
589566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(mC, MATPRODUCT_PtAP));
599566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(mC, "default"));
609566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(mC, PETSC_DEFAULT));
619566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(mC));
629566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(mC));
639566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(mC));
649566063dSJacob Faibussowitsch   PetscCall(MatProductNumeric(mC));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Check mC = C */
679566063dSJacob Faibussowitsch   PetscCall(MatEqual(C, mC, &flg));
6828b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatProduct C != mC");
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mC));
70c20d7725SJed Brown 
71c20d7725SJed Brown   /* User API */
729566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1., &mC));
739566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, 1., &mC));
74c20d7725SJed Brown 
75c20d7725SJed Brown   /* Check mC = C */
769566063dSJacob Faibussowitsch   PetscCall(MatEqual(C, mC, &flg));
7728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatPtAP C != mC");
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mC));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Cleanup */
819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aijP));
839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
86b122ec5aSJacob Faibussowitsch   return 0;
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown /*TEST
90c4762a1bSJed Brown 
91c4762a1bSJed Brown    test:
92c4762a1bSJed Brown       output_file: output/ex101.out
93c4762a1bSJed Brown 
94c4762a1bSJed Brown TEST*/
95