Actual source code: ex88.c

petsc-3.9.1 2018-04-29
Report Typos and Errors

  2: static char help[] = "Tests MatShift(), MatScale(), and MatDiagonalScale() for SHELL and NEST matrices\n\n";

  4:  #include <petscmat.h>

  6: typedef struct _n_User *User;
  7: struct _n_User {
  8:   Mat B;
  9: };

 11: static PetscErrorCode MatView_User(Mat A,PetscViewer viewer)
 12: {
 13:   User           user;

 17:   MatShellGetContext(A,&user);
 18:   MatView(user->B,viewer);
 19:   return(0);
 20: }

 22: static PetscErrorCode MatMult_User(Mat A,Vec X,Vec Y)
 23: {
 24:   User           user;

 28:   MatShellGetContext(A,&user);
 29:   MatMult(user->B,X,Y);
 30:   return(0);
 31: }

 33: static PetscErrorCode MatMultTranspose_User(Mat A,Vec X,Vec Y)
 34: {
 35:   User           user;

 39:   MatShellGetContext(A,&user);
 40:   MatMultTranspose(user->B,X,Y);
 41:   return(0);
 42: }

 44: static PetscErrorCode TestMatrix(Mat A,Vec X,Vec Y,Vec Z)
 45: {
 47:   Vec            W1,W2;
 48:   Mat            E;
 49:   const char     *mattypename;
 50:   PetscViewer    viewer = PETSC_VIEWER_STDOUT_WORLD;

 53:   VecDuplicate(X,&W1);
 54:   VecDuplicate(X,&W2);
 55:   MatScale(A,31);
 56:   MatShift(A,37);
 57:   MatDiagonalScale(A,X,Y);
 58:   MatScale(A,41);
 59:   MatDiagonalScale(A,Y,Z);
 60:   MatComputeExplicitOperator(A,&E);

 62:   PetscObjectGetType((PetscObject)A,&mattypename);
 63:   PetscViewerASCIIPrintf(viewer,"Matrix of type: %s\n",mattypename);
 64:   MatView(E,viewer);
 65:   MatMult(A,Z,W1);
 66:   MatMultTranspose(A,W1,W2);
 67:   VecView(W2,viewer);
 68:   MatGetDiagonal(A,W2);
 69:   VecView(W2,viewer);
 70:   /* MATSHELL does not support MatDiagonalSet after MatScale */
 71:   if (strncmp(mattypename, "shell", 5)) {
 72:     MatDiagonalSet(A,X,INSERT_VALUES);
 73:     MatGetDiagonal(A,W1);
 74:     VecView(W1,viewer);
 75:   } else {
 76:     PetscViewerASCIIPrintf(viewer,"MatDiagonalSet not tested on MATSHELL\n");
 77:   }
 78:   MatDestroy(&E);
 79:   VecDestroy(&W1);
 80:   VecDestroy(&W2);
 81:   return(0);
 82: }

 84: int main(int argc,char **args)
 85: {
 86:   const PetscScalar xvals[] = {11,13},yvals[] = {17,19},zvals[] = {23,29};
 87:   const PetscInt    inds[]  = {0,1};
 88:   PetscScalar       avals[] = {2,3,5,7};
 89:   Mat               A,S,D[4],N;
 90:   Vec               X,Y,Z;
 91:   User              user;
 92:   PetscInt          i;
 93:   PetscErrorCode    ierr;

 95:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 96:   MatCreateSeqAIJ(PETSC_COMM_WORLD,2,2,2,NULL,&A);
 97:   MatSetUp(A);
 98:   MatSetValues(A,2,inds,2,inds,avals,INSERT_VALUES);
 99:   VecCreateSeq(PETSC_COMM_WORLD,2,&X);
100:   VecDuplicate(X,&Y);
101:   VecDuplicate(X,&Z);
102:   VecSetValues(X,2,inds,xvals,INSERT_VALUES);
103:   VecSetValues(Y,2,inds,yvals,INSERT_VALUES);
104:   VecSetValues(Z,2,inds,zvals,INSERT_VALUES);
105:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
106:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
107:   VecAssemblyBegin(X);
108:   VecAssemblyBegin(Y);
109:   VecAssemblyBegin(Z);
110:   VecAssemblyEnd(X);
111:   VecAssemblyEnd(Y);
112:   VecAssemblyEnd(Z);

114:   PetscNew(&user);
115:   user->B = A;

117:   MatCreateShell(PETSC_COMM_WORLD,2,2,2,2,user,&S);
118:   MatSetUp(S);
119:   MatShellSetOperation(S,MATOP_VIEW,(void (*)(void))MatView_User);
120:   MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MatMult_User);
121:   MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_User);

123:   for (i=0; i<4; i++) {
124:     MatCreateSeqDense(PETSC_COMM_WORLD,1,1,&avals[i],&D[i]);
125:   }
126:   MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,D,&N);
127:   MatSetUp(N);

129:   TestMatrix(S,X,Y,Z);
130:   TestMatrix(A,X,Y,Z);
131:   TestMatrix(N,X,Y,Z);

133:   for (i=0; i<4; i++) {MatDestroy(&D[i]);}
134:   MatDestroy(&A);
135:   MatDestroy(&S);
136:   MatDestroy(&N);
137:   VecDestroy(&X);
138:   VecDestroy(&Y);
139:   VecDestroy(&Z);
140:   PetscFree(user);
141:   PetscFinalize();
142:   return ierr;
143: }


146: /*TEST

148:    test:

150: TEST*/