xref: /petsc/src/mat/tests/ex88.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&user));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetType((PetscObject)A,&mattypename));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nMatrix of type: %s\n",mattypename));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&W1));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&W2));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,31));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,37));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(A,X,Y));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,41));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(A,Y,Z));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeOperator(A,MATDENSE,&E));
73c4762a1bSJed Brown 
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(E,viewer));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Testing MatMult + MatMultTranspose\n"));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,Z,W1));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(A,W1,W2));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
79c4762a1bSJed Brown 
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Testing MatMultAdd\n"));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multadd,&diff));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(W1,-1.0));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A,W1,W1,W2));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(W2,-1.0,diff));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(W2,-1.0));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(A,W1,W2,W2));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(W2,-1.0,diff));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(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
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&diff));
100c4762a1bSJed Brown 
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Testing MatMultTransposeAdd\n"));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,multtadd,&diff));
103c4762a1bSJed Brown 
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(W1,-1.0));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAdd(A,W1,W1,W2));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(W2,-1.0,diff));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(W2,-1.0));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAdd(A,W1,W2,W2));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(W2,-1.0,diff));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(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
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&diff));
122c4762a1bSJed Brown 
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Testing MatGetDiagonal\n"));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetDiagonal(A,W2));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(W2,viewer));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,2,diag,&diff));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(diff,-1.0,W2));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(diff,NORM_2,&nrm));
1292c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nrm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() produces incorrect result");
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&diff));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* MATSHELL does not support MatDiagonalSet after MatScale */
133c4762a1bSJed Brown   if (strncmp(mattypename, "shell", 5)) {
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDiagonalSet(A,X,INSERT_VALUES));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetDiagonal(A,W1));
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(W1,viewer));
137c4762a1bSJed Brown   } else {
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"MatDiagonalSet not tested on MATSHELL\n"));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&E));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&W1));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,2,&X));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&Y));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&Z));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(X,2,inds,xvals,INSERT_VALUES));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(Y,2,inds,yvals,INSERT_VALUES));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(Z,2,inds,zvals,INSERT_VALUES));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(X));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(Y));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(Z));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(X));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(Y));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(Z));
175c4762a1bSJed Brown 
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&user));
177c4762a1bSJed Brown   user->B = A;
178c4762a1bSJed Brown 
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(S));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_VIEW,(void (*)(void))MatView_User));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_User));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   for (i=0; i<4; i++) {
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSeqDense(PETSC_COMM_WORLD,1,1,&avals[i],&D[i]));
188c4762a1bSJed Brown   }
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,D,&N));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(N));
191c4762a1bSJed Brown 
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(TestMatrix(S,X,Y,Z));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(TestMatrix(A,X,Y,Z));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(TestMatrix(N,X,Y,Z));
195c4762a1bSJed Brown 
1965f80ce2aSJacob Faibussowitsch   for (i=0; i<4; i++) CHKERRQ(MatDestroy(&D[i]));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&S));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&N));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Y));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Z));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user));
204*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
205*b122ec5aSJacob Faibussowitsch   return 0;
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown /*TEST
209c4762a1bSJed Brown 
210c4762a1bSJed Brown    test:
211c4762a1bSJed Brown 
212c4762a1bSJed Brown TEST*/
213