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