Actual source code: ex127.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: static char help[] = "Test MatMult() for Hermitian matrix.\n\n";
  2: /*
  3:   Example of usage
  4:     ./ex127 -check_Hermitian -display_mat -display_vec
  5:     mpiexec -n 2 ./ex127
  6: */

  8: #include <petscmat.h>

 12: int main(int argc,char **args)
 13: {
 14:   Mat            A,As;
 15:   Vec            x,y,ys;
 16:   PetscBool      flg,disp_mat=PETSC_FALSE,disp_vec=PETSC_FALSE;
 18:   PetscMPIInt    size,rank;
 19:   PetscInt       m,i,j;
 20:   PetscScalar    v,sigma2;
 21:   PetscRandom    rctx;
 22:   PetscReal      h2,sigma1=100.0,norm;
 23:   PetscInt       dim,Ii,J,n = 3,use_random,rstart,rend;

 25:   PetscInitialize(&argc,&args,(char*)0,help);
 26: #if !defined(PETSC_USE_COMPLEX)
 27:   SETERRQ(PETSC_COMM_WORLD,1,"This example requires complex numbers");
 28: #endif
 29:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 30:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 31:   PetscOptionsHasName(NULL, "-display_mat", &disp_mat);
 32:   PetscOptionsHasName(NULL, "-display_vec", &disp_vec);

 34:   PetscOptionsGetReal(NULL,"-sigma1",&sigma1,NULL);
 35:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 36:   dim  = n*n;

 38:   MatCreate(PETSC_COMM_WORLD,&A);
 39:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
 40:   MatSetType(A,MATAIJ);
 41:   MatSetFromOptions(A);
 42:   MatSetUp(A);

 44:   PetscOptionsHasName(NULL,"-norandom",&flg);
 45:   if (flg) use_random = 0;
 46:   else use_random = 1;
 47:   if (use_random) {
 48:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 49:     PetscRandomSetFromOptions(rctx);
 50:     PetscRandomSetInterval(rctx,0.0,PETSC_i);
 51:     PetscRandomGetValue(rctx,&sigma2); /* RealPart(sigma2) == 0.0 */
 52:   } else {
 53:     sigma2 = 10.0*PETSC_i;
 54:   }
 55:   h2 = 1.0/((n+1)*(n+1));

 57:   MatGetOwnershipRange(A,&rstart,&rend);
 58:   for (Ii=rstart; Ii<rend; Ii++) {
 59:     v = -1.0; i = Ii/n; j = Ii - i*n;
 60:     if (i>0) {
 61:       J = Ii-n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 62:     }
 63:     if (i<n-1) {
 64:       J = Ii+n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 65:     }
 66:     if (j>0) {
 67:       J = Ii-1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 68:     }
 69:     if (j<n-1) {
 70:       J = Ii+1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 71:     }
 72:     v    = 4.0 - sigma1*h2;
 73:     MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
 74:   }
 75:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 76:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 78:   /* Check whether A is symmetric */
 79:   PetscOptionsHasName(NULL, "-check_symmetric", &flg);
 80:   if (flg) {
 81:     Mat Trans;
 82:     MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
 83:     MatEqual(A, Trans, &flg);
 84:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A is not symmetric");
 85:     MatDestroy(&Trans);
 86:   }
 87:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);

 89:   /* make A complex Hermitian */
 90:   Ii = 0; J = dim-1;
 91:   if (Ii >= rstart && Ii < rend) {
 92:     v    = sigma2*h2; /* RealPart(v) = 0.0 */
 93:     MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 94:     v    = -sigma2*h2;
 95:     MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);
 96:   }

 98:   Ii = dim-2; J = dim-1;
 99:   if (Ii >= rstart && Ii < rend) {
100:     v    = sigma2*h2; /* RealPart(v) = 0.0 */
101:     MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
102:     v    = -sigma2*h2;
103:     MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);
104:   }

106:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
107:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

109:   /* Check whether A is Hermitian */
110:   PetscOptionsHasName(NULL, "-check_Hermitian", &flg);
111:   if (flg) {
112:     Mat Hermit;
113:     if (disp_mat) {
114:       if (!rank) printf(" A:\n");
115:       MatView(A,PETSC_VIEWER_STDOUT_WORLD);
116:     }
117:     MatHermitianTranspose(A,MAT_INITIAL_MATRIX, &Hermit);
118:     if (disp_mat) {
119:       if (!rank) printf(" A_Hermitian:\n");
120:       MatView(Hermit,PETSC_VIEWER_STDOUT_WORLD);
121:     }
122:     MatEqual(A, Hermit, &flg);
123:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A is not Hermitian");
124:     MatDestroy(&Hermit);
125:   }
126:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

128:   /* Create a Hermitian matrix As in sbaij format */
129:   MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&As);
130:   if (disp_mat) {
131:     if (!rank) {PetscPrintf(PETSC_COMM_SELF," As:\n");}
132:     MatView(As,PETSC_VIEWER_STDOUT_WORLD);
133:   }

135:   MatGetLocalSize(A,&m,&n);
136:   VecCreate(PETSC_COMM_WORLD,&x);
137:   VecSetSizes(x,n,PETSC_DECIDE);
138:   VecSetFromOptions(x);
139:   if (use_random) {
140:     VecSetRandom(x,rctx);
141:   } else {
142:     VecSet(x,1.0);
143:   }

145:   /* Create vectors */
146:   VecCreate(PETSC_COMM_WORLD,&y);
147:   VecSetSizes(y,m,PETSC_DECIDE);
148:   VecSetFromOptions(y);
149:   VecDuplicate(y,&ys);

151:   /* Test MatMult */
152:   MatMult(A,x,y);
153:   MatMult(As,x,ys);
154:   if (disp_vec) {
155:     printf("y = A*x:\n");
156:     VecView(y,PETSC_VIEWER_STDOUT_WORLD);
157:     PetscPrintf(PETSC_COMM_WORLD,"ys = As*x:\n");
158:     VecView(ys,PETSC_VIEWER_STDOUT_WORLD);
159:   }
160:   VecAXPY(y,-1.0,ys);
161:   VecNorm(y,NORM_INFINITY,&norm);
162:   if (norm > 1.e-12 || disp_vec) {
163:     printf("|| A*x - As*x || = %g\n",(double)norm);
164:   }

166:   /* Free spaces */
167:   if (use_random) {PetscRandomDestroy(&rctx);}
168:   MatDestroy(&A);
169:   MatDestroy(&As);

171:   VecDestroy(&x);
172:   VecDestroy(&y);
173:   VecDestroy(&ys);
174:   PetscFinalize();
175:   return 0;
176: }