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