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
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat pA, P, aijP;
8c4762a1bSJed Brown PetscScalar pa[] = {1., -1., 0., 0., 1., -1., 0., 0., 1.};
9c4762a1bSJed Brown PetscInt i, pij[] = {0, 1, 2};
109371c9d4SSatish Balay PetscInt aij[3][3] = {
119371c9d4SSatish Balay {0, 1, 2},
129371c9d4SSatish Balay {3, 4, 5},
139371c9d4SSatish Balay {6, 7, 8}
149371c9d4SSatish Balay };
15c4762a1bSJed Brown Mat A, mC, C;
16c4762a1bSJed Brown PetscScalar one = 1.;
17c4762a1bSJed Brown PetscMPIInt size;
18c20d7725SJed Brown PetscBool flg;
198fcac130SJose E. Roman MatProductAlgorithm alg;
20c4762a1bSJed Brown
21327415f7SBarry Smith PetscFunctionBeginUser;
22*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
25c4762a1bSJed Brown
26c4762a1bSJed Brown /* Create MAIJ matrix, P */
279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &pA));
289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pA, 3, 3, 3, 3));
299566063dSJacob Faibussowitsch PetscCall(MatSetType(pA, MATSEQAIJ));
309566063dSJacob Faibussowitsch PetscCall(MatSetUp(pA));
319566063dSJacob Faibussowitsch PetscCall(MatSetOption(pA, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
329566063dSJacob Faibussowitsch PetscCall(MatSetValues(pA, 3, pij, 3, pij, pa, ADD_VALUES));
339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pA, MAT_FINAL_ASSEMBLY));
349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pA, MAT_FINAL_ASSEMBLY));
359566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(pA, 3, &P));
369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pA));
37c4762a1bSJed Brown
38c4762a1bSJed Brown /* Create AIJ equivalent matrix, aijP, for comparison testing */
399566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATSEQAIJ, MAT_INITIAL_MATRIX, &aijP));
40c4762a1bSJed Brown
41c20d7725SJed Brown /* Create AIJ matrix A */
429566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A));
439566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 9, 9, 9, 9));
449566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ));
459566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
469566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
479566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[0], 3, aij[0], pa, ADD_VALUES));
489566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[1], 3, aij[1], pa, ADD_VALUES));
499566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[2], 3, aij[2], pa, ADD_VALUES));
5048a46eb9SPierre Jolivet for (i = 0; i < 9; i++) PetscCall(MatSetValue(A, i, i, one, ADD_VALUES));
519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
53c4762a1bSJed Brown
54c4762a1bSJed Brown /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */
559566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, aijP, MAT_INITIAL_MATRIX, 1., &C));
56c20d7725SJed Brown
57c20d7725SJed Brown /* Perform PtAP_SeqAIJ_SeqMAIJ */
58c20d7725SJed Brown /* Developer API */
599566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, P, NULL, &mC));
609566063dSJacob Faibussowitsch PetscCall(MatProductSetType(mC, MATPRODUCT_PtAP));
618fcac130SJose E. Roman PetscCall(MatProductSetAlgorithm(mC, MATPRODUCTALGORITHMDEFAULT));
629566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(mC, PETSC_DEFAULT));
639566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mC));
649566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mC));
659566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC));
669566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC));
678fcac130SJose E. Roman PetscCall(MatProductGetAlgorithm(mC, &alg));
688fcac130SJose E. Roman PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatProduct algorithm: %s\n", alg));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* Check mC = C */
719566063dSJacob Faibussowitsch PetscCall(MatEqual(C, mC, &flg));
7228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatProduct C != mC");
739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC));
74c20d7725SJed Brown
75c20d7725SJed Brown /* User API */
769566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1., &mC));
779566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, 1., &mC));
78c20d7725SJed Brown
79c20d7725SJed Brown /* Check mC = C */
809566063dSJacob Faibussowitsch PetscCall(MatEqual(C, mC, &flg));
8128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatPtAP C != mC");
829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC));
83c4762a1bSJed Brown
84c4762a1bSJed Brown /* Cleanup */
859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P));
869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aijP));
879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
899566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
90b122ec5aSJacob Faibussowitsch return 0;
91c4762a1bSJed Brown }
92c4762a1bSJed Brown
93c4762a1bSJed Brown /*TEST
94c4762a1bSJed Brown
95c4762a1bSJed Brown test:
96c4762a1bSJed Brown output_file: output/ex101.out
97c4762a1bSJed Brown
98c4762a1bSJed Brown TEST*/
99