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 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 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}; 10c4762a1bSJed Brown PetscInt aij[3][3]={{0,1,2},{3,4,5},{6,7,8}}; 11c4762a1bSJed Brown Mat A,mC,C; 12c4762a1bSJed Brown PetscScalar one=1.; 13c4762a1bSJed Brown PetscMPIInt size; 14c20d7725SJed Brown PetscBool flg; 15c4762a1bSJed Brown 16*327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* Create MAIJ matrix, P */ 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&pA)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(pA,3,3,3,3)); 249566063dSJacob Faibussowitsch PetscCall(MatSetType(pA,MATSEQAIJ)); 259566063dSJacob Faibussowitsch PetscCall(MatSetUp(pA)); 269566063dSJacob Faibussowitsch PetscCall(MatSetOption(pA,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 279566063dSJacob Faibussowitsch PetscCall(MatSetValues(pA,3,pij,3,pij,pa,ADD_VALUES)); 289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(pA,MAT_FINAL_ASSEMBLY)); 299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(pA,MAT_FINAL_ASSEMBLY)); 309566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(pA,3,&P)); 319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pA)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create AIJ equivalent matrix, aijP, for comparison testing */ 349566063dSJacob Faibussowitsch PetscCall(MatConvert(P,MATSEQAIJ,MAT_INITIAL_MATRIX,&aijP)); 35c4762a1bSJed Brown 36c20d7725SJed Brown /* Create AIJ matrix A */ 379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,9,9,9,9)); 399566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 409566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 419566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 429566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,aij[0],3,aij[0],pa,ADD_VALUES)); 439566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,aij[1],3,aij[1],pa,ADD_VALUES)); 449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,aij[2],3,aij[2],pa,ADD_VALUES)); 45c4762a1bSJed Brown for (i=0; i<9; i++) { 469566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,i,one,ADD_VALUES)); 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */ 529566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,aijP,MAT_INITIAL_MATRIX,1.,&C)); 53c20d7725SJed Brown 54c20d7725SJed Brown /* Perform PtAP_SeqAIJ_SeqMAIJ */ 55c20d7725SJed Brown /* Developer API */ 569566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A,P,NULL,&mC)); 579566063dSJacob Faibussowitsch PetscCall(MatProductSetType(mC,MATPRODUCT_PtAP)); 589566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(mC,"default")); 599566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(mC,PETSC_DEFAULT)); 609566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(mC)); 619566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(mC)); 629566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC)); 639566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(mC)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Check mC = C */ 669566063dSJacob Faibussowitsch PetscCall(MatEqual(C,mC,&flg)); 6728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"MatProduct C != mC"); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC)); 69c20d7725SJed Brown 70c20d7725SJed Brown /* User API */ 719566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,1.,&mC)); 729566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,P,MAT_REUSE_MATRIX,1.,&mC)); 73c20d7725SJed Brown 74c20d7725SJed Brown /* Check mC = C */ 759566063dSJacob Faibussowitsch PetscCall(MatEqual(C,mC,&flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"MatPtAP C != mC"); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mC)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Cleanup */ 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aijP)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 85b122ec5aSJacob Faibussowitsch return 0; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /*TEST 89c4762a1bSJed Brown 90c4762a1bSJed Brown test: 91c4762a1bSJed Brown output_file: output/ex101.out 92c4762a1bSJed Brown 93c4762a1bSJed Brown TEST*/ 94