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 diamond_xxs_A -fB diamond_xxs_B \n\n";
11: #include petscmat.h
12: #include petscblaslapack.h
13: #include src/mat/impls/sbaij/seq/sbaij.h
19: PetscInt main(PetscInt argc,char **args)
20: {
21: Mat A,B,A_dense,B_dense,mats[2],A_sp;
22: Vec *evecs;
23: PetscViewer fd; /* viewer */
24: char file[3][PETSC_MAX_PATH_LEN]; /* input file name */
25: PetscTruth flg,flgA=PETSC_FALSE,flgB=PETSC_FALSE,TestSYGVX=PETSC_TRUE;
27: PetscTruth preload=PETSC_TRUE,isSymmetric;
28: PetscScalar sigma,one=1.0,*arrayA,*arrayB,*evecs_array,*work,*evals;
29: PetscMPIInt size;
30: PetscInt m,n,i,j,nevs,il,iu,stages[2];
31: PetscReal vl,vu,abstol=1.e-8;
32: PetscBLASInt *iwork,*ifail,lone=1,lwork,lierr,bn;
33: PetscInt ievbd_loc[2],offset=0,cklvl=2;
34: PetscReal tols[2];
35: Mat_SeqSBAIJ *sbaij;
36: PetscScalar *aa;
37: PetscInt *ai,*aj;
38: PetscInt nzeros[2],nz;
39: PetscReal ratio;
40:
41: PetscInitialize(&argc,&args,(char *)0,help);
42: MPI_Comm_size(PETSC_COMM_WORLD,&size);
43: if (size != 1) SETERRQ(PETSC_ERR_SUP,"This is a uniprocessor example only!");
44: PetscLogStageRegister(&stages[0],"EigSolve");
45: PetscLogStageRegister(&stages[1],"EigCheck");
47: /* Determine files from which we read the two matrices */
48: PetscOptionsGetString(PETSC_NULL,"-f0",file[0],PETSC_MAX_PATH_LEN-1,&flg);
49: if (!flg) {
50: PetscOptionsGetString(PETSC_NULL,"-fA",file[0],PETSC_MAX_PATH_LEN-1,&flgA);
51: if (!flgA) SETERRQ(PETSC_ERR_USER,"Must indicate binary file with the -fA or -fB options");
52: PetscOptionsGetString(PETSC_NULL,"-fB",file[1],PETSC_MAX_PATH_LEN-1,&flgB);
53: preload = PETSC_FALSE;
54: } else {
55: PetscOptionsGetString(PETSC_NULL,"-fA",file[1],PETSC_MAX_PATH_LEN-1,&flgA);
56: if (!flgA) {preload = PETSC_FALSE;} /* don't bother with second system */
57: PetscOptionsGetString(PETSC_NULL,"-fB",file[2],PETSC_MAX_PATH_LEN-1,&flgB);
58: }
60: PreLoadBegin(preload,"Load system");
61: /* Load matrices */
62: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt],FILE_MODE_READ,&fd);
63: MatLoad(fd,MATSBAIJ,&A);
64: PetscViewerDestroy(fd);
65: MatGetSize(A,&m,&n);
66: if ((flgB && PreLoadIt) || (flgB && !preload)){
67: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt+1],FILE_MODE_READ,&fd);
68: MatLoad(fd,MATSBAIJ,&B);
69: PetscViewerDestroy(fd);
70: } else { /* create B=I */
71: MatCreate(PETSC_COMM_WORLD,&B);
72: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
73: MatSetType(B,MATSEQSBAIJ);
74: MatSetFromOptions(B);
75: for (i=0; i<m; i++) {
76: MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES);
77: }
78: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
80: }
81:
82: /* Add a shift to A */
83: PetscOptionsGetScalar(PETSC_NULL,"-mat_sigma",&sigma,&flg);
84: if(flg) {
85: MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN); /* A <- sigma*B + A */
86: }
88: /* Check whether A is symmetric */
89: PetscOptionsHasName(PETSC_NULL, "-check_symmetry", &flg);
90: if (flg) {
91: Mat Trans;
92: MatTranspose(A, &Trans);
93: MatEqual(A, Trans, &isSymmetric);
94: if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"A must be symmetric");
95: MatDestroy(Trans);
96: if (flgB && PreLoadIt){
97: MatTranspose(B, &Trans);
98: MatEqual(B, Trans, &isSymmetric);
99: if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"B must be symmetric");
100: MatDestroy(Trans);
101: }
102: }
104: /* View small entries of A */
105: PetscOptionsHasName(PETSC_NULL, "-Asp_view", &flg);
106: if (flg){
107: MatCreate(PETSC_COMM_SELF,&A_sp);
108: MatSetSizes(A_sp,PETSC_DECIDE,PETSC_DECIDE,m,n);
109: MatSetType(A_sp,MATSEQSBAIJ);
111: tols[0] = 1.e-6, tols[1] = 1.e-9;
112: sbaij = (Mat_SeqSBAIJ*)A->data;
113: ai = sbaij->i;
114: aj = sbaij->j;
115: aa = sbaij->a;
116: nzeros[0] = nzeros[1] = 0;
117: for (i=0; i<m; i++) {
118: nz = ai[i+1] - ai[i];
119: for (j=0; j<nz; j++){
120: if (PetscAbsScalar(*aa)<tols[0]) {
121: MatSetValues(A_sp,1,&i,1,aj,aa,INSERT_VALUES);
122: nzeros[0]++;
123: }
124: if (PetscAbsScalar(*aa)<tols[1]) nzeros[1]++;
125: aa++; aj++;
126: }
127: }
128: MatAssemblyBegin(A_sp,MAT_FINAL_ASSEMBLY);
129: MatAssemblyEnd(A_sp,MAT_FINAL_ASSEMBLY);
131: MatDestroy(A_sp);
133: ratio = (PetscReal)nzeros[0]/sbaij->nz;
134: PetscPrintf(PETSC_COMM_SELF," %d matrix entries < %e, ratio %G of %d nonzeros\n",nzeros[0],tols[0],ratio,sbaij->nz);
135: PetscPrintf(PETSC_COMM_SELF," %d matrix entries < %e\n",nzeros[1],tols[1]);
136: }
138: /* Convert aij matrix to MatSeqDense for LAPACK */
139: PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
140: if (!flg) {
141: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
142: }
143: PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
144: if (!flg) {MatConvert(B,MATSEQDENSE,MAT_INITIAL_MATRIX,&B_dense);}
146: /* Solve eigenvalue problem: A*x = lambda*B*x */
147: /*============================================*/
148: lwork = 8*n;
149: bn = (PetscBLASInt)n;
150: PetscMalloc(n*sizeof(PetscScalar),&evals);
151: PetscMalloc(lwork*sizeof(PetscScalar),&work);
152: MatGetArray(A_dense,&arrayA);
153: MatGetArray(B_dense,&arrayB);
155: if (!TestSYGVX){ /* test sygv() */
156: evecs_array = arrayA;
157: LAPACKsygv_(&lone,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,&lierr);
158: nevs = m;
159: il=1;
160: } else { /* test sygvx() */
161: il = 1; iu=(PetscBLASInt)(.6*m); /* request 1 to 60%m evalues */
162: PetscMalloc((m*n+1)*sizeof(PetscScalar),&evecs_array);
163: PetscMalloc((6*n+1)*sizeof(PetscBLASInt),&iwork);
164: ifail = iwork + 5*n;
165: if(PreLoadIt){PetscLogStagePush(stages[0]);}
166: /* in the case "I", vl and vu are not referenced */
167: 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);
168: if(PreLoadIt){PetscLogStagePop();}
169: PetscFree(iwork);
170: }
171: MatRestoreArray(A,&arrayA);
172: MatRestoreArray(B,&arrayB);
174: if (nevs <= 0 ) SETERRQ1(PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
175: /* View evals */
176: PetscOptionsHasName(PETSC_NULL, "-eig_view", &flg);
177: if (flg){
178: printf(" %d evals: \n",nevs);
179: for (i=0; i<nevs; i++) printf("%d %G\n",i+il,evals[i]);
180: }
182: /* Check residuals and orthogonality */
183: if(PreLoadIt){
184: mats[0] = A; mats[1] = B;
185: one = (PetscInt)one;
186: PetscMalloc((nevs+1)*sizeof(Vec),&evecs);
187: for (i=0; i<nevs; i++){
188: VecCreate(PETSC_COMM_SELF,&evecs[i]);
189: VecSetSizes(evecs[i],PETSC_DECIDE,n);
190: VecSetFromOptions(evecs[i]);
191: VecPlaceArray(evecs[i],evecs_array+i*n);
192: }
193:
194: ievbd_loc[0] = 0; ievbd_loc[1] = nevs-1;
195: tols[0] = 1.e-8; tols[1] = 1.e-8;
196: PetscLogStagePush(stages[1]);
197: CkEigenSolutions(&cklvl,mats,evals,evecs,ievbd_loc,&offset,tols);
198: PetscLogStagePop();
199: for (i=0; i<nevs; i++){ VecDestroy(evecs[i]);}
200: PetscFree(evecs);
201: }
202:
203: /* Free work space. */
204: if (TestSYGVX){PetscFree(evecs_array);}
205:
206: PetscFree(evals);
207: PetscFree(work);
209: MatDestroy(A_dense);
210: MatDestroy(B_dense);
211: MatDestroy(B);
212: MatDestroy(A);
214: PreLoadEnd();
215: PetscFinalize();
216: return 0;
217: }
218: /*------------------------------------------------
219: Check the accuracy of the eigen solution
220: ----------------------------------------------- */
221: /*
222: input:
223: cklvl - check level:
224: 1: check residual
225: 2: 1 and check B-orthogonality locally
226: fA, fB - matrix pencil
227: eval, evec - eigenvalues and eigenvectors stored in this process
228: ievbd_loc - local eigenvalue bounds, see eigc()
229: offset - see eigc()
230: tols[0] - reporting tol_res: || A evec[i] - eval[i] B evec[i]||
231: tols[1] - reporting tol_orth: evec[i] B evec[j] - delta_ij
232: */
233: #undef DEBUG_CkEigenSolutions
236: PetscErrorCode CkEigenSolutions(PetscInt *fcklvl,Mat *mats,
237: PetscReal *eval,Vec *evec,PetscInt *ievbd_loc,PetscInt *offset,
238: PetscReal *tols)
239: {
240: PetscInt ierr,cklvl=*fcklvl,nev_loc,i,j;
241: Mat A=mats[0], B=mats[1];
242: Vec vt1,vt2; /* tmp vectors */
243: PetscReal norm,tmp,dot,norm_max,dot_max;
246: nev_loc = ievbd_loc[1] - ievbd_loc[0];
247: if (nev_loc == 0) return(0);
249: nev_loc += (*offset);
250: VecDuplicate(evec[*offset],&vt1);
251: VecDuplicate(evec[*offset],&vt2);
253: switch (cklvl){
254: case 2:
255: dot_max = 0.0;
256: for (i = *offset; i<nev_loc; i++){
257: MatMult(B, evec[i], vt1);
258: for (j=i; j<nev_loc; j++){
259: VecDot(evec[j],vt1,&dot);
260: if (j == i){
261: dot = PetscAbsScalar(dot - 1.0);
262: } else {
263: dot = PetscAbsScalar(dot);
264: }
265: if (dot > dot_max) dot_max = dot;
266: #ifdef DEBUG_CkEigenSolutions
267: if (dot > tols[1] ) {
268: VecNorm(evec[i],NORM_INFINITY,&norm);
269: PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %G, norm: %G\n",i,j,dot,norm);
270: }
271: #endif
272: } /* for (j=i; j<nev_loc; j++) */
273: }
274: PetscPrintf(PETSC_COMM_SELF," max|(x_j*B*x_i) - delta_ji|: %G\n",dot_max);
276: case 1:
277: norm_max = 0.0;
278: for (i = *offset; i< nev_loc; i++){
279: MatMult(A, evec[i], vt1);
280: MatMult(B, evec[i], vt2);
281: tmp = -eval[i];
282: VecAXPY(vt1,tmp,vt2);
283: VecNorm(vt1, NORM_INFINITY, &norm);
284: norm = PetscAbsScalar(norm);
285: if (norm > norm_max) norm_max = norm;
286: #ifdef DEBUG_CkEigenSolutions
287: /* sniff, and bark if necessary */
288: if (norm > tols[0]){
289: printf( " residual violation: %d, resi: %g\n",i, norm);
290: }
291: #endif
292: }
293:
294: PetscPrintf(PETSC_COMM_SELF," max_resi: %G\n", norm_max);
295:
296: break;
297: default:
298: PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl);
299: }
300: VecDestroy(vt2);
301: VecDestroy(vt1);
302: return(0);
303: }