1c4762a1bSJed Brown static char help[] = "Tests incorrect use of MatDiagonalSet() 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 B; 8c4762a1bSJed Brown }; 9c4762a1bSJed Brown 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_User(Mat A, Vec X) 11d71ae5a4SJacob Faibussowitsch { 12c4762a1bSJed Brown User user; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user)); 169566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(user->B, X)); 173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 20d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 21d71ae5a4SJacob Faibussowitsch { 22c4762a1bSJed Brown const PetscScalar xvals[] = {11, 13}; 23c4762a1bSJed Brown const PetscInt inds[] = {0, 1}; 24c4762a1bSJed Brown PetscScalar avals[] = {2, 3, 5, 7}; 25c4762a1bSJed Brown Mat A, S; 26c4762a1bSJed Brown Vec X, Y; 27c4762a1bSJed Brown User user; 28c4762a1bSJed Brown 29327415f7SBarry Smith PetscFunctionBeginUser; 30c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 319566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2, 2, 2, NULL, &A)); 329566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 339566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, inds, 2, inds, avals, INSERT_VALUES)); 349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 369566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 2, &X)); 379566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 2, inds, xvals, INSERT_VALUES)); 389566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 399566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 43c4762a1bSJed Brown user->B = A; 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, 2, 2, 2, 2, user, &S)); 46*57d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_User)); 479566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(MatShift(S, 42)); 509566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S, Y)); 519566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD)); 529566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(S, X, ADD_VALUES)); 539566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S, Y)); 549566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD)); 559566063dSJacob Faibussowitsch PetscCall(MatScale(S, 42)); 569566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(S, Y)); 579566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD)); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 649566063dSJacob 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