Actual source code: ex99.c

  1: static char help[] = "Test LAPACK routine DSYGV() or DSYGVX(). \n\
  2: Reads PETSc matrix A and B (or create B=I), \n\
  3: then computes selected eigenvalues, and optionally, eigenvectors of \n\
  4: a real generalized symmetric-definite eigenproblem \n\
  5:  A*x = lambda*B*x \n\
  6: Input parameters include\n\
  7:   -f0 <input_file> : first file to load (small system)\n\
  8:   -fA <input_file> -fB <input_file>: second files to load (larger system) \n\
  9: e.g. ./ex99 -f0 $D/small -fA $D/Eigdftb/dftb_bin/diamond_xxs_A -fB $D/Eigdftb/dftb_bin/diamond_xxs_B -mat_getrow_uppertriangular,\n\
 10:      where $D = /home/petsc/datafiles/matrices/Eigdftb/dftb_bin\n\n";

 12:  #include petscmat.h
 13:  #include ../src/mat/impls/sbaij/seq/sbaij.h
 14:  #include petscblaslapack.h


 20: PetscInt main(PetscInt argc,char **args)
 21: {
 22:   Mat            A,B,A_dense,B_dense,mats[2],A_sp;
 23:   Vec            *evecs;
 24:   PetscViewer    fd;                /* viewer */
 25:   char           file[3][PETSC_MAX_PATH_LEN];     /* input file name */
 26:   PetscTruth     flg,flgA=PETSC_FALSE,flgB=PETSC_FALSE,TestSYGVX=PETSC_TRUE;
 28:   PetscTruth     preload=PETSC_TRUE,isSymmetric;
 29:   PetscScalar    sigma,one=1.0,*arrayA,*arrayB,*evecs_array,*work,*evals;
 30:   PetscMPIInt    size;
 31:   PetscInt       m,n,i,j,nevs,il,iu;
 32:   PetscLogStage  stages[2];
 33:   PetscReal      vl,vu,abstol=1.e-8;
 34:   PetscBLASInt   *iwork,*ifail,lone=1,lwork,lierr,bn;
 35:   PetscInt       ievbd_loc[2],offset=0,cklvl=2;
 36:   PetscReal      tols[2];
 37:   Mat_SeqSBAIJ   *sbaij;
 38:   PetscScalar    *aa;
 39:   PetscInt       *ai,*aj;
 40:   PetscInt       nzeros[2],nz;
 41:   PetscReal      ratio;
 42: 
 43:   PetscInitialize(&argc,&args,(char *)0,help);
 44:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 45:   if (size != 1) SETERRQ(PETSC_ERR_SUP,"This is a uniprocessor example only!");
 46:   PetscLogStageRegister("EigSolve",&stages[0]);
 47:   PetscLogStageRegister("EigCheck",&stages[1]);

 49:   /* Determine files from which we read the two matrices */
 50:   PetscOptionsGetString(PETSC_NULL,"-f0",file[0],PETSC_MAX_PATH_LEN-1,&flg);
 51:   if (!flg) {
 52:     PetscOptionsGetString(PETSC_NULL,"-fA",file[0],PETSC_MAX_PATH_LEN-1,&flgA);
 53:     if (!flgA) SETERRQ(PETSC_ERR_USER,"Must indicate binary file with the -fA or -fB options");
 54:     PetscOptionsGetString(PETSC_NULL,"-fB",file[1],PETSC_MAX_PATH_LEN-1,&flgB);
 55:     preload = PETSC_FALSE;
 56:   } else {
 57:     PetscOptionsGetString(PETSC_NULL,"-fA",file[1],PETSC_MAX_PATH_LEN-1,&flgA);
 58:     if (!flgA) {preload = PETSC_FALSE;} /* don't bother with second system */
 59:     PetscOptionsGetString(PETSC_NULL,"-fB",file[2],PETSC_MAX_PATH_LEN-1,&flgB);
 60:   }

 62:   PreLoadBegin(preload,"Load system");
 63:     /* Load matrices */
 64:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt],FILE_MODE_READ,&fd);
 65:     MatLoad(fd,MATSBAIJ,&A);
 66:     PetscViewerDestroy(fd);
 67:     MatGetSize(A,&m,&n);
 68:     if ((flgB && PreLoadIt) || (flgB && !preload)){
 69:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt+1],FILE_MODE_READ,&fd);
 70:       MatLoad(fd,MATSBAIJ,&B);
 71:       PetscViewerDestroy(fd);
 72:     } else { /* create B=I */
 73:       MatCreate(PETSC_COMM_WORLD,&B);
 74:       MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
 75:       MatSetType(B,MATSEQSBAIJ);
 76:       MatSetFromOptions(B);
 77:       for (i=0; i<m; i++) {
 78:         MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES);
 79:       }
 80:       MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 81:       MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 82:     }
 83: 
 84:     /* Add a shift to A */
 85:     PetscOptionsGetScalar(PETSC_NULL,"-mat_sigma",&sigma,&flg);
 86:     if(flg) {
 87:       MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN); /* A <- sigma*B + A */
 88:     }

 90:     /* Check whether A is symmetric */
 91:     PetscOptionsHasName(PETSC_NULL, "-check_symmetry", &flg);
 92:     if (flg) {
 93:       Mat Trans;
 94:       MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
 95:       MatEqual(A, Trans, &isSymmetric);
 96:       if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"A must be symmetric");
 97:       MatDestroy(Trans);
 98:       if (flgB && PreLoadIt){
 99:         MatTranspose(B,MAT_INITIAL_MATRIX, &Trans);
100:         MatEqual(B, Trans, &isSymmetric);
101:         if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"B must be symmetric");
102:         MatDestroy(Trans);
103:       }
104:     }

