1c4762a1bSJed Brown static char help[] = "Tests MatShift(), MatScale(), and MatDiagonalScale() for SHELL and NEST 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
MatView_User(Mat A,PetscViewer viewer)10d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_User(Mat A, PetscViewer viewer)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown User user;
13c4762a1bSJed Brown
14c4762a1bSJed Brown PetscFunctionBegin;
159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user));
169566063dSJacob Faibussowitsch PetscCall(MatView(user->B, viewer));
173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18c4762a1bSJed Brown }
19c4762a1bSJed Brown
MatMult_User(Mat A,Vec X,Vec Y)20d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_User(Mat A, Vec X, Vec Y)
21d71ae5a4SJacob Faibussowitsch {
22c4762a1bSJed Brown User user;
23c4762a1bSJed Brown
24c4762a1bSJed Brown PetscFunctionBegin;
259566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user));
269566063dSJacob Faibussowitsch PetscCall(MatMult(user->B, X, Y));
273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown
MatMultTranspose_User(Mat A,Vec X,Vec Y)30d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_User(Mat A, Vec X, Vec Y)
31d71ae5a4SJacob Faibussowitsch {
32c4762a1bSJed Brown User user;
33c4762a1bSJed Brown
34c4762a1bSJed Brown PetscFunctionBegin;
359566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user));
369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(user->B, X, Y));
373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
38c4762a1bSJed Brown }
39c4762a1bSJed Brown
MatGetDiagonal_User(Mat A,Vec X)40d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_User(Mat A, Vec X)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown User user;
43c4762a1bSJed Brown
44c4762a1bSJed Brown PetscFunctionBegin;
459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &user));
469566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(user->B, X));
473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown
TestMatrix(Mat A,Vec X,Vec Y,Vec Z)50d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestMatrix(Mat A, Vec X, Vec Y, Vec Z)
51d71ae5a4SJacob Faibussowitsch {
52cf242a8cSBlanca Mellado Pinto Vec W1, W2, W3, diff;
53c4762a1bSJed Brown Mat E;
54c4762a1bSJed Brown const char *mattypename;
55c4762a1bSJed Brown PetscViewer viewer = PETSC_VIEWER_STDOUT_WORLD;
56c4762a1bSJed Brown PetscReal nrm;
57cf242a8cSBlanca Mellado Pinto #if defined(PETSC_USE_COMPLEX)
58cf242a8cSBlanca Mellado Pinto const PetscScalar diag[2] = {PetscCMPLX(-6.2902938000000000e+07, 4.5741953400000000e+08), PetscCMPLX(1.0828994620000000e+09, 1.2955916360000000e+09)};
59cf242a8cSBlanca Mellado Pinto const PetscScalar multadd[2] = {PetscCMPLX(1.4926230300000000e+08, -1.2811063360000000e+09), PetscCMPLX(-1.2985220710000000e+09, -2.1029893020000000e+09)};
60cf242a8cSBlanca Mellado Pinto const PetscScalar multtadd[2] = {PetscCMPLX(-1.5271967100000000e+08, -1.2648172000000000e+09), PetscCMPLX(-9.9654009700000000e+08, -2.1192784380000000e+09)};
61cf242a8cSBlanca Mellado Pinto #else
62cf242a8cSBlanca Mellado Pinto const PetscScalar diag[2] = {2.9678190300000000e+08, 1.4173141580000000e+09};
63cf242a8cSBlanca Mellado Pinto const PetscScalar multadd[2] = {-6.8966198500000000e+08, -2.0310609940000000e+09};
64cf242a8cSBlanca Mellado Pinto const PetscScalar multtadd[2] = {-9.1052873900000000e+08, -1.8101942400000000e+09};
65cf242a8cSBlanca Mellado Pinto #endif
66c4762a1bSJed Brown
67c4762a1bSJed Brown PetscFunctionBegin;
689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject)A, &mattypename));
699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nMatrix of type: %s\n", mattypename));
709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &W1));
719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &W2));
72cf242a8cSBlanca Mellado Pinto PetscCall(VecDuplicate(X, &W3));
739566063dSJacob Faibussowitsch PetscCall(MatScale(A, 31));
749566063dSJacob Faibussowitsch PetscCall(MatShift(A, 37));
759566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, X, Y));
769566063dSJacob Faibussowitsch PetscCall(MatScale(A, 41));
779566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, Y, Z));
789566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(A, MATDENSE, &E));
79c4762a1bSJed Brown
809566063dSJacob Faibussowitsch PetscCall(MatView(E, viewer));
819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMult + MatMultTranspose\n"));
829566063dSJacob Faibussowitsch PetscCall(MatMult(A, Z, W1));
839566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, W1, W2));
849566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
85cf242a8cSBlanca Mellado Pinto PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultHermitianTranspose\n"));
86cf242a8cSBlanca Mellado Pinto PetscCall(VecConjugate(W1));
87cf242a8cSBlanca Mellado Pinto PetscCall(MatMultHermitianTranspose(A, W1, W2));
88cf242a8cSBlanca Mellado Pinto PetscCall(VecConjugate(W2));
89cf242a8cSBlanca Mellado Pinto PetscCall(VecView(W2, viewer));
90c4762a1bSJed Brown
919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultAdd\n"));
929566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, multadd, &diff));
939566063dSJacob Faibussowitsch PetscCall(VecSet(W1, -1.0));
949566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, W1, W1, W2));
959566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
969566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff));
979566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm));
98c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
9908401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultAdd(A,x,x,y) produces incorrect result");
100c4762a1bSJed Brown #endif
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(VecSet(W2, -1.0));
1039566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, W1, W2, W2));
1049566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
1059566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff));
1069566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm));
107c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
10808401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultAdd(A,x,y,y) produces incorrect result");
109c4762a1bSJed Brown #endif
1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff));
111c4762a1bSJed Brown
1129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultTransposeAdd\n"));
1139566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, multtadd, &diff));
114c4762a1bSJed Brown
1159566063dSJacob Faibussowitsch PetscCall(VecSet(W1, -1.0));
1169566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, W1, W1, W2));
1179566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
1189566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff));
1199566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm));
120c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
12108401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultTransposeAdd(A,x,x,y) produces incorrect result");
122c4762a1bSJed Brown #endif
123c4762a1bSJed Brown
1249566063dSJacob Faibussowitsch PetscCall(VecSet(W2, -1.0));
1259566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, W1, W2, W2));
1269566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
1279566063dSJacob Faibussowitsch PetscCall(VecAXPY(W2, -1.0, diff));
1289566063dSJacob Faibussowitsch PetscCall(VecNorm(W2, NORM_2, &nrm));
129c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
13008401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultTransposeAdd(A,x,y,y) produces incorrect result");
131c4762a1bSJed Brown #endif
1329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff));
133c4762a1bSJed Brown
134cf242a8cSBlanca Mellado Pinto PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatMultHermitianTransposeAdd\n"));
135cf242a8cSBlanca Mellado Pinto PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, multtadd, &diff));
136cf242a8cSBlanca Mellado Pinto
137cf242a8cSBlanca Mellado Pinto PetscCall(VecSet(W1, -1.0));
138cf242a8cSBlanca Mellado Pinto PetscCall(MatMultHermitianTransposeAdd(A, W1, W1, W3));
139cf242a8cSBlanca Mellado Pinto PetscCall(VecConjugate(W3));
140cf242a8cSBlanca Mellado Pinto PetscCall(VecView(W3, viewer));
141cf242a8cSBlanca Mellado Pinto PetscCall(VecAXPY(W3, -1.0, diff));
142cf242a8cSBlanca Mellado Pinto PetscCall(VecNorm(W3, NORM_2, &nrm));
143cf242a8cSBlanca Mellado Pinto #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
144cf242a8cSBlanca Mellado Pinto PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultHermitianTransposeAdd(A,x,x,y) produces incorrect result");
145cf242a8cSBlanca Mellado Pinto #endif
146cf242a8cSBlanca Mellado Pinto
147cf242a8cSBlanca Mellado Pinto PetscCall(VecSet(W3, -1.0));
148cf242a8cSBlanca Mellado Pinto PetscCall(MatMultHermitianTransposeAdd(A, W1, W3, W3));
149cf242a8cSBlanca Mellado Pinto PetscCall(VecConjugate(W3));
150cf242a8cSBlanca Mellado Pinto PetscCall(VecView(W3, viewer));
151cf242a8cSBlanca Mellado Pinto PetscCall(VecAXPY(W3, -1.0, diff));
152cf242a8cSBlanca Mellado Pinto PetscCall(VecNorm(W3, NORM_2, &nrm));
153cf242a8cSBlanca Mellado Pinto #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
154cf242a8cSBlanca Mellado Pinto PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatMultHermitianTransposeAdd(A,x,y,y) produces incorrect result");
155cf242a8cSBlanca Mellado Pinto #endif
156cf242a8cSBlanca Mellado Pinto PetscCall(VecDestroy(&diff));
157cf242a8cSBlanca Mellado Pinto
1589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Testing MatGetDiagonal\n"));
1599566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, W2));
1609566063dSJacob Faibussowitsch PetscCall(VecView(W2, viewer));
1619566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, 2, diag, &diff));
1629566063dSJacob Faibussowitsch PetscCall(VecAXPY(diff, -1.0, W2));
1639566063dSJacob Faibussowitsch PetscCall(VecNorm(diff, NORM_2, &nrm));
164cf242a8cSBlanca Mellado Pinto #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
16508401ef6SPierre Jolivet PetscCheck(nrm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetDiagonal() produces incorrect result");
166cf242a8cSBlanca Mellado Pinto #endif
1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diff));
168c4762a1bSJed Brown
169c4762a1bSJed Brown /* MATSHELL does not support MatDiagonalSet after MatScale */
170c4762a1bSJed Brown if (strncmp(mattypename, "shell", 5)) {
1719566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, X, INSERT_VALUES));
1729566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, W1));
1739566063dSJacob Faibussowitsch PetscCall(VecView(W1, viewer));
174c4762a1bSJed Brown } else {
1759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MatDiagonalSet not tested on MATSHELL\n"));
176c4762a1bSJed Brown }
177c4762a1bSJed Brown
1789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E));
1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W1));
1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&W2));
181cf242a8cSBlanca Mellado Pinto PetscCall(VecDestroy(&W3));
1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown
main(int argc,char ** args)185d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
186d71ae5a4SJacob Faibussowitsch {
187c4762a1bSJed Brown const PetscInt inds[] = {0, 1};
188cf242a8cSBlanca Mellado Pinto #if defined(PETSC_USE_COMPLEX)
189cf242a8cSBlanca Mellado Pinto const PetscScalar xvals[] = {PetscCMPLX(11, 4), PetscCMPLX(13, 2)}, yvals[] = {PetscCMPLX(17, 3), PetscCMPLX(19, 1)}, zvals[] = {PetscCMPLX(23, 6), PetscCMPLX(29, 2)};
190cf242a8cSBlanca Mellado Pinto PetscScalar avals[] = {PetscCMPLX(2, 3), PetscCMPLX(3, 5), PetscCMPLX(5, 4), PetscCMPLX(7, 5)};
191cf242a8cSBlanca Mellado Pinto #else
192cf242a8cSBlanca Mellado Pinto const PetscScalar xvals[] = {11, 13}, yvals[] = {17, 19}, zvals[] = {23, 29};
193c4762a1bSJed Brown PetscScalar avals[] = {2, 3, 5, 7};
194cf242a8cSBlanca Mellado Pinto #endif
195c4762a1bSJed Brown Mat A, S, D[4], N;
196c4762a1bSJed Brown Vec X, Y, Z;
197c4762a1bSJed Brown User user;
198c4762a1bSJed Brown PetscInt i;
199c4762a1bSJed Brown
200327415f7SBarry Smith PetscFunctionBeginUser;
201c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
2029566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2, 2, 2, NULL, &A));
2039566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
2049566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 2, &X));
2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y));
2069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Z));
207cf242a8cSBlanca Mellado Pinto PetscCall(MatSetValues(A, 2, inds, 2, inds, avals, INSERT_VALUES));
2089566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 2, inds, xvals, INSERT_VALUES));
2099566063dSJacob Faibussowitsch PetscCall(VecSetValues(Y, 2, inds, yvals, INSERT_VALUES));
2109566063dSJacob Faibussowitsch PetscCall(VecSetValues(Z, 2, inds, zvals, INSERT_VALUES));
2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2139566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X));
2149566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y));
2159566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Z));
2169566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X));
2179566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y));
2189566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Z));
219c4762a1bSJed Brown
2209566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
221c4762a1bSJed Brown user->B = A;
222c4762a1bSJed Brown
2239566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, 2, 2, 2, 2, user, &S));
2249566063dSJacob Faibussowitsch PetscCall(MatSetUp(S));
225*57d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_VIEW, (PetscErrorCodeFn *)MatView_User));
226*57d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_MULT, (PetscErrorCodeFn *)MatMult_User));
227*57d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_User));
228*57d50842SBarry Smith PetscCall(MatShellSetOperation(S, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_User));
229c4762a1bSJed Brown
23048a46eb9SPierre Jolivet for (i = 0; i < 4; i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, 1, 1, &avals[i], &D[i]));
2319566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, D, &N));
2329566063dSJacob Faibussowitsch PetscCall(MatSetUp(N));
233c4762a1bSJed Brown
2349566063dSJacob Faibussowitsch PetscCall(TestMatrix(S, X, Y, Z));
2359566063dSJacob Faibussowitsch PetscCall(TestMatrix(A, X, Y, Z));
2369566063dSJacob Faibussowitsch PetscCall(TestMatrix(N, X, Y, Z));
237c4762a1bSJed Brown
2389566063dSJacob Faibussowitsch for (i = 0; i < 4; i++) PetscCall(MatDestroy(&D[i]));
2399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S));
2419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&N));
2429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Z));
2459566063dSJacob Faibussowitsch PetscCall(PetscFree(user));
2469566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
247b122ec5aSJacob Faibussowitsch return 0;
248c4762a1bSJed Brown }
249c4762a1bSJed Brown
250c4762a1bSJed Brown /*TEST
251c4762a1bSJed Brown
252cf242a8cSBlanca Mellado Pinto testset:
253c4762a1bSJed Brown test:
254cf242a8cSBlanca Mellado Pinto suffix: 1
255cf242a8cSBlanca Mellado Pinto requires:!complex
256cf242a8cSBlanca Mellado Pinto test:
257cf242a8cSBlanca Mellado Pinto suffix: 2
258cf242a8cSBlanca Mellado Pinto requires: complex
259c4762a1bSJed Brown
260c4762a1bSJed Brown TEST*/
261