Actual source code: ex120.c

  1: static char help[] = "Test LAPACK routine ZHEEV, ZHEEVX, ZHEGV and ZHEGVX. \n\
  2: ZHEEV computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. \n\n";

 4:  #include petscmat.h
 5:  #include petscblaslapack.h


 11: PetscInt main(PetscInt argc,char **args)
 12: {
 13:   Mat            A,A_dense,B;
 14:   Vec            *evecs;
 15:   PetscTruth     flg,TestZHEEV=PETSC_TRUE,TestZHEEVX=PETSC_FALSE,TestZHEGV=PETSC_FALSE,TestZHEGVX=PETSC_FALSE;
 17:   PetscTruth     isSymmetric;
 18:   PetscScalar    sigma,*arrayA,*arrayB,*evecs_array=PETSC_NULL,*work;
 19:   PetscReal      *evals,*rwork;
 20:   PetscMPIInt    size;
 21:   PetscInt       m,i,j,nevs,il,iu,cklvl=2;
 22:   PetscReal      vl,vu,abstol=1.e-8;
 23:   PetscBLASInt   *iwork,*ifail,lwork,lierr,bn;
 24:   PetscReal      tols[2];
 25:   PetscInt       nzeros[2],nz;
 26:   PetscReal      ratio;
 27:   PetscScalar    v,none = -1.0,sigma2,pfive = 0.5,*xa;
 28:   PetscRandom    rctx;
 29:   PetscReal      h2,sigma1 = 100.0;
 30:   PetscInt       dim,Ii,J,Istart,Iend,n = 6,its,use_random,one=1;
 31: 
 32:   PetscInitialize(&argc,&args,(char *)0,help);
 33: #if !defined(PETSC_USE_COMPLEX)
 34:   SETERRQ(1,"This example requires complex numbers");
 35: #endif
 36:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 37:   if (size != 1) SETERRQ(PETSC_ERR_SUP,"This is a uniprocessor example only!");

 39:   PetscOptionsHasName(PETSC_NULL, "-test_zheevx", &flg);
 40:   if (flg){
 41:     TestZHEEV  = PETSC_FALSE;
 42:     TestZHEEVX = PETSC_TRUE;
 43:   }
 44:   PetscOptionsHasName(PETSC_NULL, "-test_zhegv", &flg);
 45:   if (flg){
 46:     TestZHEEV  = PETSC_FALSE;
 47:     TestZHEGV= PETSC_TRUE;
 48:   }
 49:   PetscOptionsHasName(PETSC_NULL, "-test_zhegvx", &flg);
 50:   if (flg){
 51:     TestZHEEV  = PETSC_FALSE;
 52:     TestZHEGVX = PETSC_TRUE;
 53:   }

 55:   PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
 56:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 57:   dim = n*n;

 59:   MatCreate(PETSC_COMM_SELF,&A);
 60:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
 61:   MatSetType(A,MATSEQDENSE);
 62:   MatSetFromOptions(A);

 64:   PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
 65:   if (flg) use_random = 0;
 66:   else     use_random = 1;
 67:   if (use_random) {
 68:     PetscRandomCreate(PETSC_COMM_SELF,&rctx);
 69:     PetscRandomSetFromOptions(rctx);
 70:   } else {
 71:     sigma2 = 10.0*PETSC_i;
 72:   }
 73:   h2 = 1.0/((n+1)*(n+1));
 74:   for (Ii=0; Ii<dim; Ii++) {
 75:     v = -1.0; i = Ii/n; j = Ii - i*n;
 76:     if (i>0) {
 77:       J = Ii-n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 78:     if (i<n-1) {
 79:       J = Ii+n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 80:     if (j>0) {
 81:       J = Ii-1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 82:     if (j<n-1) {
 83:       J = Ii+1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 84:     if (use_random) {PetscRandomGetValueImaginary(rctx,&sigma2);}
 85:     v = 4.0 - sigma1*h2;
 86:     MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
 87:   }
 88:   /* make A complex Hermitian */
 89:   v = sigma2*h2;
 90:   Ii = 0; J = 1;
 91:   MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
 92:   v = -sigma2*h2;
 93:   MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);
 94:   if (use_random) {PetscRandomDestroy(rctx);}
 95:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 96:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 97:   m = n = dim;

 99:   /* Check whether A is symmetric */
100:   PetscOptionsHasName(PETSC_NULL, "-check_symmetry", &flg);
101:   if (flg) {
102:     Mat Trans;
103:     MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
104:     MatEqual(A, Trans, &isSymmetric);
105:     if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"A must be symmetric");
106:     MatDestroy(Trans);
107:   }

109:   /* Convert aij matrix to MatSeqDense for LAPACK */
110:   PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
111:   if (flg) {
112:     MatDuplicate(A,MAT_COPY_VALUES,&A_dense);
113:   } else {
114:     MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
115:   }

117:   MatCreate(PETSC_COMM_SELF,&B);
118:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
119:   MatSetType(B,MATSEQDENSE);
120:   MatSetFromOptions(B);
121:   v    = 1.0;
122:   for (Ii=0; Ii<dim; Ii++) {
123:     MatSetValues(B,1,&Ii,1,&Ii,&v,ADD_VALUES);
124:   }

126:   /* Solve standard eigenvalue problem: A*x = lambda*x */
127:   /*===================================================*/
128:   lwork = PetscBLASIntCast(2*n);
129:   bn    = PetscBLASIntCast(n);
130:   PetscMalloc(n*sizeof(PetscReal),&evals);
131:   PetscMalloc(lwork*sizeof(PetscScalar),&work);
132:   MatGetArray(A_dense,&arrayA);

134:   if (TestZHEEV){ /* test zheev() */
135:     printf(" LAPACKsyev: compute all %d eigensolutions...\n",m);
136:     PetscMalloc((3*n-2)*sizeof(PetscReal),&rwork);
137:     LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,rwork,&lierr);
138:     PetscFree(rwork);
139:     evecs_array = arrayA;
140:     nevs = m;
141:     il=1; iu=m;
142:   }
143:   if (TestZHEEVX){
144:     il = 1; iu=PetscBLASIntCast((0.2*m)); /* request 1 to 20%m evalues */
145:     printf(" LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu);
146:     PetscMalloc((m*n+1)*sizeof(PetscScalar),&evecs_array);
147:     PetscMalloc((7*n+1)*sizeof(PetscReal),&rwork);
148:     PetscMalloc((5*n+1)*sizeof(PetscBLASInt),&iwork);
149:     PetscMalloc((n+1)*sizeof(PetscBLASInt),&ifail);
150: 
151:     /* in the case "I", vl and vu are not referenced */
152:     vl = 0.0; vu = 8.0;
153:     LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&n,work,&lwork,rwork,iwork,ifail,&lierr);
154:     PetscFree(iwork);
155:     PetscFree(ifail);
156:     PetscFree(rwork);
157:   }
158:   if (TestZHEGV){
159:     printf(" LAPACKsygv: compute all %d eigensolutions...\n",m);
160:     PetscMalloc((3*n+1)*sizeof(PetscReal),&rwork);
161:     MatGetArray(B,&arrayB);
162:     LAPACKsygv_(&one,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,rwork,&lierr);
163:     evecs_array = arrayA;
164:     nevs = m;
165:     il=1; iu=m;
166:     MatRestoreArray(B,&arrayB);
167:     PetscFree(rwork);
168:   }
169:   if (TestZHEGVX){
170:     il = 1; iu=PetscBLASIntCast((0.2*m)); /* request 1 to 20%m evalues */
171:     printf(" LAPACKsygv: compute %d to %d-th eigensolutions...\n",il,iu);
172:     PetscMalloc((m*n+1)*sizeof(PetscScalar),&evecs_array);
173:     PetscMalloc((6*n+1)*sizeof(PetscBLASInt),&iwork);
174:     ifail = iwork + 5*n;
175:     PetscMalloc((7*n+1)*sizeof(PetscReal),&rwork);
176:     MatGetArray(B,&arrayB);
177:     vl = 0.0; vu = 8.0;
178:     LAPACKsygvx_(&one,"V","I","U",&bn,arrayA,&bn,arrayB,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&n,work,&lwork,rwork,iwork,ifail,&lierr);
179:     MatRestoreArray(B,&arrayB);
180:     PetscFree(iwork);
181:     PetscFree(rwork);
182:   }
183:   MatRestoreArray(A_dense,&arrayA);
184:   if (nevs <= 0 ) SETERRQ1(PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);

186:   /* View evals */
187:   PetscOptionsHasName(PETSC_NULL, "-eig_view", &flg);
188:   if (flg){
189:     printf(" %d evals: \n",nevs);
190:     for (i=0; i<nevs; i++) printf("%d  %G\n",i+il,evals[i]);
191:   }

193:   /* Check residuals and orthogonality */
194:   PetscMalloc((nevs+1)*sizeof(Vec),&evecs);
195:   for (i=0; i<nevs; i++){
196:     VecCreate(PETSC_COMM_SELF,&evecs[i]);
197:     VecSetSizes(evecs[i],PETSC_DECIDE,n);
198:     VecSetFromOptions(evecs[i]);
199:     VecPlaceArray(evecs[i],evecs_array+i*n);
200:   }
201: 
202:   tols[0] = 1.e-8;  tols[1] = 1.e-8;
203:   CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols);
204:   for (i=0; i<nevs; i++){ VecDestroy(evecs[i]);}
205:   PetscFree(evecs);
206: 
207:   /* Free work space. */
208:   if (TestZHEEVX || TestZHEGVX){
209:     PetscFree(evecs_array);
210:   }
211:   PetscFree(evals);
212:   PetscFree(work);
213:   MatDestroy(A_dense);
214:   MatDestroy(A);
215:   MatDestroy(B);
216:   PetscFinalize();
217:   return 0;
218: }
219: /*------------------------------------------------
220:   Check the accuracy of the eigen solution
221:   ----------------------------------------------- */
222: /*
223:   input: 
224:      cklvl      - check level: 
225:                     1: check residual
226:                     2: 1 and check B-orthogonality locally 
227:      A          - matrix 
228:      il,iu      - lower and upper index bound of eigenvalues 
229:      eval, evec - eigenvalues and eigenvectors stored in this process
230:      tols[0]    - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
231:      tols[1]    - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
232: */
233: #undef DEBUG_CkEigenSolutions
236: PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
237: {
238:   PetscInt     ierr,i,j,nev;
239:   Vec          vt1,vt2; /* tmp vectors */
240:   PetscReal    norm,tmp,norm_max,dot_max,rdot;
241:   PetscScalar  dot;

244:   nev = iu - il;
245:   if (nev <= 0) return(0);

247:   //VecView(evec[0],PETSC_VIEWER_STDOUT_SELF);
248:   VecDuplicate(evec[0],&vt1);
249:   VecDuplicate(evec[0],&vt2);

251:   switch (cklvl){
252:   case 2:
253:     dot_max = 0.0;
254:     for (i = il; i<iu; i++){
255:       //printf("ck %d-th\n",i);
256:       VecCopy(evec[i], vt1);
257:       for (j=il; j<iu; j++){
258:         VecDot(evec[j],vt1,&dot);
259:         if (j == i){
260:           rdot = PetscAbsScalar(dot - 1.0);
261:         } else {
262:           rdot = PetscAbsScalar(dot);
263:         }
264:         if (rdot > dot_max) dot_max = rdot;
265: #ifdef DEBUG_CkEigenSolutions
266:         if (rdot > tols[1] ) {
267:           VecNorm(evec[i],NORM_INFINITY,&norm);
268:           PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %G, norm: %G\n",i,j,dot,norm);
269:         }
270: #endif
271:       }
272:     }
273:     PetscPrintf(PETSC_COMM_SELF,"    max|(x_j^T*x_i) - delta_ji|: %G\n",dot_max);

275:   case 1:
276:     norm_max = 0.0;
277:     for (i = il; i< iu; i++){
278:       MatMult(A, evec[i], vt1);
279:       VecCopy(evec[i], vt2);
280:       tmp  = -eval[i];
281:       VecAXPY(vt1,tmp,vt2);
282:       VecNorm(vt1, NORM_INFINITY, &norm);
283:       norm = PetscAbsScalar(norm);
284:       if (norm > norm_max) norm_max = norm;
285: #ifdef DEBUG_CkEigenSolutions
286:       /* sniff, and bark if necessary */
287:       if (norm > tols[0]){
288:         printf( "  residual violation: %d, resi: %g\n",i, norm);
289:       }
290: #endif
291:     }
292:     PetscPrintf(PETSC_COMM_SELF,"    max_resi:                    %G\n", norm_max);
293:    break;
294:   default:
295:     PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl);
296:   }
297:   VecDestroy(vt2);
298:   VecDestroy(vt1);
299:   return(0);
300: }