1c4762a1bSJed Brown static char help[] = "Tests MatCopy() for SHELL matrices\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct _n_User *User; 6c4762a1bSJed Brown struct _n_User { 7c4762a1bSJed Brown Mat A; 8c4762a1bSJed Brown }; 9c4762a1bSJed Brown 10c4762a1bSJed Brown static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown User user; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscFunctionBegin; 155f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(user->A,X,Y)); 17c4762a1bSJed Brown PetscFunctionReturn(0); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 20c4762a1bSJed Brown static PetscErrorCode MatCopy_User(Mat A,Mat B,MatStructure str) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown User userA,userB; 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionBegin; 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&userA)); 26c4762a1bSJed Brown if (userA) { 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&userB)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(userA->A,MAT_COPY_VALUES,&userB->A)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(B, userB)); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown static PetscErrorCode MatDestroy_User(Mat A) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown User user; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBegin; 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A, &user)); 40c4762a1bSJed Brown if (user) { 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->A)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown int main(int argc,char **args) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown const PetscScalar xvals[] = {11,13},yvals[] = {17,19}; 50c4762a1bSJed Brown const PetscInt inds[] = {0,1}; 51c4762a1bSJed Brown PetscScalar avals[] = {2,3,5,7}; 52c4762a1bSJed Brown Mat S1,S2; 53c4762a1bSJed Brown Vec X,Y; 54c4762a1bSJed Brown User user; 55c4762a1bSJed Brown 56*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 57c4762a1bSJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&user->A)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user->A)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(user->A,2,inds,2,inds,avals,INSERT_VALUES)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->A,MAT_FINAL_ASSEMBLY)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->A,MAT_FINAL_ASSEMBLY)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,2,&X)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(X,2,inds,xvals,INSERT_VALUES)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(X)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(X)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&Y)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(Y,2,inds,yvals,INSERT_VALUES)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(Y)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(Y)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S1)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(S1)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_MULT,(void (*)(void))MatMult_User)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_COPY,(void (*)(void))MatCopy_User)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_DESTROY,(void (*)(void))MatDestroy_User)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,NULL,&S2)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(S2)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_MULT,(void (*)(void))MatMult_User)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_COPY,(void (*)(void))MatCopy_User)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_DESTROY,(void (*)(void))MatDestroy_User)); 83c4762a1bSJed Brown 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(S1,31)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(S1,37)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(S1,X,Y)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(S1,S2,SAME_NONZERO_PATTERN)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S1,X,Y)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S2,X,Y)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 92c4762a1bSJed Brown 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S1)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S2)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 97*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 98*b122ec5aSJacob Faibussowitsch return 0; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /*TEST 102c4762a1bSJed Brown 103c4762a1bSJed Brown test: 104c4762a1bSJed Brown args: -malloc_dump 105c4762a1bSJed Brown 106c4762a1bSJed Brown TEST*/ 107