106:     /* View small entries of A */
107:     PetscOptionsHasName(PETSC_NULL, "-Asp_view", &flg);
108:     if (flg){
109:       MatCreate(PETSC_COMM_SELF,&A_sp);
110:       MatSetSizes(A_sp,PETSC_DECIDE,PETSC_DECIDE,m,n);
111:       MatSetType(A_sp,MATSEQSBAIJ);

113:       tols[0] = 1.e-6, tols[1] = 1.e-9;
114:       sbaij = (Mat_SeqSBAIJ*)A->data;
115:       ai    = sbaij->i;
116:       aj    = sbaij->j;
117:       aa    = sbaij->a;
118:       nzeros[0] = nzeros[1] = 0;
119:       for (i=0; i<m; i++) {
120:         nz = ai[i+1] - ai[i];
121:         for (j=0; j<nz; j++){
122:           if (PetscAbsScalar(*aa)<tols[0]) {
123:             MatSetValues(A_sp,1,&i,1,aj,aa,INSERT_VALUES);
124:             nzeros[0]++;
125:           }
126:           if (PetscAbsScalar(*aa)<tols[1]) nzeros[1]++;
127:           aa++; aj++;
128:         }
129:       }
130:       MatAssemblyBegin(A_sp,MAT_FINAL_ASSEMBLY);
131:       MatAssemblyEnd(A_sp,MAT_FINAL_ASSEMBLY);

133:       MatDestroy(A_sp);

135:       ratio = (PetscReal)nzeros[0]/sbaij->nz;
136:       PetscPrintf(PETSC_COMM_SELF," %d matrix entries < %e, ratio %G of %d nonzeros\n",nzeros[0],tols[0],ratio,sbaij->nz);
137:       PetscPrintf(PETSC_COMM_SELF," %d matrix entries < %e\n",nzeros[1],tols[1]);
138:     }

140:     /* Convert aij matrix to MatSeqDense for LAPACK */
141:     PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
142:     if (!flg) {
143:       MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
144:     }
145:     PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
146:     if (!flg) {MatConvert(B,MATSEQDENSE,MAT_INITIAL_MATRIX,&B_dense);}

148:     /* Solve eigenvalue problem: A*x = lambda*B*x */
149:     /*============================================*/
150:     lwork = PetscBLASIntCast(8*n);
151:     bn    = PetscBLASIntCast(n);
152:     PetscMalloc(n*sizeof(PetscScalar),&evals);
153:     PetscMalloc(lwork*sizeof(PetscScalar),&work);
154:     MatGetArray(A_dense,&arrayA);
155:     MatGetArray(B_dense,&arrayB);

157:     if (!TestSYGVX){ /* test sygv()  */
158:       evecs_array = arrayA;
159:       LAPACKsygv_(&lone,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,&lierr);
160:       nevs = m;
161:       il=1;
162:     } else { /* test sygvx()  */
163:       il = 1; iu=PetscBLASIntCast(.6*m); /* request 1 to 60%m evalues */
164:       PetscMalloc((m*n+1)*sizeof(PetscScalar),&evecs_array);
165:       PetscMalloc((6*n+1)*sizeof(PetscBLASInt),&iwork);
166:       ifail = iwork + 5*n;
167:       if(PreLoadIt){PetscLogStagePush(stages[0]);}
168:       /* in the case "I", vl and vu are not referenced */
169:       LAPACKsygvx_(&lone,"V","I","U",&bn,arrayA,&bn,arrayB,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&n,work,&lwork,iwork,ifail,&lierr);
170:       if(PreLoadIt){PetscLogStagePop();}
171:       PetscFree(iwork);
172:     }
173:     MatRestoreArray(A,&arrayA);
174:     MatRestoreArray(B,&arrayB);

