Actual source code: ex116.c
petsc-3.13.6 2020-09-29
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 $DATAFILESPATH/matrices/small \n\n";
10: #include <petscmat.h>
11: #include <petscblaslapack.h>
13: extern PetscErrorCode CkEigenSolutions(PetscInt,Mat,PetscInt,PetscInt,PetscReal*,Vec*,PetscReal*);
15: int main(int argc,char **args)
16: {
17: Mat A,A_dense;
18: Vec *evecs;
19: PetscViewer fd; /* viewer */
20: char file[1][PETSC_MAX_PATH_LEN]; /* input file name */
21: PetscBool flg,TestSYEVX=PETSC_TRUE;
23: PetscBool isSymmetric;
24: PetscScalar *arrayA,*evecs_array,*work,*evals;
25: PetscMPIInt size;
26: PetscInt m,n,i,cklvl=2;
27: PetscBLASInt nevs,il,iu,in;
28: PetscReal vl,vu,abstol=1.e-8;
29: PetscBLASInt *iwork,*ifail,lwork,lierr,bn;
30: PetscReal tols[2];
32: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
33: MPI_Comm_size(PETSC_COMM_WORLD,&size);
34: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
36: PetscOptionsHasName(NULL,NULL, "-test_syev", &flg);
37: if (flg) {
38: TestSYEVX = PETSC_FALSE;
39: }
41: /* Determine files from which we read the two matrices */
42: PetscOptionsGetString(NULL,NULL,"-f",file[0],PETSC_MAX_PATH_LEN,&flg);
44: /* Load matrix A */
45: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd);
46: MatCreate(PETSC_COMM_WORLD,&A);
47: MatSetType(A,MATSEQAIJ);
48: MatLoad(A,fd);
49: PetscViewerDestroy(&fd);
50: MatGetSize(A,&m,&n);
52: /* Check whether A is symmetric */
53: PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg);
54: if (flg) {
55: Mat Trans;
56: MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);
57: MatEqual(A, Trans, &isSymmetric);
58: if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
59: MatDestroy(&Trans);
60: }
62: /* Solve eigenvalue problem: A_dense*x = lambda*B*x */
63: /*==================================================*/
64: /* Convert aij matrix to MatSeqDense for LAPACK */
65: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
67: PetscBLASIntCast(8*n,&lwork);
68: PetscBLASIntCast(n,&bn);
69: PetscMalloc1(n,&evals);
70: PetscMalloc1(lwork,&work);
71: MatDenseGetArray(A_dense,&arrayA);
73: if (!TestSYEVX) { /* test syev() */
74: PetscPrintf(PETSC_COMM_SELF," LAPACKsyev: compute all %D eigensolutions...\n",m);
75: LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,&lierr);
76: evecs_array = arrayA;
77: PetscBLASIntCast(m,&nevs);
78: il = 1;
79: PetscBLASIntCast(m,&iu);
80: } else { /* test syevx() */
81: il = 1;
82: PetscBLASIntCast(0.2*m,&iu);
83: PetscBLASIntCast(n,&in);
84: PetscPrintf(PETSC_COMM_SELF," LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu);
85: PetscMalloc1(m*n+1,&evecs_array);
86: PetscMalloc1(6*n+1,&iwork);
87: ifail = iwork + 5*n;
89: /* in the case "I", vl and vu are not referenced */
90: vl = 0.0; vu = 8.0;
91: LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&in,work,&lwork,iwork,ifail,&lierr);
92: PetscFree(iwork);
93: }
94: MatDenseRestoreArray(A_dense,&arrayA);
95: if (nevs <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%D, no eigensolution has found", nevs);
97: /* View eigenvalues */
98: PetscOptionsHasName(NULL,NULL, "-eig_view", &flg);
99: if (flg) {
100: PetscPrintf(PETSC_COMM_SELF," %D evals: \n",nevs);
101: for (i=0; i<nevs; i++) {PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i+il,(double)evals[i]);}
102: }
104: /* Check residuals and orthogonality */
105: PetscMalloc1(nevs+1,&evecs);
106: for (i=0; i<nevs; i++) {
107: VecCreate(PETSC_COMM_SELF,&evecs[i]);
108: VecSetSizes(evecs[i],PETSC_DECIDE,n);
109: VecSetFromOptions(evecs[i]);
110: VecPlaceArray(evecs[i],evecs_array+i*n);
111: }
113: tols[0] = tols[1] = PETSC_SQRT_MACHINE_EPSILON;
114: CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols);
116: /* Free work space. */
117: for (i=0; i<nevs; i++) { VecDestroy(&evecs[i]);}
118: PetscFree(evecs);
119: MatDestroy(&A_dense);
120: PetscFree(work);
121: if (TestSYEVX) {PetscFree(evecs_array);}
123: /* Compute SVD: A_dense = U*SIGMA*transpose(V),
124: JOBU=JOBV='S': the first min(m,n) columns of U and V are returned in the arrayU and arrayV; */
125: /*==============================================================================================*/
126: {
127: /* Convert aij matrix to MatSeqDense for LAPACK */
128: PetscScalar *arrayU,*arrayVT,*arrayErr,alpha=1.0,beta=-1.0;
129: Mat Err;
130: PetscBLASInt minMN,maxMN,im,in;
131: PetscInt j;
132: PetscReal norm;
134: MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);
136: minMN = PetscMin(m,n);
137: maxMN = PetscMax(m,n);
138: lwork = 5*minMN + maxMN;
139: PetscMalloc4(m*minMN,&arrayU,m*minMN,&arrayVT,m*minMN,&arrayErr,lwork,&work);
141: /* Create matrix Err for checking error */
142: MatCreate(PETSC_COMM_WORLD,&Err);
143: MatSetSizes(Err,PETSC_DECIDE,PETSC_DECIDE,m,minMN);
144: MatSetType(Err,MATSEQDENSE);
145: MatSeqDenseSetPreallocation(Err,(PetscScalar*)arrayErr);
147: /* Save A to arrayErr for checking accuracy later. arrayA will be destroyed by LAPACKgesvd_() */
148: MatDenseGetArray(A_dense,&arrayA);
149: PetscArraycpy(arrayErr,arrayA,m*minMN);
151: PetscBLASIntCast(m,&im);
152: PetscBLASIntCast(n,&in);
153: /* Compute A = U*SIGMA*VT */
154: LAPACKgesvd_("S","S",&im,&in,arrayA,&im,evals,arrayU,&minMN,arrayVT,&minMN,work,&lwork,&lierr);
155: MatDenseRestoreArray(A_dense,&arrayA);
156: if (!lierr) {
157: PetscPrintf(PETSC_COMM_SELF," 1st 10 of %d singular values: \n",minMN);
158: for (i=0; i<10; i++) {PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i,(double)evals[i]);}
159: } else {
160: PetscPrintf(PETSC_COMM_SELF,"LAPACKgesvd_ fails!");
161: }
163: /* Check Err = (U*Sigma*V^T - A) using BLASgemm() */
164: /* U = U*Sigma */
165: for (j=0; j<minMN; j++) { /* U[:,j] = sigma[j]*U[:,j] */
166: for (i=0; i<m; i++) arrayU[j*m+i] *= evals[j];
167: }
168: /* Err = U*VT - A = alpha*U*VT + beta*Err */
169: BLASgemm_("N","N",&im,&minMN,&minMN,&alpha,arrayU,&im,arrayVT,&minMN,&beta,arrayErr,&im);
170: MatNorm(Err,NORM_FROBENIUS,&norm);
171: PetscPrintf(PETSC_COMM_SELF," || U*Sigma*VT - A || = %g\n",(double)norm);
173: PetscFree4(arrayU,arrayVT,arrayErr,work);
174: PetscFree(evals);
175: MatDestroy(&A_dense);
176: MatDestroy(&Err);
177: }
179: MatDestroy(&A);
180: PetscFinalize();
181: return ierr;
182: }
183: /*------------------------------------------------
184: Check the accuracy of the eigen solution
185: ----------------------------------------------- */
186: /*
187: input:
188: cklvl - check level:
189: 1: check residual
190: 2: 1 and check B-orthogonality locally
191: A - matrix
192: il,iu - lower and upper index bound of eigenvalues
193: eval, evec - eigenvalues and eigenvectors stored in this process
194: tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
195: tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
196: */
197: PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
198: {
199: PetscInt ierr,i,j,nev;
200: Vec vt1,vt2; /* tmp vectors */
201: PetscReal norm,tmp,dot,norm_max,dot_max;
204: nev = iu - il;
205: if (nev <= 0) return(0);
207: /*VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);*/
208: VecDuplicate(evec[0],&vt1);
209: VecDuplicate(evec[0],&vt2);
211: switch (cklvl) {
212: case 2:
213: dot_max = 0.0;
214: for (i = il; i<iu; i++) {
215: VecCopy(evec[i], vt1);
216: for (j=il; j<iu; j++) {
217: VecDot(evec[j],vt1,&dot);
218: if (j == i) {
219: dot = PetscAbsScalar(dot - 1);
220: } else {
221: dot = PetscAbsScalar(dot);
222: }
223: if (dot > dot_max) dot_max = dot;
224: if (dot > tols[1]) {
225: VecNorm(evec[i],NORM_INFINITY,&norm);
226: PetscPrintf(PETSC_COMM_SELF,"|delta(%D,%D)|: %g, norm: %g\n",i,j,(double)dot,(double)norm);
227: }
228: }
229: }
230: PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max);
232: case 1:
233: norm_max = 0.0;
234: for (i = il; i< iu; i++) {
235: MatMult(A, evec[i], vt1);
236: VecCopy(evec[i], vt2);
237: tmp = -eval[i];
238: VecAXPY(vt1,tmp,vt2);
239: VecNorm(vt1, NORM_INFINITY, &norm);
240: norm = PetscAbsScalar(norm);
241: if (norm > norm_max) norm_max = norm;
242: /* sniff, and bark if necessary */
243: if (norm > tols[0]) {
244: PetscPrintf(PETSC_COMM_SELF," residual violation: %D, resi: %g\n",i, norm);
245: }
246: }
247: PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max);
248: break;
249: default:
250: PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%D is not supported \n",cklvl);
251: }
252: VecDestroy(&vt2);
253: VecDestroy(&vt1);
254: return(0);
255: }
258: /*TEST
260: build:
261: requires: !complex
263: test:
264: requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
265: args: -f ${DATAFILESPATH}/matrices/small
266: output_file: output/ex116_1.out
268: test:
269: suffix: 2
270: requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
271: args: -f ${DATAFILESPATH}/matrices/small -test_syev -check_symmetry
273: TEST*/