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 5*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6*d71ae5a4SJacob 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; 19c4762a1bSJed Brown 20327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Create MAIJ matrix, P */ 269566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &pA)); 279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pA, 3, 3, 3, 3)); 289566063dSJacob Faibussowitsch PetscCall(MatSetType(pA, MATSEQAIJ)); 299566063dSJacob Faibussowitsch PetscCall(MatSetUp(pA)); 309566063dSJacob Faibussowitsch PetscCall(MatSetOption(pA, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 319566063dSJacob Faibussowitsch PetscCall(MatSetValues(pA, 3, pij, 3, pij, pa, ADD_VALUES)); 329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pA, MAT_FINAL_ASSEMBLY)); 339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pA, MAT_FINAL_ASSEMBLY)); 349566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(pA, 3, &P)); 359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pA)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Create AIJ equivalent matrix, aijP, for comparison testing */ 389566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATSEQAIJ, MAT_INITIAL_MATRIX, &aijP)); 39c4762a1bSJed Brown 40c20d7725SJed Brown /* Create AIJ matrix A */ 419566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 9, 9, 9, 9)); 439566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ)); 449566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 459566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[0], 3, aij[0], pa, ADD_VALUES)); 479566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[1], 3, aij[1], pa, ADD_VALUES)); 489566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, aij[2], 3, aij[2], pa, ADD_VALUES)); 4948a46eb9SPierre Jolivet for (i = 0; i < 9; i++) PetscCall(MatSetValue(A, i, i, one, ADD_VALUES)); 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */ 549566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, aijP, MAT_INITIAL_MATRIX, 1., &C)); 55c20d7725SJed Brown 56c20d7725SJed Brown /* Perform PtAP_SeqAIJ_SeqMAIJ */ 57c20d7725SJed Brown /* Developer API */ 589566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, P, NULL, &mC)); 599566063dSJacob Faibussowitsch PetscCall(MatProductSetType(mC, MATPRODUCT_PtAP)); 609566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(mC, "default")); 619566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(mC, PETSC_DEFAULT)); 629566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mC)); 639566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mC)); 649566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC)); 659566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Check mC = C */ 689566063dSJacob Faibussowitsch PetscCall(MatEqual(C, mC, &flg)); 6928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatProduct C != mC"); 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC)); 71c20d7725SJed Brown 72c20d7725SJed Brown /* User API */ 739566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1., &mC)); 749566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, 1., &mC)); 75c20d7725SJed Brown 76c20d7725SJed Brown /* Check mC = C */ 779566063dSJacob Faibussowitsch PetscCall(MatEqual(C, mC, &flg)); 7828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "MatPtAP C != mC"); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Cleanup */ 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aijP)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 87b122ec5aSJacob Faibussowitsch return 0; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /*TEST 91c4762a1bSJed Brown 92c4762a1bSJed Brown test: 93c4762a1bSJed Brown output_file: output/ex101.out 94c4762a1bSJed Brown 95c4762a1bSJed Brown TEST*/ 96