1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatShift(), MatScale(), and MatDiagonalScale() for SHELL and NEST 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 MatView_User(Mat A,PetscViewer viewer) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown User user; 14c4762a1bSJed Brown PetscErrorCode ierr; 15c4762a1bSJed Brown 16c4762a1bSJed Brown PetscFunctionBegin; 17c4762a1bSJed Brown ierr = MatShellGetContext(A,&user);CHKERRQ(ierr); 18c4762a1bSJed Brown ierr = MatView(user->B,viewer);CHKERRQ(ierr); 19c4762a1bSJed Brown PetscFunctionReturn(0); 20c4762a1bSJed Brown } 21c4762a1bSJed Brown 22c4762a1bSJed Brown static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown User user; 25c4762a1bSJed Brown PetscErrorCode ierr; 26c4762a1bSJed Brown 27c4762a1bSJed Brown PetscFunctionBegin; 28c4762a1bSJed Brown ierr = MatShellGetContext(A,&user);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatMult(user->B,X,Y);CHKERRQ(ierr); 30c4762a1bSJed Brown PetscFunctionReturn(0); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown static PetscErrorCode MatMultTranspose_User(Mat A,Vec X,Vec Y) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown User user; 36c4762a1bSJed Brown PetscErrorCode ierr; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBegin; 39c4762a1bSJed Brown ierr = MatShellGetContext(A,&user);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatMultTranspose(user->B,X,Y);CHKERRQ(ierr); 41c4762a1bSJed Brown PetscFunctionReturn(0); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown static PetscErrorCode MatGetDiagonal_User(Mat A,Vec X) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown User user; 47c4762a1bSJed Brown PetscErrorCode ierr; 48c4762a1bSJed Brown 49c4762a1bSJed Brown PetscFunctionBegin; 50c4762a1bSJed Brown ierr = MatShellGetContext(A,&user);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatGetDiagonal(user->B,X);CHKERRQ(ierr); 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown static PetscErrorCode TestMatrix(Mat A,Vec X,Vec Y,Vec Z) 56c4762a1bSJed Brown { 57c4762a1bSJed Brown PetscErrorCode ierr; 58c4762a1bSJed Brown Vec W1,W2,diff; 59c4762a1bSJed Brown Mat E; 60c4762a1bSJed Brown const char *mattypename; 61c4762a1bSJed Brown PetscViewer viewer = PETSC_VIEWER_STDOUT_WORLD; 62c4762a1bSJed Brown PetscScalar diag[2] = { 2.9678190300000000e+08, 1.4173141580000000e+09}; 63c4762a1bSJed Brown PetscScalar multadd[2] = {-6.8966198500000000e+08,-2.0310609940000000e+09}; 64c4762a1bSJed Brown PetscScalar multtadd[2] = {-9.1052873900000000e+08,-1.8101942400000000e+09}; 65c4762a1bSJed Brown PetscReal nrm; 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscFunctionBegin; 68c4762a1bSJed Brown ierr = PetscObjectGetType((PetscObject)A,&mattypename);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"\nMatrix of type: %s\n",mattypename);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = VecDuplicate(X,&W1);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = VecDuplicate(X,&W2);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = MatScale(A,31);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = MatShift(A,37);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = MatDiagonalScale(A,X,Y);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = MatScale(A,41);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatDiagonalScale(A,Y,Z);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = MatComputeOperator(A,MATDENSE,&E);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = MatView(E,viewer);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Testing MatMult + MatMultTranspose\n");CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = MatMult(A,Z,W1);CHKERRQ(ierr); 82c4762a1bSJed Brown ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Testing MatMultAdd\n");CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multadd,&diff);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = VecSet(W1,-1.0);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatMultAdd(A,W1,W1,W2);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = VecAXPY(W2,-1.0,diff);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = VecNorm(W2,NORM_2,&nrm);CHKERRQ(ierr); 92c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 93c4762a1bSJed Brown if (nrm > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultAdd(A,x,x,y) produces incorrect result"); 94c4762a1bSJed Brown #endif 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = VecSet(W2,-1.0);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = MatMultAdd(A,W1,W2,W2);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = VecAXPY(W2,-1.0,diff);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = VecNorm(W2,NORM_2,&nrm);CHKERRQ(ierr); 101c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 102c4762a1bSJed Brown if (nrm > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultAdd(A,x,y,y) produces incorrect result"); 103c4762a1bSJed Brown #endif 104c4762a1bSJed Brown ierr = VecDestroy(&diff);CHKERRQ(ierr); 105c4762a1bSJed Brown 106*2b1bc0daSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer,"Testing MatMultTransposeAdd\n");CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multtadd,&diff);CHKERRQ(ierr); 108c4762a1bSJed Brown 109c4762a1bSJed Brown ierr = VecSet(W1,-1.0);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = MatMultTransposeAdd(A,W1,W1,W2);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = VecAXPY(W2,-1.0,diff);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecNorm(W2,NORM_2,&nrm);CHKERRQ(ierr); 114c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 115*2b1bc0daSPierre Jolivet if (nrm > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultTransposeAdd(A,x,x,y) produces incorrect result"); 116c4762a1bSJed Brown #endif 117c4762a1bSJed Brown 118c4762a1bSJed Brown ierr = VecSet(W2,-1.0);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatMultTransposeAdd(A,W1,W2,W2);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecAXPY(W2,-1.0,diff);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecNorm(W2,NORM_2,&nrm);CHKERRQ(ierr); 123c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 124*2b1bc0daSPierre Jolivet if (nrm > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultTransposeAdd(A,x,y,y) produces incorrect result"); 125c4762a1bSJed Brown #endif 126c4762a1bSJed Brown ierr = VecDestroy(&diff);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Testing MatGetDiagonal\n");CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = MatGetDiagonal(A,W2);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = VecView(W2,viewer);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,diag,&diff);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = VecAXPY(diff,-1.0,W2);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = VecNorm(diff,NORM_2,&nrm);CHKERRQ(ierr); 134c4762a1bSJed Brown if (nrm > PETSC_SMALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() produces incorrect result"); 135c4762a1bSJed Brown ierr = VecDestroy(&diff);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* MATSHELL does not support MatDiagonalSet after MatScale */ 138c4762a1bSJed Brown if (strncmp(mattypename, "shell", 5)) { 139c4762a1bSJed Brown ierr = MatDiagonalSet(A,X,INSERT_VALUES);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = MatGetDiagonal(A,W1);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = VecView(W1,viewer);CHKERRQ(ierr); 142c4762a1bSJed Brown } else { 143c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"MatDiagonalSet not tested on MATSHELL\n");CHKERRQ(ierr); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = MatDestroy(&E);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = VecDestroy(&W1);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = VecDestroy(&W2);CHKERRQ(ierr); 149c4762a1bSJed Brown PetscFunctionReturn(0); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown int main(int argc,char **args) 153c4762a1bSJed Brown { 154c4762a1bSJed Brown const PetscScalar xvals[] = {11,13},yvals[] = {17,19},zvals[] = {23,29}; 155c4762a1bSJed Brown const PetscInt inds[] = {0,1}; 156c4762a1bSJed Brown PetscScalar avals[] = {2,3,5,7}; 157c4762a1bSJed Brown Mat A,S,D[4],N; 158c4762a1bSJed Brown Vec X,Y,Z; 159c4762a1bSJed Brown User user; 160c4762a1bSJed Brown PetscInt i; 161c4762a1bSJed Brown PetscErrorCode ierr; 162c4762a1bSJed Brown 163c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 164c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_WORLD,2,&X);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecDuplicate(X,&Z);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = VecSetValues(X,2,inds,xvals,INSERT_VALUES);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = VecSetValues(Y,2,inds,yvals,INSERT_VALUES);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = VecSetValues(Z,2,inds,zvals,INSERT_VALUES);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = VecAssemblyBegin(X);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = VecAssemblyBegin(Z);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = VecAssemblyEnd(X);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = VecAssemblyEnd(Z);CHKERRQ(ierr); 181c4762a1bSJed Brown 182c4762a1bSJed Brown ierr = PetscNew(&user);CHKERRQ(ierr); 183c4762a1bSJed Brown user->B = A; 184c4762a1bSJed Brown 185c4762a1bSJed Brown ierr = MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = MatSetUp(S);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = MatShellSetOperation(S,MATOP_VIEW,(void (*)(void))MatView_User);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User);CHKERRQ(ierr); 191c4762a1bSJed Brown 192c4762a1bSJed Brown for (i=0; i<4; i++) { 193c4762a1bSJed Brown ierr = MatCreateSeqDense(PETSC_COMM_WORLD,1,1,&avals[i],&D[i]);CHKERRQ(ierr); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,D,&N);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = MatSetUp(N);CHKERRQ(ierr); 197c4762a1bSJed Brown 198c4762a1bSJed Brown ierr = TestMatrix(S,X,Y,Z);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = TestMatrix(A,X,Y,Z);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = TestMatrix(N,X,Y,Z);CHKERRQ(ierr); 201c4762a1bSJed Brown 202c4762a1bSJed Brown for (i=0; i<4; i++) {ierr = MatDestroy(&D[i]);CHKERRQ(ierr);} 203c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatDestroy(&S);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = MatDestroy(&N);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = VecDestroy(&Z);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = PetscFree(user);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = PetscFinalize(); 211c4762a1bSJed Brown return ierr; 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown /*TEST 215c4762a1bSJed Brown 216c4762a1bSJed Brown test: 217c4762a1bSJed Brown 218c4762a1bSJed Brown TEST*/ 219