xref: /petsc/src/mat/tests/ex120.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Test LAPACK routine ZHEEV, ZHEEVX, ZHEGV and ZHEGVX. \n\
2c4762a1bSJed Brown ZHEEV computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #include <petscblaslapack.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown extern PetscErrorCode CkEigenSolutions(PetscInt, Mat, PetscInt, PetscInt, PetscReal *, Vec *, PetscReal *);
8c4762a1bSJed Brown 
99371c9d4SSatish Balay int main(int argc, char **args) {
10c4762a1bSJed Brown   Mat          A, A_dense, B;
11c4762a1bSJed Brown   Vec         *evecs;
12c4762a1bSJed Brown   PetscBool    flg, TestZHEEV = PETSC_TRUE, TestZHEEVX = PETSC_FALSE, TestZHEGV = PETSC_FALSE, TestZHEGVX = PETSC_FALSE;
13c4762a1bSJed Brown   PetscBool    isSymmetric;
14c4762a1bSJed Brown   PetscScalar *arrayA, *arrayB, *evecs_array = NULL, *work;
15c4762a1bSJed Brown   PetscReal   *evals, *rwork;
16c4762a1bSJed Brown   PetscMPIInt  size;
17c4762a1bSJed Brown   PetscInt     m, i, j, cklvl = 2;
18c4762a1bSJed Brown   PetscReal    vl, vu, abstol = 1.e-8;
19c4762a1bSJed Brown   PetscBLASInt nn, nevs, il, iu, *iwork, *ifail, lwork, lierr, bn, one = 1;
20c4762a1bSJed Brown   PetscReal    tols[2];
21c4762a1bSJed Brown   PetscScalar  v, sigma2;
22c4762a1bSJed Brown   PetscRandom  rctx;
23c4762a1bSJed Brown   PetscReal    h2, sigma1 = 100.0;
24c4762a1bSJed Brown   PetscInt     dim, Ii, J, n = 6, use_random;
25c4762a1bSJed Brown 
26327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
29be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_zheevx", &flg));
32c4762a1bSJed Brown   if (flg) {
33c4762a1bSJed Brown     TestZHEEV  = PETSC_FALSE;
34c4762a1bSJed Brown     TestZHEEVX = PETSC_TRUE;
35c4762a1bSJed Brown   }
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_zhegv", &flg));
37c4762a1bSJed Brown   if (flg) {
38c4762a1bSJed Brown     TestZHEEV = PETSC_FALSE;
39c4762a1bSJed Brown     TestZHEGV = PETSC_TRUE;
40c4762a1bSJed Brown   }
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_zhegvx", &flg));
42c4762a1bSJed Brown   if (flg) {
43c4762a1bSJed Brown     TestZHEEV  = PETSC_FALSE;
44c4762a1bSJed Brown     TestZHEGVX = PETSC_TRUE;
45c4762a1bSJed Brown   }
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-sigma1", &sigma1, NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
49c4762a1bSJed Brown   dim = n * n;
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
529566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, dim, dim));
539566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSEQDENSE));
549566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-norandom", &flg));
58c4762a1bSJed Brown   if (flg) use_random = 0;
59c4762a1bSJed Brown   else use_random = 1;
60c4762a1bSJed Brown   if (use_random) {
619566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rctx));
629566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(rctx));
639566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(rctx, 0.0, PETSC_i));
64c4762a1bSJed Brown   } else {
65c4762a1bSJed Brown     sigma2 = 10.0 * PETSC_i;
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown   h2 = 1.0 / ((n + 1) * (n + 1));
68c4762a1bSJed Brown   for (Ii = 0; Ii < dim; Ii++) {
699371c9d4SSatish Balay     v = -1.0;
709371c9d4SSatish Balay     i = Ii / n;
719371c9d4SSatish Balay     j = Ii - i * n;
72c4762a1bSJed Brown     if (i > 0) {
739371c9d4SSatish Balay       J = Ii - n;
749371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown     if (i < n - 1) {
779371c9d4SSatish Balay       J = Ii + n;
789371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown     if (j > 0) {
819371c9d4SSatish Balay       J = Ii - 1;
829371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
83c4762a1bSJed Brown     }
84c4762a1bSJed Brown     if (j < n - 1) {
859371c9d4SSatish Balay       J = Ii + 1;
869371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
87c4762a1bSJed Brown     }
889566063dSJacob Faibussowitsch     if (use_random) PetscCall(PetscRandomGetValue(rctx, &sigma2));
89c4762a1bSJed Brown     v = 4.0 - sigma1 * h2;
909566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   /* make A complex Hermitian */
93c4762a1bSJed Brown   v  = sigma2 * h2;
949371c9d4SSatish Balay   Ii = 0;
959371c9d4SSatish Balay   J  = 1;
969566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
97c4762a1bSJed Brown   v = -sigma2 * h2;
989566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &J, 1, &Ii, &v, ADD_VALUES));
999566063dSJacob Faibussowitsch   if (use_random) PetscCall(PetscRandomDestroy(&rctx));
1009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
102c4762a1bSJed Brown   m = n = dim;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Check whether A is symmetric */
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetry", &flg));
106c4762a1bSJed Brown   if (flg) {
107c4762a1bSJed Brown     Mat Trans;
1089566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Trans));
1099566063dSJacob Faibussowitsch     PetscCall(MatEqual(A, Trans, &isSymmetric));
1105f80ce2aSJacob Faibussowitsch     PetscCheck(isSymmetric, PETSC_COMM_SELF, PETSC_ERR_USER, "A must be symmetric");
1119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Trans));
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Convert aij matrix to MatSeqDense for LAPACK */
1159566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &flg));
116c4762a1bSJed Brown   if (flg) {
1179566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A_dense));
118c4762a1bSJed Brown   } else {
1199566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &A_dense));
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &B));
1239566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, dim, dim));
1249566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATSEQDENSE));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
1269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
127c4762a1bSJed Brown   v = 1.0;
128*48a46eb9SPierre Jolivet   for (Ii = 0; Ii < dim; Ii++) PetscCall(MatSetValues(B, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Solve standard eigenvalue problem: A*x = lambda*x */
131c4762a1bSJed Brown   /*===================================================*/
1329566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(2 * n, &lwork));
1339566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &bn));
1349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &evals));
1359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lwork, &work));
1369566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A_dense, &arrayA));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   if (TestZHEEV) { /* test zheev() */
1399566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n", m));
1409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(3 * n - 2, &rwork));
141c4762a1bSJed Brown     LAPACKsyev_("V", "U", &bn, arrayA, &bn, evals, work, &lwork, rwork, &lierr);
1429566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown     evecs_array = arrayA;
145c4762a1bSJed Brown     nevs        = m;
1469371c9d4SSatish Balay     il          = 1;
1479371c9d4SSatish Balay     iu          = m;
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown   if (TestZHEEVX) {
150c4762a1bSJed Brown     il = 1;
1519566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast((0.2 * m), &iu));
1529566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " LAPACKsyevx: compute %d to %d-th eigensolutions...\n", il, iu));
1539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m * n + 1, &evecs_array));
1549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(7 * n + 1, &rwork));
1559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * n + 1, &iwork));
1569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n + 1, &ifail));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     /* in the case "I", vl and vu are not referenced */
1599371c9d4SSatish Balay     vl = 0.0;
1609371c9d4SSatish Balay     vu = 8.0;
1619566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &nn));
162c4762a1bSJed Brown     LAPACKsyevx_("V", "I", "U", &bn, arrayA, &bn, &vl, &vu, &il, &iu, &abstol, &nevs, evals, evecs_array, &nn, work, &lwork, rwork, iwork, ifail, &lierr);
1639566063dSJacob Faibussowitsch     PetscCall(PetscFree(iwork));
1649566063dSJacob Faibussowitsch     PetscCall(PetscFree(ifail));
1659566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
166c4762a1bSJed Brown   }
167c4762a1bSJed Brown   if (TestZHEGV) {
1689566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " LAPACKsygv: compute all %" PetscInt_FMT " eigensolutions...\n", m));
1699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(3 * n + 1, &rwork));
1709566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(B, &arrayB));
171c4762a1bSJed Brown     LAPACKsygv_(&one, "V", "U", &bn, arrayA, &bn, arrayB, &bn, evals, work, &lwork, rwork, &lierr);
172c4762a1bSJed Brown     evecs_array = arrayA;
173c4762a1bSJed Brown     nevs        = m;
1749371c9d4SSatish Balay     il          = 1;
1759371c9d4SSatish Balay     iu          = m;
1769566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(B, &arrayB));
1779566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   if (TestZHEGVX) {
180c4762a1bSJed Brown     il = 1;
1819566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast((0.2 * m), &iu));
1829566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " LAPACKsygv: compute %d to %d-th eigensolutions...\n", il, iu));
1839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m * n + 1, &evecs_array));
1849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(6 * n + 1, &iwork));
185c4762a1bSJed Brown     ifail = iwork + 5 * n;
1869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(7 * n + 1, &rwork));
1879566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(B, &arrayB));
1889371c9d4SSatish Balay     vl = 0.0;
1899371c9d4SSatish Balay     vu = 8.0;
1909566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &nn));
191c4762a1bSJed Brown     LAPACKsygvx_(&one, "V", "I", "U", &bn, arrayA, &bn, arrayB, &bn, &vl, &vu, &il, &iu, &abstol, &nevs, evals, evecs_array, &nn, work, &lwork, rwork, iwork, ifail, &lierr);
1929566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(B, &arrayB));
1939566063dSJacob Faibussowitsch     PetscCall(PetscFree(iwork));
1949566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
195c4762a1bSJed Brown   }
1969566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(A_dense, &arrayA));
1975f80ce2aSJacob Faibussowitsch   PetscCheck(nevs > 0, PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* View evals */
2009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-eig_view", &flg));
201c4762a1bSJed Brown   if (flg) {
2029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %d evals: \n", nevs));
2039566063dSJacob Faibussowitsch     for (i = 0; i < nevs; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "  %g\n", i + il, (double)evals[i]));
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* Check residuals and orthogonality */
2079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevs + 1, &evecs));
208c4762a1bSJed Brown   for (i = 0; i < nevs; i++) {
2099566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &evecs[i]));
2109566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(evecs[i], PETSC_DECIDE, n));
2119566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(evecs[i]));
2129566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(evecs[i], evecs_array + i * n));
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown 
2159371c9d4SSatish Balay   tols[0] = PETSC_SQRT_MACHINE_EPSILON;
2169371c9d4SSatish Balay   tols[1] = PETSC_SQRT_MACHINE_EPSILON;
2179566063dSJacob Faibussowitsch   PetscCall(CkEigenSolutions(cklvl, A, il - 1, iu - 1, evals, evecs, tols));
2189566063dSJacob Faibussowitsch   for (i = 0; i < nevs; i++) PetscCall(VecDestroy(&evecs[i]));
2199566063dSJacob Faibussowitsch   PetscCall(PetscFree(evecs));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Free work space. */
222*48a46eb9SPierre Jolivet   if (TestZHEEVX || TestZHEGVX) PetscCall(PetscFree(evecs_array));
2239566063dSJacob Faibussowitsch   PetscCall(PetscFree(evals));
2249566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
2259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_dense));
2269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
229b122ec5aSJacob Faibussowitsch   return 0;
230c4762a1bSJed Brown }
231c4762a1bSJed Brown /*------------------------------------------------
232c4762a1bSJed Brown   Check the accuracy of the eigen solution
233c4762a1bSJed Brown   ----------------------------------------------- */
234c4762a1bSJed Brown /*
235c4762a1bSJed Brown   input:
236c4762a1bSJed Brown      cklvl      - check level:
237c4762a1bSJed Brown                     1: check residual
238c4762a1bSJed Brown                     2: 1 and check B-orthogonality locally
239c4762a1bSJed Brown      A          - matrix
240c4762a1bSJed Brown      il,iu      - lower and upper index bound of eigenvalues
241c4762a1bSJed Brown      eval, evec - eigenvalues and eigenvectors stored in this process
242c4762a1bSJed Brown      tols[0]    - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
243c4762a1bSJed Brown      tols[1]    - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
244c4762a1bSJed Brown */
2459371c9d4SSatish Balay PetscErrorCode CkEigenSolutions(PetscInt cklvl, Mat A, PetscInt il, PetscInt iu, PetscReal *eval, Vec *evec, PetscReal *tols) {
2465f80ce2aSJacob Faibussowitsch   PetscInt    i, j, nev;
247c4762a1bSJed Brown   Vec         vt1, vt2; /* tmp vectors */
248c4762a1bSJed Brown   PetscReal   norm, tmp, norm_max, dot_max, rdot;
249c4762a1bSJed Brown   PetscScalar dot;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   PetscFunctionBegin;
252c4762a1bSJed Brown   nev = iu - il;
253c4762a1bSJed Brown   if (nev <= 0) PetscFunctionReturn(0);
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(evec[0], &vt1));
2569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(evec[0], &vt2));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   switch (cklvl) {
259c4762a1bSJed Brown   case 2:
260c4762a1bSJed Brown     dot_max = 0.0;
261c4762a1bSJed Brown     for (i = il; i < iu; i++) {
2629566063dSJacob Faibussowitsch       PetscCall(VecCopy(evec[i], vt1));
263c4762a1bSJed Brown       for (j = il; j < iu; j++) {
2649566063dSJacob Faibussowitsch         PetscCall(VecDot(evec[j], vt1, &dot));
265c4762a1bSJed Brown         if (j == i) {
266c4762a1bSJed Brown           rdot = PetscAbsScalar(dot - (PetscScalar)1.0);
267c4762a1bSJed Brown         } else {
268c4762a1bSJed Brown           rdot = PetscAbsScalar(dot);
269c4762a1bSJed Brown         }
270c4762a1bSJed Brown         if (rdot > dot_max) dot_max = rdot;
271c4762a1bSJed Brown         if (rdot > tols[1]) {
2729566063dSJacob Faibussowitsch           PetscCall(VecNorm(evec[i], NORM_INFINITY, &norm));
2739566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n", i, j, (double)rdot, (double)norm));
274c4762a1bSJed Brown         }
275c4762a1bSJed Brown       }
276c4762a1bSJed Brown     }
2779566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "    max|(x_j^T*x_i) - delta_ji|: %g\n", (double)dot_max));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   case 1:
280c4762a1bSJed Brown     norm_max = 0.0;
281c4762a1bSJed Brown     for (i = il; i < iu; i++) {
2829566063dSJacob Faibussowitsch       PetscCall(MatMult(A, evec[i], vt1));
2839566063dSJacob Faibussowitsch       PetscCall(VecCopy(evec[i], vt2));
284c4762a1bSJed Brown       tmp = -eval[i];
2859566063dSJacob Faibussowitsch       PetscCall(VecAXPY(vt1, tmp, vt2));
2869566063dSJacob Faibussowitsch       PetscCall(VecNorm(vt1, NORM_INFINITY, &norm));
287c4762a1bSJed Brown       norm = PetscAbs(norm);
288c4762a1bSJed Brown       if (norm > norm_max) norm_max = norm;
289c4762a1bSJed Brown       /* sniff, and bark if necessary */
290*48a46eb9SPierre Jolivet       if (norm > tols[0]) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  residual violation: %" PetscInt_FMT ", resi: %g\n", i, norm));
291c4762a1bSJed Brown     }
2929566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "    max_resi:                    %g\n", (double)norm_max));
293c4762a1bSJed Brown     break;
2949371c9d4SSatish Balay   default: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: cklvl=%" PetscInt_FMT " is not supported \n", cklvl));
295c4762a1bSJed Brown   }
2969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vt2));
2979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vt1));
298c4762a1bSJed Brown   PetscFunctionReturn(0);
299c4762a1bSJed Brown }
300c4762a1bSJed Brown 
301c4762a1bSJed Brown /*TEST
302c4762a1bSJed Brown 
303c4762a1bSJed Brown    build:
304c4762a1bSJed Brown       requires: complex
305c4762a1bSJed Brown 
306c4762a1bSJed Brown    test:
307c4762a1bSJed Brown 
308c4762a1bSJed Brown    test:
309c4762a1bSJed Brown       suffix: 2
310c4762a1bSJed Brown       args: -test_zheevx
311c4762a1bSJed Brown 
312c4762a1bSJed Brown    test:
313c4762a1bSJed Brown       suffix: 3
314c4762a1bSJed Brown       args: -test_zhegv
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    test:
317c4762a1bSJed Brown       suffix: 4
318c4762a1bSJed Brown       args: -test_zhegvx
319c4762a1bSJed Brown 
320c4762a1bSJed Brown TEST*/
321