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; 15*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&user)); 16*5f80ce2aSJacob 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; 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&userA)); 26c4762a1bSJed Brown if (userA) { 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&userB)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(userA->A,MAT_COPY_VALUES,&userB->A)); 29*5f80ce2aSJacob 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; 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A, &user)); 40c4762a1bSJed Brown if (user) { 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user->A)); 42*5f80ce2aSJacob 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 PetscErrorCode ierr; 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 58c4762a1bSJed Brown 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&user)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&user->A)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user->A)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(user->A,2,inds,2,inds,avals,INSERT_VALUES)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(user->A,MAT_FINAL_ASSEMBLY)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(user->A,MAT_FINAL_ASSEMBLY)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,2,&X)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(X,2,inds,xvals,INSERT_VALUES)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(X)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(X)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&Y)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(Y,2,inds,yvals,INSERT_VALUES)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(Y)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(Y)); 73c4762a1bSJed Brown 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S1)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(S1)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_MULT,(void (*)(void))MatMult_User)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_COPY,(void (*)(void))MatCopy_User)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S1,MATOP_DESTROY,(void (*)(void))MatDestroy_User)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,NULL,&S2)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(S2)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_MULT,(void (*)(void))MatMult_User)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_COPY,(void (*)(void))MatCopy_User)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S2,MATOP_DESTROY,(void (*)(void))MatDestroy_User)); 84c4762a1bSJed Brown 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(S1,31)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(S1,37)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(S1,X,Y)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(S1,S2,SAME_NONZERO_PATTERN)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S1,X,Y)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(S2,X,Y)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 93c4762a1bSJed Brown 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S1)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S2)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 98c4762a1bSJed Brown ierr = PetscFinalize(); 99c4762a1bSJed Brown return ierr; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /*TEST 103c4762a1bSJed Brown 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown args: -malloc_dump 106c4762a1bSJed Brown 107c4762a1bSJed Brown TEST*/ 108