Actual source code: ex99.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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: /* This example only works with real numbers */

 14: #include <petscmat.h>
 15: #include <../src/mat/impls/sbaij/seq/sbaij.h>
 16: #include <petscblaslapack.h>

 18: extern PetscErrorCode CkEigenSolutions(PetscInt*,Mat*,PetscReal*,Vec*,PetscInt*,PetscInt*,PetscReal*);

 22: PetscInt main(PetscInt argc,char **args)
 23: {
 24:   Mat            A,B,A_dense,B_dense,mats[2],A_sp;
 25:   Vec            *evecs;
 26:   PetscViewer    fd;                /* viewer */
 27:   char           file[3][PETSC_MAX_PATH_LEN];     /* input file name */
 28:   PetscBool      flg,flgA=PETSC_FALSE,flgB=PETSC_FALSE,TestSYGVX=PETSC_TRUE;
 30:   PetscBool      preload=PETSC_TRUE,isSymmetric;
 31:   PetscScalar    sigma,one=1.0,*arrayA,*arrayB,*evecs_array,*work,*evals;
 32:   PetscMPIInt    size;
 33:   PetscInt       m,n,i,j,nevs,il,iu;
 34:   PetscLogStage  stages[2];
 35:   PetscReal      vl,vu,abstol=1.e-8;
 36:   PetscBLASInt   *iwork,*ifail,lone=1,lwork,lierr,bn;
 37:   PetscInt       ievbd_loc[2],offset=0,cklvl=2;
 38:   PetscReal      tols[2];
 39:   Mat_SeqSBAIJ   *sbaij;
 40:   PetscScalar    *aa;
 41:   PetscInt       *ai,*aj;
 42:   PetscInt       nzeros[2],nz;
 43:   PetscReal      ratio;

 45:   PetscInitialize(&argc,&args,(char*)0,help);
 46:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 47:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
 48:   PetscLogStageRegister("EigSolve",&stages[0]);
 49:   PetscLogStageRegister("EigCheck",&stages[1]);

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

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

 91:   /* Add a shift to A */
 92:   PetscOptionsGetScalar(NULL,"-mat_sigma",&sigma,&flg);
 93:   if (flg) {
 94:     MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN);   /* A <- sigma*B + A */
 95:   }

 97:   /* Check whether A is symmetric */
 98:   PetscOptionsHasName(NULL, "-check_symmetry", &flg);
 99:   if (flg) {
100:     Mat Trans;
101:     MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
102:     MatEqual(A, Trans, &isSymmetric);
103:     if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
104:     MatDestroy(&Trans);
105:     if (flgB && PetscPreLoadIt) {
106:       MatTranspose(B,MAT_INITIAL_MATRIX, &Trans);
107:       MatEqual(B, Trans, &isSymmetric);
108:       if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"B must be symmetric");
109:       MatDestroy(&Trans);
110:     }
111:   }

113:   /* View small entries of A */
114:   PetscOptionsHasName(NULL, "-Asp_view", &flg);
115:   if (flg) {
116:     MatCreate(PETSC_COMM_SELF,&A_sp);
117:     MatSetSizes(A_sp,PETSC_DECIDE,PETSC_DECIDE,m,n);
118:     MatSetType(A_sp,MATSEQSBAIJ);

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

140:     MatDestroy(&A_sp);

142:     ratio = (PetscReal)nzeros[0]/sbaij->nz;
143:     PetscPrintf(PETSC_COMM_SELF," %D matrix entries < %g, ratio %g of %d nonzeros\n",nzeros[0],(double)tols[0],(double)ratio,sbaij->nz);
144:     PetscPrintf(PETSC_COMM_SELF," %D matrix entries < %g\n",nzeros[1],(double)tols[1]);
145:   }

147:   /* Convert aij matrix to MATSEQDENSE for LAPACK */
148:   PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
149:   if (!flg) {
150:     MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
151:   }
152:   PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
153:   if (!flg) {MatConvert(B,MATSEQDENSE,MAT_INITIAL_MATRIX,&B_dense);}

155:   /* Solve eigenvalue problem: A*x = lambda*B*x */
156:   /*============================================*/
157:   PetscBLASIntCast(8*n,&lwork);
158:   PetscBLASIntCast(n,&bn);
159:   PetscMalloc1(n,&evals);
160:   PetscMalloc1(lwork,&work);
161:   MatDenseGetArray(A_dense,&arrayA);
162:   MatDenseGetArray(B_dense,&arrayB);

