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; 169566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A,&user)); 179566063dSJacob 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*327415f7SBarry Smith PetscFunctionBeginUser; 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 329566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A)); 339566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 349566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES)); 359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 379566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD,2,&X)); 389566063dSJacob Faibussowitsch PetscCall(VecSetValues(X,2,inds,xvals,INSERT_VALUES)); 399566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 409566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&Y)); 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 44c4762a1bSJed Brown user->B = A; 45c4762a1bSJed Brown 469566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S)); 479566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User)); 489566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(MatShift(S,42)); 519566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 529566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 539566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(S,X,ADD_VALUES)); 549566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 559566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 569566063dSJacob Faibussowitsch PetscCall(MatScale(S,42)); 579566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S,Y)); 589566063dSJacob Faibussowitsch PetscCall(VecView(Y,PETSC_VIEWER_STDOUT_WORLD)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 649566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 66b122ec5aSJacob Faibussowitsch return 0; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown /*TEST 70c4762a1bSJed Brown 71c4762a1bSJed Brown test: 72c4762a1bSJed Brown args: -malloc_dump 73c4762a1bSJed Brown 74c4762a1bSJed Brown TEST*/ 75