Actual source code: ex99.c
petsc-3.6.1 2015-08-06
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: int main(int 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;
34: PetscBLASInt il,iu,nevs,nn;
35: PetscLogStage stages[2];
36: PetscReal vl,vu,abstol=1.e-8;
37: PetscBLASInt *iwork,*ifail,lone=1,lwork,lierr,bn;
38: PetscInt ievbd_loc[2],offset=0,cklvl=2;
39: PetscReal tols[2];
40: Mat_SeqSBAIJ *sbaij;
41: PetscScalar *aa;
42: PetscInt *ai,*aj;
43: PetscInt nzeros[2],nz;
44: PetscReal ratio;
46: PetscInitialize(&argc,&args,(char*)0,help);
47: MPI_Comm_size(PETSC_COMM_WORLD,&size);
48: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
49: PetscLogStageRegister("EigSolve",&stages[0]);
50: PetscLogStageRegister("EigCheck",&stages[1]);
52: /* Determine files from which we read the two matrices */
53: PetscOptionsGetString(NULL,"-f0",file[0],PETSC_MAX_PATH_LEN,&flg);
54: if (!flg) {
55: PetscOptionsGetString(NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,&flgA);
56: if (!flgA) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -fA or -fB options");
57: PetscOptionsGetString(NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flgB);
58: preload = PETSC_FALSE;
59: } else {
60: PetscOptionsGetString(NULL,"-fA",file[1],PETSC_MAX_PATH_LEN,&flgA);
61: if (!flgA) preload = PETSC_FALSE; /* don't bother with second system */
62: PetscOptionsGetString(NULL,"-fB",file[2],PETSC_MAX_PATH_LEN,&flgB);
63: }
65: PetscPreLoadBegin(preload,"Load system");
66: /* Load matrices */
67: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PetscPreLoadIt],FILE_MODE_READ,&fd);
68: MatCreate(PETSC_COMM_WORLD,&A);
69: MatSetType(A,MATSBAIJ);
70: MatLoad(A,fd);
71: PetscViewerDestroy(&fd);
72: MatGetSize(A,&m,&n);
73: if ((flgB && PetscPreLoadIt) || (flgB && !preload)) {
74: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PetscPreLoadIt+1],FILE_MODE_READ,&fd);
75: MatCreate(PETSC_COMM_WORLD,&B);
76: MatSetType(B,MATSBAIJ);
77: MatLoad(B,fd);
78: PetscViewerDestroy(&fd);
79: } else { /* create B=I */
80: MatCreate(PETSC_COMM_WORLD,&B);
81: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);
82: MatSetType(B,MATSEQSBAIJ);
83: MatSetFromOptions(B);
84: MatSetUp(B);
85: for (i=0; i<m; i++) {
86: MatSetValues(B,1,&i,1,&i,&one,INSERT_VALUES);
87: }
88: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
90: }
92: /* Add a shift to A */
93: PetscOptionsGetScalar(NULL,"-mat_sigma",&sigma,&flg);
94: if (flg) {
95: MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN); /* A <- sigma*B + A */
96: }
98: /* Check whether A is symmetric */
99: PetscOptionsHasName(NULL, "-check_symmetry", &flg);
100: if (flg) {
101: Mat Trans;
102: MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
103: MatEqual(A, Trans, &isSymmetric);
104: if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
105: MatDestroy(&Trans);
106: if (flgB && PetscPreLoadIt) {
107: MatTranspose(B,MAT_INITIAL_MATRIX, &Trans);
108: MatEqual(B, Trans, &isSymmetric);
109: if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"B must be symmetric");
110: MatDestroy(&Trans);
111: }
112: }
114: /* View small entries of A */
115: PetscOptionsHasName(NULL, "-Asp_view", &flg);
116: if (flg) {
117: MatCreate(PETSC_COMM_SELF,&A_sp);
118: MatSetSizes(A_sp,PETSC_DECIDE,PETSC_DECIDE,m,n);
119: MatSetType(A_sp,MATSEQSBAIJ);
121: tols[0] = 1.e-6, tols[1] = 1.e-9;
122: sbaij = (Mat_SeqSBAIJ*)A->data;
123: ai = sbaij->i;
124: aj = sbaij->j;
125: aa = sbaij->a;
126: nzeros[0] = nzeros[1] = 0;
127: for (i=0; i<m; i++) {
128: nz = ai[i+1] - ai[i];
129: for (j=0; j<nz; j++) {
130: if (PetscAbsScalar(*aa)<tols[0]) {
131: MatSetValues(A_sp,1,&i,1,aj,aa,INSERT_VALUES);
132: nzeros[0]++;
133: }
134: if (PetscAbsScalar(*aa)<tols[1]) nzeros[1]++;
135: aa++; aj++;
136: }
137: }
138: MatAssemblyBegin(A_sp,MAT_FINAL_ASSEMBLY);
139: MatAssemblyEnd(A_sp,MAT_FINAL_ASSEMBLY);
141: MatDestroy(&A_sp);
143: ratio = (PetscReal)nzeros[0]/sbaij->nz;
144: PetscPrintf(PETSC_COMM_SELF," %D matrix entries < %g, ratio %g of %d nonzeros\n",nzeros[0],(double)tols[0],(double)ratio,sbaij->nz);
145: PetscPrintf(PETSC_COMM_SELF," %D matrix entries < %g\n",nzeros[1],(double)tols[1]);
146: }
148: /* Convert aij matrix to MATSEQDENSE for LAPACK */
149: PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg);
150: if (!flg) {
151: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
152: }
153: PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
154: if (!flg) {MatConvert(B,MATSEQDENSE,MAT_INITIAL_MATRIX,&B_dense);}
156: /* Solve eigenvalue problem: A*x = lambda*B*x */
157: /*============================================*/
158: PetscBLASIntCast(8*n,&lwork);
159: PetscBLASIntCast(n,&bn);
160: PetscMalloc1(n,&evals);
161: PetscMalloc1(lwork,&work);
162: MatDenseGetArray(A_dense,&arrayA);
163: MatDenseGetArray(B_dense,&arrayB);
165: if (!TestSYGVX) { /* test sygv() */
166: evecs_array = arrayA;
167: LAPACKsygv_(&lone,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,&lierr);
168: nevs = m;
169: il =1;
170: } else { /* test sygvx() */
171: il = 1;
172: PetscBLASIntCast((PetscInt).6*m,&iu);
173: PetscBLASIntCast(n,&nn);
174: PetscMalloc1(m*n+1,&evecs_array);
175: PetscMalloc1(6*n+1,&iwork);
176: ifail = iwork + 5*n;
177: if (PetscPreLoadIt) {PetscLogStagePush(stages[0]);}
178: /* in the case "I", vl and vu are not referenced */
179: LAPACKsygvx_(&lone,"V","I","U",&bn,arrayA,&bn,arrayB,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&nn,work,&lwork,iwork,ifail,&lierr);
180: if (PetscPreLoadIt) PetscLogStagePop();
181: PetscFree(iwork);
182: }
183: MatDenseRestoreArray(A_dense,&arrayA);
184: MatDenseRestoreArray(B_dense,&arrayB);
186: if (nevs <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
187: /* View evals */
188: PetscOptionsHasName(NULL, "-eig_view", &flg);
189: if (flg) {
190: PetscPrintf(PETSC_COMM_SELF," %D evals: \n",nevs);
191: for (i=0; i<nevs; i++) {
192: PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i+il,(double)evals[i]);
193: }
194: }
196: /* Check residuals and orthogonality */
197: if (PetscPreLoadIt) {
198: mats[0] = A; mats[1] = B;
199: one = (PetscInt)one;
200: PetscMalloc1(nevs+1,&evecs);
201: for (i=0; i<nevs; i++) {
202: VecCreate(PETSC_COMM_SELF,&evecs[i]);
203: VecSetSizes(evecs[i],PETSC_DECIDE,n);
204: VecSetFromOptions(evecs[i]);
205: VecPlaceArray(evecs[i],evecs_array+i*n);
206: }
208: ievbd_loc[0] = 0; ievbd_loc[1] = nevs-1;
209: tols[0] = 1.e-8; tols[1] = 1.e-8;
211: PetscLogStagePush(stages[1]);
212: CkEigenSolutions(&cklvl,mats,evals,evecs,ievbd_loc,&offset,tols);
213: PetscLogStagePop();
214: for (i=0; i<nevs; i++) { VecDestroy(&evecs[i]);}
215: PetscFree(evecs);
216: }
218: /* Free work space. */
219: if (TestSYGVX) {PetscFree(evecs_array);}
221: PetscFree(evals);
222: PetscFree(work);
224: MatDestroy(&A_dense);
225: MatDestroy(&B_dense);
226: MatDestroy(&B);
227: MatDestroy(&A);
229: PetscPreLoadEnd();
230: PetscFinalize();
231: return 0;
232: }
233: /*------------------------------------------------
234: Check the accuracy of the eigen solution
235: ----------------------------------------------- */
236: /*
237: input:
238: cklvl - check level:
239: 1: check residual
240: 2: 1 and check B-orthogonality locally
241: mats - matrix pencil
242: eval, evec - eigenvalues and eigenvectors stored in this process
243: ievbd_loc - local eigenvalue bounds, see eigc()
244: offset - see eigc()
245: tols[0] - reporting tol_res: || A evec[i] - eval[i] B evec[i]||
246: tols[1] - reporting tol_orth: evec[i] B evec[j] - delta_ij
247: */
248: #undef DEBUG_CkEigenSolutions
251: PetscErrorCode CkEigenSolutions(PetscInt *fcklvl,Mat *mats,PetscReal *eval,Vec *evec,PetscInt *ievbd_loc,PetscInt *offset,PetscReal *tols)
252: {
253: PetscInt ierr,cklvl=*fcklvl,nev_loc,i,j;
254: Mat A=mats[0], B=mats[1];
255: Vec vt1,vt2; /* tmp vectors */
256: PetscReal norm,tmp,dot,norm_max,dot_max;
259: nev_loc = ievbd_loc[1] - ievbd_loc[0];
260: if (nev_loc == 0) return(0);
262: nev_loc += (*offset);
263: VecDuplicate(evec[*offset],&vt1);
264: VecDuplicate(evec[*offset],&vt2);
266: switch (cklvl) {
267: case 2:
268: dot_max = 0.0;
269: for (i = *offset; i<nev_loc; i++) {
270: MatMult(B, evec[i], vt1);
271: for (j=i; j<nev_loc; j++) {
272: VecDot(evec[j],vt1,&dot);
273: if (j == i) {
274: dot = PetscAbsScalar(dot - 1.0);
275: } else {
276: dot = PetscAbsScalar(dot);
277: }
278: if (dot > dot_max) dot_max = dot;
279: #if defined(DEBUG_CkEigenSolutions)
280: if (dot > tols[1]) {
281: VecNorm(evec[i],NORM_INFINITY,&norm);
282: PetscPrintf(PETSC_COMM_SELF,"|delta(%D,%D)|: %g, norm: %g\n",i,j,(double)ndot,(double)nnorm);
283: }
284: #endif
285: } /* for (j=i; j<nev_loc; j++) */
286: }
287: PetscPrintf(PETSC_COMM_SELF," max|(x_j*B*x_i) - delta_ji|: %g\n",(double)dot_max);
289: case 1:
290: norm_max = 0.0;
291: for (i = *offset; i< nev_loc; i++) {
292: MatMult(A, evec[i], vt1);
293: MatMult(B, evec[i], vt2);
294: tmp = -eval[i];
295: VecAXPY(vt1,tmp,vt2);
296: VecNorm(vt1, NORM_INFINITY, &norm);
297: norm = PetscAbsScalar(norm);
298: if (norm > norm_max) norm_max = norm;
299: #if defined(DEBUG_CkEigenSolutions)
300: /* sniff, and bark if necessary */
301: if (norm > tols[0]) {
302: PetscPrintf(PETSC_COMM_SELF," residual violation: %D, resi: %g\n",i, (double)nnorm);
303: }
304: #endif
305: }
307: PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max);
309: break;
310: default:
311: PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%D is not supported \n",cklvl);
312: }
313: VecDestroy(&vt2);
314: VecDestroy(&vt1);
315: return(0);
316: }