164:   if (!TestSYGVX) {   /* test sygv()  */
165:     evecs_array = arrayA;
166:     LAPACKsygv_(&lone,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,&lierr);
167:     nevs = m;
168:     il   =1;
169:   } else {   /* test sygvx()  */
170:     il    = 1;
171:     PetscBLASIntCast(.6*m,&iu);
172:     PetscMalloc1((m*n+1),&evecs_array);
173:     PetscMalloc1((6*n+1),&iwork);
174:     ifail = iwork + 5*n;
175:     if (PetscPreLoadIt) {PetscLogStagePush(stages[0]);}
176:     /* in the case "I", vl and vu are not referenced */
177:     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);
178:     if (PetscPreLoadIt) PetscLogStagePop();
179:     PetscFree(iwork);
180:   }
181:   MatDenseRestoreArray(A_dense,&arrayA);
182:   MatDenseRestoreArray(B_dense,&arrayB);

184:   if (nevs <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
185:   /* View evals */
186:   PetscOptionsHasName(NULL, "-eig_view", &flg);
187:   if (flg) {
188:     printf(" %d evals: \n",nevs);
189:     for (i=0; i<nevs; i++) printf("%D  %g\n",i+il,(double)evals[i]);
190:   }

192:   /* Check residuals and orthogonality */
193:   if (PetscPreLoadIt) {
194:     mats[0] = A; mats[1] = B;
195:     one     = (PetscInt)one;
196:     PetscMalloc1((nevs+1),&evecs);
197:     for (i=0; i<nevs; i++) {
198:       VecCreate(PETSC_COMM_SELF,&evecs[i]);
199:       VecSetSizes(evecs[i],PETSC_DECIDE,n);
200:       VecSetFromOptions(evecs[i]);
201:       VecPlaceArray(evecs[i],evecs_array+i*n);
202:     }

204:     ievbd_loc[0] = 0; ievbd_loc[1] = nevs-1;
205:     tols[0]      = 1.e-8;  tols[1] = 1.e-8;

207:     PetscLogStagePush(stages[1]);
208:     CkEigenSolutions(&cklvl,mats,evals,evecs,ievbd_loc,&offset,tols);
209:     PetscLogStagePop();
210:     for (i=0; i<nevs; i++) { VecDestroy(&evecs[i]);}
211:     PetscFree(evecs);
212:   }

214:   /* Free work space. */
215:   if (TestSYGVX) {PetscFree(evecs_array);}

217:   PetscFree(evals);
218:   PetscFree(work);

220:   MatDestroy(&A_dense);
221:   MatDestroy(&B_dense);
222:   MatDestroy(&B);
223:   MatDestroy(&A);

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

255:   nev_loc = ievbd_loc[1] - ievbd_loc[0];
256:   if (nev_loc == 0) return(0);

258:   nev_loc += (*offset);
259:   VecDuplicate(evec[*offset],&vt1);
260:   VecDuplicate(evec[*offset],&vt2);

262:   switch (cklvl) {
263:   case 2:
264:     dot_max = 0.0;
265:     for (i = *offset; i<nev_loc; i++) {
266:       MatMult(B, evec[i], vt1);
267:       for (j=i; j<nev_loc; j++) {
268:         VecDot(evec[j],vt1,&dot);
269:         if (j == i) {
270:           dot = PetscAbsScalar(dot - 1.0);
271:         } else {
272:           dot = PetscAbsScalar(dot);
273:         }
274:         if (dot > dot_max) dot_max = dot;
275: #if defined(DEBUG_CkEigenSolutions)
276:         if (dot > tols[1]) {
277:           VecNorm(evec[i],NORM_INFINITY,&norm);
278:           PetscPrintf(PETSC_COMM_SELF,"|delta(%D,%D)|: %g, norm: %g\n",i,j,(double)ndot,(double)nnorm);
279:         }
280: #endif
281:       } /* for (j=i; j<nev_loc; j++) */
282:     }
283:     PetscPrintf(PETSC_COMM_SELF,"    max|(x_j*B*x_i) - delta_ji|: %g\n",(double)dot_max);

285:   case 1:
286:     norm_max = 0.0;
287:     for (i = *offset; i< nev_loc; i++) {
288:       MatMult(A, evec[i], vt1);
289:       MatMult(B, evec[i], vt2);
290:       tmp  = -eval[i];
291:       VecAXPY(vt1,tmp,vt2);
292:       VecNorm(vt1, NORM_INFINITY, &norm);
293:       norm = PetscAbsScalar(norm);
294:       if (norm > norm_max) norm_max = norm;
295: #if defined(DEBUG_CkEigenSolutions)
296:       /* sniff, and bark if necessary */
297:       if (norm > tols[0]) {
298:         printf("  residual violation: %D, resi: %g\n",i, (double)nnorm);
299:       }
300: #endif
301:     }

303:     PetscPrintf(PETSC_COMM_SELF,"    max_resi:                    %g\n", (double)norm_max);

305:     break;
306:   default:
307:     PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%D is not supported \n",cklvl);
308:   }
309:   VecDestroy(&vt2);
310:   VecDestroy(&vt1);
311:   return(0);
312: }