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 119371c9d4SSatish Balay static PetscErrorCode MatView_User(Mat A, PetscViewer viewer) { 12c4762a1bSJed Brown User user; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user)); 169566063dSJacob Faibussowitsch PetscCall(MatView(user->B, viewer)); 17c4762a1bSJed Brown PetscFunctionReturn(0); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 209371c9d4SSatish Balay static PetscErrorCode MatMult_User(Mat A, Vec X, Vec Y) { 21c4762a1bSJed Brown User user; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user)); 259566063dSJacob Faibussowitsch PetscCall(MatMult(user->B, X, Y)); 26c4762a1bSJed Brown PetscFunctionReturn(0); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown 299371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_User(Mat A, Vec X, Vec Y) { 30c4762a1bSJed Brown User user; 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user)); 349566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->B, X, Y)); 35c4762a1bSJed Brown PetscFunctionReturn(0); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 389371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_User(Mat A, Vec X) { 39c4762a1bSJed Brown User user; 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscFunctionBegin; 429566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user)); 439566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(user->B, X)); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 479371c9d4SSatish Balay static PetscErrorCode TestMatrix(Mat A, Vec X, Vec Y, Vec Z) { 48c4762a1bSJed Brown Vec W1, W2, diff; 49c4762a1bSJed Brown Mat E; 50c4762a1bSJed Brown const char *mattypename; 51c4762a1bSJed Brown PetscViewer viewer = PETSC_VIEWER_STDOUT_WORLD; 52c4762a1bSJed Brown PetscScalar diag[2] = {2.9678190300000000e+08, 1.4173141580000000e+09}; 53c4762a1bSJed Brown PetscScalar multadd[2] = {-6.8966198500000000e+08, -2.0310609940000000e+09}; 54c4762a1bSJed Brown PetscScalar multtadd[2] = {-9.1052873900000000e+08, -1.8101942400000000e+09}; 55c4762a1bSJed Brown PetscReal nrm; 56c4762a1bSJed Brown 57c4762a1bSJed Brown PetscFunctionBegin; 589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject)A, &mattypename)); 599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nMatrix of type: %s\n", mattypename)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &W1)); 619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &W2)); 629566063dSJacob Faibussowitsch PetscCall(MatScale(A, 31)); 639566063dSJacob Faibussowitsch PetscCall(MatShift(A, 37)); 649566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, X, Y)); 659566063dSJacob Faibussowitsch PetscCall(MatScale(A, 41)); 669566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, Y, Z)); 679566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(A, MATDENSE, &E)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatView(E, viewer)); 709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMult + MatMultTranspose\n")); 719566063dSJacob Faibussowitsch PetscCall(MatMult(A, Z, W1)); 729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2)); 739566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultAdd\n")); 769566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, multadd, &diff)); 779566063dSJacob Faibussowitsch PetscCall(VecSet(W1, -1.0)); 789566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, W1, W1, W2)); 799566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 809566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff)); 819566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm)); 82c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 8308401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultAdd(A,x,x,y) produces incorrect result"); 84c4762a1bSJed Brown #endif 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(VecSet(W2, -1.0)); 879566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, W1, W2, W2)); 889566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 899566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff)); 909566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm)); 91c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 9208401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultAdd(A,x,y,y) produces incorrect result"); 93c4762a1bSJed Brown #endif 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff)); 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultTransposeAdd\n")); 979566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, multtadd, &diff)); 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(VecSet(W1, -1.0)); 1009566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, W1, W1, W2)); 1019566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 1029566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff)); 1039566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm)); 104c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 10508401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultTransposeAdd(A,x,x,y) produces incorrect result"); 106c4762a1bSJed Brown #endif 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(VecSet(W2, -1.0)); 1099566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, W1, W2, W2)); 1109566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 1119566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff)); 1129566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm)); 113c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128) 11408401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultTransposeAdd(A,x,y,y) produces incorrect result"); 115c4762a1bSJed Brown #endif 1169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff)); 117c4762a1bSJed Brown 1189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatGetDiagonal\n")); 1199566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, W2)); 1209566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer)); 1219566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, diag, &diff)); 1229566063dSJacob Faibussowitsch PetscCall(VecAXPY(diff, -1.0, W2)); 1239566063dSJacob Faibussowitsch PetscCall(VecNorm(diff, NORM_2, &nrm)); 12408401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetDiagonal() produces incorrect result"); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* MATSHELL does not support MatDiagonalSet after MatScale */ 128c4762a1bSJed Brown if (strncmp(mattypename, "shell", 5)) { 1299566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, X, INSERT_VALUES)); 1309566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, W1)); 1319566063dSJacob Faibussowitsch PetscCall(VecView(W1, viewer)); 132c4762a1bSJed Brown } else { 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MatDiagonalSet not tested on MATSHELL\n")); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E)); 1379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W1)); 1389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 1429371c9d4SSatish Balay int main(int argc, char **args) { 143c4762a1bSJed Brown const PetscScalar xvals[] = {11, 13}, yvals[] = {17, 19}, zvals[] = {23, 29}; 144c4762a1bSJed Brown const PetscInt inds[] = {0, 1}; 145c4762a1bSJed Brown PetscScalar avals[] = {2, 3, 5, 7}; 146c4762a1bSJed Brown Mat A, S, D[4], N; 147c4762a1bSJed Brown Vec X, Y, Z; 148c4762a1bSJed Brown User user; 149c4762a1bSJed Brown PetscInt i; 150c4762a1bSJed Brown 151327415f7SBarry Smith PetscFunctionBeginUser; 1529566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 1539566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2, 2, 2, NULL, &A)); 1549566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, inds, 2, inds, avals, INSERT_VALUES)); 1569566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 2, &X)); 1579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Z)); 1599566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 2, inds, xvals, INSERT_VALUES)); 1609566063dSJacob Faibussowitsch PetscCall(VecSetValues(Y, 2, inds, yvals, INSERT_VALUES)); 1619566063dSJacob Faibussowitsch PetscCall(VecSetValues(Z, 2, inds, zvals, INSERT_VALUES)); 1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1649566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 1659566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y)); 1669566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Z)); 1679566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 1689566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y)); 1699566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Z)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 172c4762a1bSJed Brown user->B = A; 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, 2, 2, 2, 2, user, &S)); 1759566063dSJacob Faibussowitsch PetscCall(MatSetUp(S)); 1769566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_VIEW, (void (*)(void))MatView_User)); 1779566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MatMult_User)); 1789566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_User)); 1799566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_User)); 180c4762a1bSJed Brown 181*48a46eb9SPierre Jolivet for (i = 0; i < 4; i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, 1, 1, &avals[i], &D[i])); 1829566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, D, &N)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetUp(N)); 184c4762a1bSJed Brown 1859566063dSJacob Faibussowitsch PetscCall(TestMatrix(S, X, Y, Z)); 1869566063dSJacob Faibussowitsch PetscCall(TestMatrix(A, X, Y, Z)); 1879566063dSJacob Faibussowitsch PetscCall(TestMatrix(N, X, Y, Z)); 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch for (i = 0; i < 4; i++) PetscCall(MatDestroy(&D[i])); 1909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&N)); 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Z)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFree(user)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 198b122ec5aSJacob Faibussowitsch return 0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /*TEST 202c4762a1bSJed Brown 203c4762a1bSJed Brown test: 204c4762a1bSJed Brown 205c4762a1bSJed Brown TEST*/ 206