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