1*c4762a1bSJed Brown static char help[] = "Testing PtAP for SeqMAIJ matrix, P, with SeqAIJ matrix, A.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc,char **argv) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown Mat pA,P,aijP; 8*c4762a1bSJed Brown PetscScalar pa[]={1.,-1.,0.,0.,1.,-1.,0.,0.,1.}; 9*c4762a1bSJed Brown PetscInt i,pij[]={0,1,2}; 10*c4762a1bSJed Brown PetscInt aij[3][3]={{0,1,2},{3,4,5},{6,7,8}}; 11*c4762a1bSJed Brown Mat A,mC,C; 12*c4762a1bSJed Brown PetscScalar one=1.; 13*c4762a1bSJed Brown PetscErrorCode ierr; 14*c4762a1bSJed Brown PetscMPIInt size; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 17*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 18*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown /* Create MAIJ matrix, P */ 21*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&pA);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = MatSetSizes(pA,3,3,3,3);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MatSetType(pA,MATSEQAIJ);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MatSetUp(pA);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetOption(pA,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatSetValues(pA,3,pij,3,pij,pa,ADD_VALUES);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatAssemblyBegin(pA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatAssemblyEnd(pA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatCreateMAIJ(pA,3,&P);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatDestroy(&pA);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown /* Create AIJ equivalent matrix, aijP, for comparison testing */ 33*c4762a1bSJed Brown ierr = MatConvert(P,MATSEQAIJ,MAT_INITIAL_MATRIX,&aijP);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* Create AIJ matrix, A */ 36*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatSetSizes(A,9,9,9,9);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatSetValues(A,3,aij[0],3,aij[0],pa,ADD_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatSetValues(A,3,aij[1],3,aij[1],pa,ADD_VALUES);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = MatSetValues(A,3,aij[2],3,aij[2],pa,ADD_VALUES);CHKERRQ(ierr); 44*c4762a1bSJed Brown for (i=0; i<9; i++) { 45*c4762a1bSJed Brown ierr = MatSetValue(A,i,i,one,ADD_VALUES);CHKERRQ(ierr); 46*c4762a1bSJed Brown } 47*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown /* Perform PtAP_SeqAIJ_SeqMAIJ */ 51*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,1.,&mC);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,1.,&mC);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatView(mC,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown /* Perform PtAP_SeqAIJ_SeqAIJ for comparison testing */ 56*c4762a1bSJed Brown ierr = MatPtAP(A,aijP,MAT_INITIAL_MATRIX,1.,&C);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown /* Check mC = C */ 60*c4762a1bSJed Brown ierr = MatAXPY(C,-1.0,mC,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 61*c4762a1bSJed Brown /* Note: We should be able to use SAME_NONZERO_PATTERN on the line above, */ 62*c4762a1bSJed Brown /* but don't because this flag doesn't assist testing. */ 63*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown /* Cleanup */ 66*c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatDestroy(&aijP);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatDestroy(&mC);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscFinalize(); 72*c4762a1bSJed Brown return ierr; 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /*TEST 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown test: 79*c4762a1bSJed Brown output_file: output/ex101.out 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown TEST*/ 82