1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests incorrect use of MatDiagonalSet() for SHELL matrices\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown typedef struct _n_User *User; 7c4762a1bSJed Brown struct _n_User { 8c4762a1bSJed Brown Mat B; 9c4762a1bSJed Brown }; 10c4762a1bSJed Brown 11c4762a1bSJed Brown static PetscErrorCode MatGetDiagonal_User(Mat A,Vec X) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown User user; 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBegin; 16*9566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&user)); 17*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(user->B,X)); 18c4762a1bSJed Brown PetscFunctionReturn(0); 19c4762a1bSJed Brown } 20c4762a1bSJed Brown 21c4762a1bSJed Brown int main(int argc,char **args) 22c4762a1bSJed Brown { 23c4762a1bSJed Brown const PetscScalar xvals[] = {11,13}; 24c4762a1bSJed Brown const PetscInt inds[] = {0,1}; 25c4762a1bSJed Brown PetscScalar avals[] = {2,3,5,7}; 26c4762a1bSJed Brown Mat A,S; 27c4762a1bSJed Brown Vec X,Y; 28c4762a1bSJed Brown User user; 29c4762a1bSJed Brown 30*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 31*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A)); 32*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 33*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES)); 34*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 35*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 36*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD,2,&X)); 37*9566063dSJacob Faibussowitsch PetscCall(VecSetValues(X,2,inds,xvals,INSERT_VALUES)); 38*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 39*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 40*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 41c4762a1bSJed Brown 42*9566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 43c4762a1bSJed Brown user->B = A; 44c4762a1bSJed Brown 45*9566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S)); 46*9566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User)); 47*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 48c4762a1bSJed Brown 49*9566063dSJacob Faibussowitsch PetscCall(MatShift(S,42)); 50*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 51*9566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 52*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(S,X,ADD_VALUES)); 53*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 54*9566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 55*9566063dSJacob Faibussowitsch PetscCall(MatScale(S,42)); 56*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 57*9566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 58c4762a1bSJed Brown 59*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 60*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 61*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 62*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 63*9566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 64*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 65b122ec5aSJacob Faibussowitsch return 0; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /*TEST 69c4762a1bSJed Brown 70c4762a1bSJed Brown test: 71c4762a1bSJed Brown args: -malloc_dump 72c4762a1bSJed Brown 73c4762a1bSJed Brown TEST*/ 74