xref: /petsc/src/mat/tests/ex88.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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;
169566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&user));
179566063dSJacob 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;
269566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&user));
279566063dSJacob 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;
369566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&user));
379566063dSJacob 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;
469566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A,&user));
479566063dSJacob 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;
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject)A,&mattypename));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nMatrix of type: %s\n",mattypename));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&W1));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&W2));
679566063dSJacob Faibussowitsch   PetscCall(MatScale(A,31));
689566063dSJacob Faibussowitsch   PetscCall(MatShift(A,37));
699566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,X,Y));
709566063dSJacob Faibussowitsch   PetscCall(MatScale(A,41));
719566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,Y,Z));
729566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(A,MATDENSE,&E));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(MatView(E,viewer));
759566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Testing MatMult + MatMultTranspose\n"));
769566063dSJacob Faibussowitsch   PetscCall(MatMult(A,Z,W1));
779566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(A,W1,W2));
789566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Testing MatMultAdd\n"));
819566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multadd,&diff));
829566063dSJacob Faibussowitsch   PetscCall(VecSet(W1,-1.0));
839566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A,W1,W1,W2));
849566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W2,-1.0,diff));
869566063dSJacob Faibussowitsch   PetscCall(VecNorm(W2,NORM_2,&nrm));
87c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
88*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultAdd(A,x,x,y) produces incorrect result");
89c4762a1bSJed Brown #endif
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(VecSet(W2,-1.0));
929566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A,W1,W2,W2));
939566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
949566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W2,-1.0,diff));
959566063dSJacob Faibussowitsch   PetscCall(VecNorm(W2,NORM_2,&nrm));
96c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
97*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultAdd(A,x,y,y) produces incorrect result");
98c4762a1bSJed Brown #endif
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diff));
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Testing MatMultTransposeAdd\n"));
1029566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multtadd,&diff));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(VecSet(W1,-1.0));
1059566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A,W1,W1,W2));
1069566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
1079566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W2,-1.0,diff));
1089566063dSJacob Faibussowitsch   PetscCall(VecNorm(W2,NORM_2,&nrm));
109c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
110*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultTransposeAdd(A,x,x,y) produces incorrect result");
111c4762a1bSJed Brown #endif
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch   PetscCall(VecSet(W2,-1.0));
1149566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(A,W1,W2,W2));
1159566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
1169566063dSJacob Faibussowitsch   PetscCall(VecAXPY(W2,-1.0,diff));
1179566063dSJacob Faibussowitsch   PetscCall(VecNorm(W2,NORM_2,&nrm));
118c4762a1bSJed Brown #if defined(PETSC_USE_REAL_DOUBLE) || defined(PETSC_USE_REAL___FLOAT128)
119*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatMultTransposeAdd(A,x,y,y) produces incorrect result");
120c4762a1bSJed Brown #endif
1219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diff));
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Testing MatGetDiagonal\n"));
1249566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A,W2));
1259566063dSJacob Faibussowitsch   PetscCall(VecView(W2,viewer));
1269566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,diag,&diff));
1279566063dSJacob Faibussowitsch   PetscCall(VecAXPY(diff,-1.0,W2));
1289566063dSJacob Faibussowitsch   PetscCall(VecNorm(diff,NORM_2,&nrm));
129*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() produces incorrect result");
1309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diff));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* MATSHELL does not support MatDiagonalSet after MatScale */
133c4762a1bSJed Brown   if (strncmp(mattypename, "shell", 5)) {
1349566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(A,X,INSERT_VALUES));
1359566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(A,W1));
1369566063dSJacob Faibussowitsch     PetscCall(VecView(W1,viewer));
137c4762a1bSJed Brown   } else {
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"MatDiagonalSet not tested on MATSHELL\n"));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&E));
1429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W1));
1439566063dSJacob 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 
1579566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
1589566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A));
1599566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES));
1619566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,2,&X));
1629566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&Y));
1639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&Z));
1649566063dSJacob Faibussowitsch   PetscCall(VecSetValues(X,2,inds,xvals,INSERT_VALUES));
1659566063dSJacob Faibussowitsch   PetscCall(VecSetValues(Y,2,inds,yvals,INSERT_VALUES));
1669566063dSJacob Faibussowitsch   PetscCall(VecSetValues(Z,2,inds,zvals,INSERT_VALUES));
1679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1699566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
1709566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
1719566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Z));
1729566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
1739566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
1749566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Z));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
177c4762a1bSJed Brown   user->B = A;
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch   PetscCall(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S));
1809566063dSJacob Faibussowitsch   PetscCall(MatSetUp(S));
1819566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S,MATOP_VIEW,(void (*)(void))MatView_User));
1829566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User));
1839566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User));
1849566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   for (i=0; i<4; i++) {
1879566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD,1,1,&avals[i],&D[i]));
188c4762a1bSJed Brown   }
1899566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,D,&N));
1909566063dSJacob Faibussowitsch   PetscCall(MatSetUp(N));
191c4762a1bSJed Brown 
1929566063dSJacob Faibussowitsch   PetscCall(TestMatrix(S,X,Y,Z));
1939566063dSJacob Faibussowitsch   PetscCall(TestMatrix(A,X,Y,Z));
1949566063dSJacob Faibussowitsch   PetscCall(TestMatrix(N,X,Y,Z));
195c4762a1bSJed Brown 
1969566063dSJacob Faibussowitsch   for (i=0; i<4; i++) PetscCall(MatDestroy(&D[i]));
1979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&S));
1999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&N));
2009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
2019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
2029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Z));
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(user));
2049566063dSJacob 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