176:     if (nevs <= 0 ) SETERRQ1(PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
177:     /* View evals */
178:     PetscOptionsHasName(PETSC_NULL, "-eig_view", &flg);
179:     if (flg){
180:       printf(" %d evals: \n",nevs);
181:       for (i=0; i<nevs; i++) printf("%d  %G\n",i+il,evals[i]);
182:     }

184:     /* Check residuals and orthogonality */
185:     if(PreLoadIt){
186:       mats[0] = A; mats[1] = B;
187:       one = (PetscInt)one;
188:       PetscMalloc((nevs+1)*sizeof(Vec),&evecs);
189:       for (i=0; i<nevs; i++){
190:         VecCreate(PETSC_COMM_SELF,&evecs[i]);
191:         VecSetSizes(evecs[i],PETSC_DECIDE,n);
192:         VecSetFromOptions(evecs[i]);
193:         VecPlaceArray(evecs[i],evecs_array+i*n);
194:       }
195: 
196:       ievbd_loc[0] = 0; ievbd_loc[1] = nevs-1;
197:       tols[0] = 1.e-8;  tols[1] = 1.e-8;
198:       PetscLogStagePush(stages[1]);
199:       CkEigenSolutions(&cklvl,mats,evals,evecs,ievbd_loc,&offset,tols);
200:       PetscLogStagePop();
201:       for (i=0; i<nevs; i++){ VecDestroy(evecs[i]);}
202:       PetscFree(evecs);
203:     }
204: 
205:     /* Free work space. */
206:     if (TestSYGVX){PetscFree(evecs_array);}
207: 
208:     PetscFree(evals);
209:     PetscFree(work);

211:     MatDestroy(A_dense);
212:     MatDestroy(B_dense);
213:     MatDestroy(B);
214:     MatDestroy(A);

216:   PreLoadEnd();
217:   PetscFinalize();
218:   return 0;
219: }
220: /*------------------------------------------------
221:   Check the accuracy of the eigen solution
222:   ----------------------------------------------- */
223: /*
224:   input: 
225:      cklvl      - check level: 
226:                     1: check residual
227:                     2: 1 and check B-orthogonality locally 
228:      mats       - matrix pencil
229:      eval, evec - eigenvalues and eigenvectors stored in this process
230:      ievbd_loc  - local eigenvalue bounds, see eigc()
231:      offset     - see eigc()
232:      tols[0]    - reporting tol_res: || A evec[i] - eval[i] B evec[i]||
233:      tols[1]    - reporting tol_orth: evec[i] B evec[j] - delta_ij
234: */
235: #undef DEBUG_CkEigenSolutions
238: PetscErrorCode CkEigenSolutions(PetscInt *fcklvl,Mat *mats,
239:                    PetscReal *eval,Vec *evec,PetscInt *ievbd_loc,PetscInt *offset, 
240:                    PetscReal *tols)
241: {
242:   PetscInt     ierr,cklvl=*fcklvl,nev_loc,i,j;
243:   Mat          A=mats[0], B=mats[1];
244:   Vec          vt1,vt2; /* tmp vectors */
245:   PetscReal    norm,tmp,dot,norm_max,dot_max;

248:   nev_loc = ievbd_loc[1] - ievbd_loc[0];
249:   if (nev_loc == 0) return(0);

251:   nev_loc += (*offset);
252:   VecDuplicate(evec[*offset],&vt1);
253:   VecDuplicate(evec[*offset],&vt2);

255:   switch (cklvl){
256:   case 2:
257:     dot_max = 0.0;
258:     for (i = *offset; i<nev_loc; i++){
259:       MatMult(B, evec[i], vt1);
260:       for (j=i; j<nev_loc; j++){
261:         VecDot(evec[j],vt1,&dot);
262:         if (j == i){
263:           dot = PetscAbsScalar(dot - 1.0);
264:         } else {
265:           dot = PetscAbsScalar(dot);
266:         }
267:         if (dot > dot_max) dot_max = dot;
268: #ifdef DEBUG_CkEigenSolutions
269:         if (dot > tols[1] ) {
270:           VecNorm(evec[i],NORM_INFINITY,&norm);
271:           PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %G, norm: %G\n",i,j,dot,norm);
272:         }
273: #endif
274:       } /* for (j=i; j<nev_loc; j++) */
275:     }
276:     PetscPrintf(PETSC_COMM_SELF,"    max|(x_j*B*x_i) - delta_ji|: %G\n",dot_max);

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