Actual source code: ex127.c

petsc-3.4.5 2014-06-29
  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: PetscInt main(PetscInt 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);

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

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

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

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

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

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

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

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

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

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

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

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

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