Actual source code: ex116.c
1: static char help[] = "Test LAPACK routine DSYEV() or DSYEVX(). \n\
2: Reads PETSc matrix A \n\
3: then computes selected eigenvalues, and optionally, eigenvectors of \n\
4: a real generalized symmetric-definite eigenproblem \n\
5: A*x = lambda*x \n\
6: Input parameters include\n\
7: -f <input_file> : file to load\n\
8: e.g. ./ex116 -f /home/petsc/datafiles/matrices/small \n\n";
10: #include petscmat.h
11: #include petscblaslapack.h
17: PetscInt main(PetscInt argc,char **args)
18: {
19: Mat A,A_dense;
20: Vec *evecs;
21: PetscViewer fd; /* viewer */
22: char file[1][PETSC_MAX_PATH_LEN]; /* input file name */
23: PetscTruth flg,flgA=PETSC_FALSE,flgB=PETSC_FALSE,TestSYEVX=PETSC_TRUE;
25: PetscTruth isSymmetric;
26: PetscScalar sigma,*arrayA,*arrayB,*evecs_array,*work,*evals;
27: PetscMPIInt size;
28: PetscInt m,n,i,j,nevs,il,iu,cklvl=2;
29: PetscReal vl,vu,abstol=1.e-8;
30: PetscBLASInt *iwork,*ifail,lwork,lierr,bn;
31: PetscReal tols[2];
32: PetscInt nzeros[2],nz;
33: PetscReal ratio;
34:
35: PetscInitialize(&argc,&args,(char *)0,help);
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_syev", &flg);
40: if (flg){
41: TestSYEVX = PETSC_FALSE;
42: }
44: /* Determine files from which we read the two matrices */
45: PetscOptionsGetString(PETSC_NULL,"-f",file[0],PETSC_MAX_PATH_LEN-1,&flg);
47: /* Load matrix A */
48: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd);
49: MatLoad(fd,MATSEQAIJ,&A);
50: PetscViewerDestroy(fd);
51: MatGetSize(A,&m,&n);
53: /* Check whether A is symmetric */
54: PetscOptionsHasName(PETSC_NULL, "-check_symmetry", &flg);
55: if (flg) {
56: Mat Trans;
57: MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
58: MatEqual(A, Trans, &isSymmetric);
59: if (!isSymmetric) SETERRQ(PETSC_ERR_USER,"A must be symmetric");
60: MatDestroy(Trans);
61: }
63: /* Convert aij matrix to MatSeqDense for LAPACK */
64: PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
65: if (!flg) {
66: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
67: }
69: /* Solve eigenvalue problem: A*x = lambda*B*x */
70: /*============================================*/
71: lwork = PetscBLASIntCast(8*n);
72: bn = PetscBLASIntCast(n);
73: PetscMalloc(n*sizeof(PetscScalar),&evals);
74: PetscMalloc(lwork*sizeof(PetscScalar),&work);
75: MatGetArray(A_dense,&arrayA);
77: if (!TestSYEVX){ /* test syev() */
78: printf(" LAPACKsyev: compute all %d eigensolutions...\n",m);
79: LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,&lierr);
80: evecs_array = arrayA;
81: nevs = PetscBLASIntCast(m);
82: il=1; iu=PetscBLASIntCast(m);
83: } else { /* test syevx() */
84: il = 1; iu=PetscBLASIntCast((0.2*m)); /* request 1 to 20%m evalues */
85: printf(" LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu);
86: PetscMalloc((m*n+1)*sizeof(PetscScalar),&evecs_array);
87: PetscMalloc((6*n+1)*sizeof(PetscBLASInt),&iwork);
88: ifail = iwork + 5*n;
89:
90: /* in the case "I", vl and vu are not referenced */
91: vl = 0.0; vu = 8.0;
92: LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&n,work,&lwork,iwork,ifail,&lierr);
93: PetscFree(iwork);
94: }
95: MatRestoreArray(A,&arrayA);
96: if (nevs <= 0 ) SETERRQ1(PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
98: /* View evals */
99: PetscOptionsHasName(PETSC_NULL, "-eig_view", &flg);
100: if (flg){
101: printf(" %d evals: \n",nevs);
102: for (i=0; i<nevs; i++) printf("%d %G\n",i+il,evals[i]);
103: }
105: /* Check residuals and orthogonality */
106: PetscMalloc((nevs+1)*sizeof(Vec),&evecs);
107: for (i=0; i<nevs; i++){
108: VecCreate(PETSC_COMM_SELF,&evecs[i]);
109: VecSetSizes(evecs[i],PETSC_DECIDE,n);
110: VecSetFromOptions(evecs[i]);
111: VecPlaceArray(evecs[i],evecs_array+i*n);
112: }
113:
114: tols[0] = 1.e-8; tols[1] = 1.e-8;
115: CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols);
116: for (i=0; i<nevs; i++){ VecDestroy(evecs[i]);}
117: PetscFree(evecs);
118:
119: /* Free work space. */
120: if (TestSYEVX){PetscFree(evecs_array);}
121:
122: PetscFree(evals);
123: PetscFree(work);
125: MatDestroy(A_dense);
126: MatDestroy(A);
127: PetscFinalize();
128: return 0;
129: }
130: /*------------------------------------------------
131: Check the accuracy of the eigen solution
132: ----------------------------------------------- */
133: /*
134: input:
135: cklvl - check level:
136: 1: check residual
137: 2: 1 and check B-orthogonality locally
138: A - matrix
139: il,iu - lower and upper index bound of eigenvalues
140: eval, evec - eigenvalues and eigenvectors stored in this process
141: tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
142: tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
143: */
144: #undef DEBUG_CkEigenSolutions
147: PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
148: {
149: PetscInt ierr,i,j,nev;
150: Vec vt1,vt2; /* tmp vectors */
151: PetscReal norm,tmp,dot,norm_max,dot_max;
154: nev = iu - il;
155: if (nev <= 0) return(0);
157: //VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);
158: VecDuplicate(evec[0],&vt1);
159: VecDuplicate(evec[0],&vt2);
161: switch (cklvl){
162: case 2:
163: dot_max = 0.0;
164: for (i = il; i<iu; i++){
165: //printf("ck %d-th\n",i);
166: VecCopy(evec[i], vt1);
167: for (j=il; j<iu; j++){
168: VecDot(evec[j],vt1,&dot);
169: if (j == i){
170: dot = PetscAbsScalar(dot - 1.0);
171: } else {
172: dot = PetscAbsScalar(dot);
173: }
174: if (dot > dot_max) dot_max = dot;
175: #ifdef DEBUG_CkEigenSolutions
176: if (dot > tols[1] ) {
177: VecNorm(evec[i],NORM_INFINITY,&norm);
178: PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %G, norm: %G\n",i,j,dot,norm);
179: }
180: #endif
181: }
182: }
183: PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %G\n",dot_max);
185: case 1:
186: norm_max = 0.0;
187: for (i = il; i< iu; i++){
188: MatMult(A, evec[i], vt1);
189: VecCopy(evec[i], vt2);
190: tmp = -eval[i];
191: VecAXPY(vt1,tmp,vt2);
192: VecNorm(vt1, NORM_INFINITY, &norm);
193: norm = PetscAbsScalar(norm);
194: if (norm > norm_max) norm_max = norm;
195: #ifdef DEBUG_CkEigenSolutions
196: /* sniff, and bark if necessary */
197: if (norm > tols[0]){
198: printf( " residual violation: %d, resi: %g\n",i, norm);
199: }
200: #endif
201: }
202: PetscPrintf(PETSC_COMM_SELF," max_resi: %G\n", norm_max);
203: break;
204: default:
205: PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl);
206: }
207: VecDestroy(vt2);
208: VecDestroy(vt1);
209: return(0);
210: }