1c4762a1bSJed Brown static char help[] = "Test LAPACK routine DSYEV() or DSYEVX(). \n\
2c4762a1bSJed Brown Reads PETSc matrix A \n\
3c4762a1bSJed Brown then computes selected eigenvalues, and optionally, eigenvectors of \n\
4c4762a1bSJed Brown a real generalized symmetric-definite eigenproblem \n\
5c4762a1bSJed Brown A*x = lambda*x \n\
6c4762a1bSJed Brown Input parameters include\n\
7c4762a1bSJed Brown -f <input_file> : file to load\n\
8c4762a1bSJed Brown e.g. ./ex116 -f $DATAFILESPATH/matrices/small \n\n";
9c4762a1bSJed Brown
10c4762a1bSJed Brown #include <petscmat.h>
11c4762a1bSJed Brown #include <petscblaslapack.h>
12c4762a1bSJed Brown
13c4762a1bSJed Brown extern PetscErrorCode CkEigenSolutions(PetscInt, Mat, PetscInt, PetscInt, PetscReal *, Vec *, PetscReal *);
14c4762a1bSJed Brown
main(int argc,char ** args)15d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
16d71ae5a4SJacob Faibussowitsch {
17c4762a1bSJed Brown Mat A, A_dense;
18c4762a1bSJed Brown Vec *evecs;
19c4762a1bSJed Brown PetscViewer fd; /* viewer */
20c4762a1bSJed Brown char file[1][PETSC_MAX_PATH_LEN]; /* input file name */
21c4762a1bSJed Brown PetscBool flg, TestSYEVX = PETSC_TRUE;
22c4762a1bSJed Brown PetscBool isSymmetric;
23c4762a1bSJed Brown PetscScalar *arrayA, *evecs_array, *work, *evals;
24c4762a1bSJed Brown PetscMPIInt size;
25c4762a1bSJed Brown PetscInt m, n, i, cklvl = 2;
26c4762a1bSJed Brown PetscBLASInt nevs, il, iu, in;
27c4762a1bSJed Brown PetscReal vl, vu, abstol = 1.e-8;
28c4762a1bSJed Brown PetscBLASInt *iwork, *ifail, lwork, lierr, bn;
29c4762a1bSJed Brown PetscReal tols[2];
30c4762a1bSJed Brown
31327415f7SBarry Smith PetscFunctionBeginUser;
32*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
34be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
35c4762a1bSJed Brown
369566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_syev", &flg));
37ad540459SPierre Jolivet if (flg) TestSYEVX = PETSC_FALSE;
38c4762a1bSJed Brown
39c4762a1bSJed Brown /* Determine files from which we read the two matrices */
409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file[0], sizeof(file[0]), &flg));
41c4762a1bSJed Brown
42c4762a1bSJed Brown /* Load matrix A */
439566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &fd));
449566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
459566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ));
469566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd));
479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
489566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n));
49c4762a1bSJed Brown
50c4762a1bSJed Brown /* Check whether A is symmetric */
519566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetry", &flg));
52c4762a1bSJed Brown if (flg) {
53c4762a1bSJed Brown Mat Trans;
549566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Trans));
559566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Trans, &isSymmetric));
5628b400f6SJacob Faibussowitsch PetscCheck(isSymmetric, PETSC_COMM_SELF, PETSC_ERR_USER, "A must be symmetric");
579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Trans));
58c4762a1bSJed Brown }
59c4762a1bSJed Brown
60c4762a1bSJed Brown /* Solve eigenvalue problem: A_dense*x = lambda*B*x */
61c4762a1bSJed Brown /*==================================================*/
62c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */
639566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &A_dense));
64c4762a1bSJed Brown
659566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(8 * n, &lwork));
669566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &bn));
679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &evals));
689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lwork, &work));
699566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(A_dense, &arrayA));
70c4762a1bSJed Brown
71c4762a1bSJed Brown if (!TestSYEVX) { /* test syev() */
729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n", m));
73c4762a1bSJed Brown LAPACKsyev_("V", "U", &bn, arrayA, &bn, evals, work, &lwork, &lierr);
74c4762a1bSJed Brown evecs_array = arrayA;
759566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &nevs));
76c4762a1bSJed Brown il = 1;
779566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &iu));
78c4762a1bSJed Brown } else { /* test syevx() */
79c4762a1bSJed Brown il = 1;
809566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(0.2 * m, &iu));
819566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &in));
829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LAPACKsyevx: compute %" PetscBLASInt_FMT " to %" PetscBLASInt_FMT "-th eigensolutions...\n", il, iu));
839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * n + 1, &evecs_array));
849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * n + 1, &iwork));
85c4762a1bSJed Brown ifail = iwork + 5 * n;
86c4762a1bSJed Brown
87c4762a1bSJed Brown /* in the case "I", vl and vu are not referenced */
889371c9d4SSatish Balay vl = 0.0;
899371c9d4SSatish Balay vu = 8.0;
90c4762a1bSJed Brown LAPACKsyevx_("V", "I", "U", &bn, arrayA, &bn, &vl, &vu, &il, &iu, &abstol, &nevs, evals, evecs_array, &in, work, &lwork, iwork, ifail, &lierr);
919566063dSJacob Faibussowitsch PetscCall(PetscFree(iwork));
92c4762a1bSJed Brown }
939566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(A_dense, &arrayA));
94e00437b9SBarry Smith PetscCheck(nevs > 0, PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "nev=%" PetscBLASInt_FMT ", no eigensolution has found", nevs);
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* View eigenvalues */
979566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-eig_view", &flg));
98c4762a1bSJed Brown if (flg) {
999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscBLASInt_FMT " evals: \n", nevs));
1009566063dSJacob Faibussowitsch for (i = 0; i < nevs; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g\n", (PetscInt)(i + il), (double)evals[i]));
101c4762a1bSJed Brown }
102c4762a1bSJed Brown
103c4762a1bSJed Brown /* Check residuals and orthogonality */
1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevs + 1, &evecs));
105c4762a1bSJed Brown for (i = 0; i < nevs; i++) {
1069566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &evecs[i]));
1079566063dSJacob Faibussowitsch PetscCall(VecSetSizes(evecs[i], PETSC_DECIDE, n));
1089566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(evecs[i]));
1099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(evecs[i], evecs_array + i * n));
110c4762a1bSJed Brown }
111c4762a1bSJed Brown
112c4762a1bSJed Brown tols[0] = tols[1] = PETSC_SQRT_MACHINE_EPSILON;
1139566063dSJacob Faibussowitsch PetscCall(CkEigenSolutions(cklvl, A, il - 1, iu - 1, evals, evecs, tols));
114c4762a1bSJed Brown
115c4762a1bSJed Brown /* Free work space. */
1169566063dSJacob Faibussowitsch for (i = 0; i < nevs; i++) PetscCall(VecDestroy(&evecs[i]));
1179566063dSJacob Faibussowitsch PetscCall(PetscFree(evecs));
1189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_dense));
1199566063dSJacob Faibussowitsch PetscCall(PetscFree(work));
1209566063dSJacob Faibussowitsch if (TestSYEVX) PetscCall(PetscFree(evecs_array));
121c4762a1bSJed Brown
122c4762a1bSJed Brown /* Compute SVD: A_dense = U*SIGMA*transpose(V),
123c4762a1bSJed Brown JOBU=JOBV='S': the first min(m,n) columns of U and V are returned in the arrayU and arrayV; */
124c4762a1bSJed Brown /*==============================================================================================*/
125c4762a1bSJed Brown {
126c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */
127c4762a1bSJed Brown PetscScalar *arrayU, *arrayVT, *arrayErr, alpha = 1.0, beta = -1.0;
128c4762a1bSJed Brown Mat Err;
129c4762a1bSJed Brown PetscBLASInt minMN, maxMN, im, in;
130c4762a1bSJed Brown PetscInt j;
131c4762a1bSJed Brown PetscReal norm;
132c4762a1bSJed Brown
1339566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &A_dense));
134c4762a1bSJed Brown
1356497c311SBarry Smith PetscCall(PetscBLASIntCast(PetscMin(m, n), &minMN));
1366497c311SBarry Smith PetscCall(PetscBLASIntCast(PetscMax(m, n), &maxMN));
1376497c311SBarry Smith PetscCall(PetscBLASIntCast(5 * minMN + maxMN, &lwork));
1389566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(m * minMN, &arrayU, m * minMN, &arrayVT, m * minMN, &arrayErr, lwork, &work));
139c4762a1bSJed Brown
140c4762a1bSJed Brown /* Create matrix Err for checking error */
1419566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Err));
1429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Err, PETSC_DECIDE, PETSC_DECIDE, m, minMN));
1439566063dSJacob Faibussowitsch PetscCall(MatSetType(Err, MATSEQDENSE));
1449566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Err, (PetscScalar *)arrayErr));
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* Save A to arrayErr for checking accuracy later. arrayA will be destroyed by LAPACKgesvd_() */
1479566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(A_dense, &arrayA));
1489566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(arrayErr, arrayA, m * minMN));
149c4762a1bSJed Brown
1509566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &im));
1519566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &in));
152c4762a1bSJed Brown /* Compute A = U*SIGMA*VT */
153c4762a1bSJed Brown LAPACKgesvd_("S", "S", &im, &in, arrayA, &im, evals, arrayU, &minMN, arrayVT, &minMN, work, &lwork, &lierr);
1549566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(A_dense, &arrayA));
155c4762a1bSJed Brown if (!lierr) {
1569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " 1st 10 of %" PetscBLASInt_FMT " singular values: \n", minMN));
1579566063dSJacob Faibussowitsch for (i = 0; i < 10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g\n", i, (double)evals[i]));
158c4762a1bSJed Brown } else {
1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "LAPACKgesvd_ fails!"));
160c4762a1bSJed Brown }
161c4762a1bSJed Brown
162c4762a1bSJed Brown /* Check Err = (U*Sigma*V^T - A) using BLASgemm() */
163c4762a1bSJed Brown /* U = U*Sigma */
164c4762a1bSJed Brown for (j = 0; j < minMN; j++) { /* U[:,j] = sigma[j]*U[:,j] */
165c4762a1bSJed Brown for (i = 0; i < m; i++) arrayU[j * m + i] *= evals[j];
166c4762a1bSJed Brown }
167c4762a1bSJed Brown /* Err = U*VT - A = alpha*U*VT + beta*Err */
168c4762a1bSJed Brown BLASgemm_("N", "N", &im, &minMN, &minMN, &alpha, arrayU, &im, arrayVT, &minMN, &beta, arrayErr, &im);
1699566063dSJacob Faibussowitsch PetscCall(MatNorm(Err, NORM_FROBENIUS, &norm));
1709566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " || U*Sigma*VT - A || = %g\n", (double)norm));
171c4762a1bSJed Brown
1729566063dSJacob Faibussowitsch PetscCall(PetscFree4(arrayU, arrayVT, arrayErr, work));
1739566063dSJacob Faibussowitsch PetscCall(PetscFree(evals));
1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_dense));
1759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Err));
176c4762a1bSJed Brown }
177c4762a1bSJed Brown
1789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1799566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
180b122ec5aSJacob Faibussowitsch return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown /*------------------------------------------------
183c4762a1bSJed Brown Check the accuracy of the eigen solution
184c4762a1bSJed Brown ----------------------------------------------- */
185c4762a1bSJed Brown /*
186c4762a1bSJed Brown input:
187c4762a1bSJed Brown cklvl - check level:
188c4762a1bSJed Brown 1: check residual
189c4762a1bSJed Brown 2: 1 and check B-orthogonality locally
190c4762a1bSJed Brown A - matrix
191c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues
192c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process
193c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
194c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
195c4762a1bSJed Brown */
CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal * eval,Vec * evec,PetscReal * tols)196d71ae5a4SJacob Faibussowitsch PetscErrorCode CkEigenSolutions(PetscInt cklvl, Mat A, PetscInt il, PetscInt iu, PetscReal *eval, Vec *evec, PetscReal *tols)
197d71ae5a4SJacob Faibussowitsch {
1985f80ce2aSJacob Faibussowitsch PetscInt i, j, nev;
199c4762a1bSJed Brown Vec vt1, vt2; /* tmp vectors */
200c4762a1bSJed Brown PetscReal norm, tmp, dot, norm_max, dot_max;
201c4762a1bSJed Brown
202c4762a1bSJed Brown PetscFunctionBegin;
203c4762a1bSJed Brown nev = iu - il;
2043ba16761SJacob Faibussowitsch if (nev <= 0) PetscFunctionReturn(PETSC_SUCCESS);
205c4762a1bSJed Brown
206c4762a1bSJed Brown /*ierr = VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);*/
2079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0], &vt1));
2089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0], &vt2));
209c4762a1bSJed Brown
210c4762a1bSJed Brown switch (cklvl) {
211c4762a1bSJed Brown case 2:
212c4762a1bSJed Brown dot_max = 0.0;
213c4762a1bSJed Brown for (i = il; i < iu; i++) {
2149566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt1));
215c4762a1bSJed Brown for (j = il; j < iu; j++) {
2169566063dSJacob Faibussowitsch PetscCall(VecDot(evec[j], vt1, &dot));
217c4762a1bSJed Brown if (j == i) {
218c4762a1bSJed Brown dot = PetscAbsScalar(dot - 1);
219c4762a1bSJed Brown } else {
220c4762a1bSJed Brown dot = PetscAbsScalar(dot);
221c4762a1bSJed Brown }
222c4762a1bSJed Brown if (dot > dot_max) dot_max = dot;
223c4762a1bSJed Brown if (dot > tols[1]) {
2249566063dSJacob Faibussowitsch PetscCall(VecNorm(evec[i], NORM_INFINITY, &norm));
2259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n", i, j, (double)dot, (double)norm));
226c4762a1bSJed Brown }
227c4762a1bSJed Brown }
228c4762a1bSJed Brown }
2299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " max|(x_j^T*x_i) - delta_ji|: %g\n", (double)dot_max));
230c4762a1bSJed Brown
231c4762a1bSJed Brown case 1:
232c4762a1bSJed Brown norm_max = 0.0;
233c4762a1bSJed Brown for (i = il; i < iu; i++) {
2349566063dSJacob Faibussowitsch PetscCall(MatMult(A, evec[i], vt1));
2359566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt2));
236c4762a1bSJed Brown tmp = -eval[i];
2379566063dSJacob Faibussowitsch PetscCall(VecAXPY(vt1, tmp, vt2));
2389566063dSJacob Faibussowitsch PetscCall(VecNorm(vt1, NORM_INFINITY, &norm));
239c4762a1bSJed Brown norm = PetscAbsScalar(norm);
240c4762a1bSJed Brown if (norm > norm_max) norm_max = norm;
241c4762a1bSJed Brown /* sniff, and bark if necessary */
24248a46eb9SPierre Jolivet if (norm > tols[0]) PetscCall(PetscPrintf(PETSC_COMM_SELF, " residual violation: %" PetscInt_FMT ", resi: %g\n", i, (double)norm));
243c4762a1bSJed Brown }
2449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " max_resi: %g\n", (double)norm_max));
245c4762a1bSJed Brown break;
246d71ae5a4SJacob Faibussowitsch default:
247d71ae5a4SJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: cklvl=%" PetscInt_FMT " is not supported \n", cklvl));
248c4762a1bSJed Brown }
2499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt2));
2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt1));
2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown
254c4762a1bSJed Brown /*TEST
255c4762a1bSJed Brown
256c4762a1bSJed Brown build:
257c4762a1bSJed Brown requires: !complex
258c4762a1bSJed Brown
259c4762a1bSJed Brown test:
260dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
261c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small
262c4762a1bSJed Brown output_file: output/ex116_1.out
263c4762a1bSJed Brown
264c4762a1bSJed Brown test:
265c4762a1bSJed Brown suffix: 2
266dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
267c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -test_syev -check_symmetry
268c4762a1bSJed Brown
269c4762a1bSJed Brown TEST*/
270