xref: /petsc/src/mat/tests/ex120.c (revision be096a4675ce3b20dd7f9c195e6046e69d9955cf)
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 
9c4762a1bSJed Brown int main(int argc,char **args)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   Mat            A,A_dense,B;
12c4762a1bSJed Brown   Vec            *evecs;
13c4762a1bSJed Brown   PetscBool      flg,TestZHEEV=PETSC_TRUE,TestZHEEVX=PETSC_FALSE,TestZHEGV=PETSC_FALSE,TestZHEGVX=PETSC_FALSE;
14c4762a1bSJed Brown   PetscBool      isSymmetric;
15c4762a1bSJed Brown   PetscScalar    *arrayA,*arrayB,*evecs_array=NULL,*work;
16c4762a1bSJed Brown   PetscReal      *evals,*rwork;
17c4762a1bSJed Brown   PetscMPIInt    size;
18c4762a1bSJed Brown   PetscInt       m,i,j,cklvl=2;
19c4762a1bSJed Brown   PetscReal      vl,vu,abstol=1.e-8;
20c4762a1bSJed Brown   PetscBLASInt   nn,nevs,il,iu,*iwork,*ifail,lwork,lierr,bn,one=1;
21c4762a1bSJed Brown   PetscReal      tols[2];
22c4762a1bSJed Brown   PetscScalar    v,sigma2;
23c4762a1bSJed Brown   PetscRandom    rctx;
24c4762a1bSJed Brown   PetscReal      h2,sigma1 = 100.0;
25c4762a1bSJed Brown   PetscInt       dim,Ii,J,n = 6,use_random;
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
29*be096a46SBarry 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++) {
69c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
70c4762a1bSJed Brown     if (i>0) {
719566063dSJacob Faibussowitsch       J = Ii-n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
72c4762a1bSJed Brown     }
73c4762a1bSJed Brown     if (i<n-1) {
749566063dSJacob Faibussowitsch       J = Ii+n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown     if (j>0) {
779566063dSJacob Faibussowitsch       J = Ii-1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
78c4762a1bSJed Brown     }
79c4762a1bSJed Brown     if (j<n-1) {
809566063dSJacob Faibussowitsch       J = Ii+1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
81c4762a1bSJed Brown     }
829566063dSJacob Faibussowitsch     if (use_random) PetscCall(PetscRandomGetValue(rctx,&sigma2));
83c4762a1bSJed Brown     v    = 4.0 - sigma1*h2;
849566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown   /* make A complex Hermitian */
87c4762a1bSJed Brown   v    = sigma2*h2;
88c4762a1bSJed Brown   Ii   = 0; J = 1;
899566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
90c4762a1bSJed Brown   v    = -sigma2*h2;
919566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES));
929566063dSJacob Faibussowitsch   if (use_random) PetscCall(PetscRandomDestroy(&rctx));
939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
95c4762a1bSJed Brown   m    = n = dim;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Check whether A is symmetric */
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg));
99c4762a1bSJed Brown   if (flg) {
100c4762a1bSJed Brown     Mat Trans;
1019566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX, &Trans));
1029566063dSJacob Faibussowitsch     PetscCall(MatEqual(A, Trans, &isSymmetric));
1035f80ce2aSJacob Faibussowitsch     PetscCheck(isSymmetric,PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
1049566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Trans));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Convert aij matrix to MatSeqDense for LAPACK */
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg));
109c4762a1bSJed Brown   if (flg) {
1109566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A_dense));
111c4762a1bSJed Brown   } else {
1129566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&B));
1169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,dim,dim));
1179566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATSEQDENSE));
1189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
1199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
120c4762a1bSJed Brown   v    = 1.0;
121c4762a1bSJed Brown   for (Ii=0; Ii<dim; Ii++) {
1229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&Ii,1,&Ii,&v,ADD_VALUES));
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Solve standard eigenvalue problem: A*x = lambda*x */
126c4762a1bSJed Brown   /*===================================================*/
1279566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(2*n,&lwork));
1289566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n,&bn));
1299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&evals));
1309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lwork,&work));
1319566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A_dense,&arrayA));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   if (TestZHEEV) { /* test zheev() */
1349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n",m));
1359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(3*n-2,&rwork));
136c4762a1bSJed Brown     LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,rwork,&lierr);
1379566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown     evecs_array = arrayA;
140c4762a1bSJed Brown     nevs        = m;
141c4762a1bSJed Brown     il          =1; iu=m;
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown   if (TestZHEEVX) {
144c4762a1bSJed Brown     il   = 1;
1459566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast((0.2*m),&iu));
1469566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu));
1479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m*n+1,&evecs_array));
1489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(7*n+1,&rwork));
1499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*n+1,&iwork));
1509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n+1,&ifail));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     /* in the case "I", vl and vu are not referenced */
153c4762a1bSJed Brown     vl = 0.0; vu = 8.0;
1549566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n,&nn));
155c4762a1bSJed Brown     LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&nn,work,&lwork,rwork,iwork,ifail,&lierr);
1569566063dSJacob Faibussowitsch     PetscCall(PetscFree(iwork));
1579566063dSJacob Faibussowitsch     PetscCall(PetscFree(ifail));
1589566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown   if (TestZHEGV) {
1619566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute all %" PetscInt_FMT " eigensolutions...\n",m));
1629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(3*n+1,&rwork));
1639566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(B,&arrayB));
164c4762a1bSJed Brown     LAPACKsygv_(&one,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,rwork,&lierr);
165c4762a1bSJed Brown     evecs_array = arrayA;
166c4762a1bSJed Brown     nevs        = m;
167c4762a1bSJed Brown     il          = 1; iu=m;
1689566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(B,&arrayB));
1699566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown   if (TestZHEGVX) {
172c4762a1bSJed Brown     il   = 1;
1739566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast((0.2*m),&iu));
1749566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute %d to %d-th eigensolutions...\n",il,iu));
1759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m*n+1,&evecs_array));
1769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(6*n+1,&iwork));
177c4762a1bSJed Brown     ifail = iwork + 5*n;
1789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(7*n+1,&rwork));
1799566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(B,&arrayB));
180c4762a1bSJed Brown     vl    = 0.0; vu = 8.0;
1819566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n,&nn));
182c4762a1bSJed 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);
1839566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(B,&arrayB));
1849566063dSJacob Faibussowitsch     PetscCall(PetscFree(iwork));
1859566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
186c4762a1bSJed Brown   }
1879566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(A_dense,&arrayA));
1885f80ce2aSJacob Faibussowitsch   PetscCheck(nevs > 0,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* View evals */
1919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg));
192c4762a1bSJed Brown   if (flg) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," %d evals: \n",nevs));
1949566063dSJacob Faibussowitsch     for (i=0; i<nevs; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "  %g\n",i+il,(double)evals[i]));
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* Check residuals and orthogonality */
1989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevs+1,&evecs));
199c4762a1bSJed Brown   for (i=0; i<nevs; i++) {
2009566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF,&evecs[i]));
2019566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(evecs[i],PETSC_DECIDE,n));
2029566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(evecs[i]));
2039566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(evecs[i],evecs_array+i*n));
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   tols[0] = PETSC_SQRT_MACHINE_EPSILON;  tols[1] = PETSC_SQRT_MACHINE_EPSILON;
2079566063dSJacob Faibussowitsch   PetscCall(CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols));
2089566063dSJacob Faibussowitsch   for (i=0; i<nevs; i++) PetscCall(VecDestroy(&evecs[i]));
2099566063dSJacob Faibussowitsch   PetscCall(PetscFree(evecs));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Free work space. */
212c4762a1bSJed Brown   if (TestZHEEVX || TestZHEGVX) {
2139566063dSJacob Faibussowitsch     PetscCall(PetscFree(evecs_array));
214c4762a1bSJed Brown   }
2159566063dSJacob Faibussowitsch   PetscCall(PetscFree(evals));
2169566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
2179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_dense));
2189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
221b122ec5aSJacob Faibussowitsch   return 0;
222c4762a1bSJed Brown }
223c4762a1bSJed Brown /*------------------------------------------------
224c4762a1bSJed Brown   Check the accuracy of the eigen solution
225c4762a1bSJed Brown   ----------------------------------------------- */
226c4762a1bSJed Brown /*
227c4762a1bSJed Brown   input:
228c4762a1bSJed Brown      cklvl      - check level:
229c4762a1bSJed Brown                     1: check residual
230c4762a1bSJed Brown                     2: 1 and check B-orthogonality locally
231c4762a1bSJed Brown      A          - matrix
232c4762a1bSJed Brown      il,iu      - lower and upper index bound of eigenvalues
233c4762a1bSJed Brown      eval, evec - eigenvalues and eigenvectors stored in this process
234c4762a1bSJed Brown      tols[0]    - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
235c4762a1bSJed Brown      tols[1]    - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
236c4762a1bSJed Brown */
237c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
238c4762a1bSJed Brown {
2395f80ce2aSJacob Faibussowitsch   PetscInt    i,j,nev;
240c4762a1bSJed Brown   Vec         vt1,vt2;  /* tmp vectors */
241c4762a1bSJed Brown   PetscReal   norm,tmp,norm_max,dot_max,rdot;
242c4762a1bSJed Brown   PetscScalar dot;
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   PetscFunctionBegin;
245c4762a1bSJed Brown   nev = iu - il;
246c4762a1bSJed Brown   if (nev <= 0) PetscFunctionReturn(0);
247c4762a1bSJed Brown 
2489566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(evec[0],&vt1));
2499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(evec[0],&vt2));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   switch (cklvl) {
252c4762a1bSJed Brown   case 2:
253c4762a1bSJed Brown     dot_max = 0.0;
254c4762a1bSJed Brown     for (i = il; i<iu; i++) {
2559566063dSJacob Faibussowitsch       PetscCall(VecCopy(evec[i], vt1));
256c4762a1bSJed Brown       for (j=il; j<iu; j++) {
2579566063dSJacob Faibussowitsch         PetscCall(VecDot(evec[j],vt1,&dot));
258c4762a1bSJed Brown         if (j == i) {
259c4762a1bSJed Brown           rdot = PetscAbsScalar(dot - (PetscScalar)1.0);
260c4762a1bSJed Brown         } else {
261c4762a1bSJed Brown           rdot = PetscAbsScalar(dot);
262c4762a1bSJed Brown         }
263c4762a1bSJed Brown         if (rdot > dot_max) dot_max = rdot;
264c4762a1bSJed Brown         if (rdot > tols[1]) {
2659566063dSJacob Faibussowitsch           PetscCall(VecNorm(evec[i],NORM_INFINITY,&norm));
2669566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF,"|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n",i,j,(double)rdot,(double)norm));
267c4762a1bSJed Brown         }
268c4762a1bSJed Brown       }
269c4762a1bSJed Brown     }
2709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"    max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max));
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   case 1:
273c4762a1bSJed Brown     norm_max = 0.0;
274c4762a1bSJed Brown     for (i = il; i< iu; i++) {
2759566063dSJacob Faibussowitsch       PetscCall(MatMult(A, evec[i], vt1));
2769566063dSJacob Faibussowitsch       PetscCall(VecCopy(evec[i], vt2));
277c4762a1bSJed Brown       tmp  = -eval[i];
2789566063dSJacob Faibussowitsch       PetscCall(VecAXPY(vt1,tmp,vt2));
2799566063dSJacob Faibussowitsch       PetscCall(VecNorm(vt1, NORM_INFINITY, &norm));
280c4762a1bSJed Brown       norm = PetscAbs(norm);
281c4762a1bSJed Brown       if (norm > norm_max) norm_max = norm;
282c4762a1bSJed Brown       /* sniff, and bark if necessary */
283c4762a1bSJed Brown       if (norm > tols[0]) {
2849566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  residual violation: %" PetscInt_FMT ", resi: %g\n",i, norm));
285c4762a1bSJed Brown       }
286c4762a1bSJed Brown     }
2879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"    max_resi:                    %g\n", (double)norm_max));
288c4762a1bSJed Brown     break;
289c4762a1bSJed Brown   default:
2909566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%" PetscInt_FMT " is not supported \n",cklvl));
291c4762a1bSJed Brown   }
2929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vt2));
2939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vt1));
294c4762a1bSJed Brown   PetscFunctionReturn(0);
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown /*TEST
298c4762a1bSJed Brown 
299c4762a1bSJed Brown    build:
300c4762a1bSJed Brown       requires: complex
301c4762a1bSJed Brown 
302c4762a1bSJed Brown    test:
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown       suffix: 2
306c4762a1bSJed Brown       args: -test_zheevx
307c4762a1bSJed Brown 
308c4762a1bSJed Brown    test:
309c4762a1bSJed Brown       suffix: 3
310c4762a1bSJed Brown       args: -test_zhegv
311c4762a1bSJed Brown 
312c4762a1bSJed Brown    test:
313c4762a1bSJed Brown       suffix: 4
314c4762a1bSJed Brown       args: -test_zhegvx
315c4762a1bSJed Brown 
316c4762a1bSJed Brown TEST*/
317