xref: /petsc/src/mat/tests/ex120.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
15c4762a1bSJed Brown   PetscBool      isSymmetric;
16c4762a1bSJed Brown   PetscScalar    *arrayA,*arrayB,*evecs_array=NULL,*work;
17c4762a1bSJed Brown   PetscReal      *evals,*rwork;
18c4762a1bSJed Brown   PetscMPIInt    size;
19c4762a1bSJed Brown   PetscInt       m,i,j,cklvl=2;
20c4762a1bSJed Brown   PetscReal      vl,vu,abstol=1.e-8;
21c4762a1bSJed Brown   PetscBLASInt   nn,nevs,il,iu,*iwork,*ifail,lwork,lierr,bn,one=1;
22c4762a1bSJed Brown   PetscReal      tols[2];
23c4762a1bSJed Brown   PetscScalar    v,sigma2;
24c4762a1bSJed Brown   PetscRandom    rctx;
25c4762a1bSJed Brown   PetscReal      h2,sigma1 = 100.0;
26c4762a1bSJed Brown   PetscInt       dim,Ii,J,n = 6,use_random;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
30*5f80ce2aSJacob Faibussowitsch   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
31c4762a1bSJed Brown 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zheevx", &flg));
33c4762a1bSJed Brown   if (flg) {
34c4762a1bSJed Brown     TestZHEEV  = PETSC_FALSE;
35c4762a1bSJed Brown     TestZHEEVX = PETSC_TRUE;
36c4762a1bSJed Brown   }
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zhegv", &flg));
38c4762a1bSJed Brown   if (flg) {
39c4762a1bSJed Brown     TestZHEEV = PETSC_FALSE;
40c4762a1bSJed Brown     TestZHEGV = PETSC_TRUE;
41c4762a1bSJed Brown   }
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zhegvx", &flg));
43c4762a1bSJed Brown   if (flg) {
44c4762a1bSJed Brown     TestZHEEV  = PETSC_FALSE;
45c4762a1bSJed Brown     TestZHEGVX = PETSC_TRUE;
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown 
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
50c4762a1bSJed Brown   dim  = n*n;
51c4762a1bSJed Brown 
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQDENSE));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
57c4762a1bSJed Brown 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-norandom",&flg));
59c4762a1bSJed Brown   if (flg) use_random = 0;
60c4762a1bSJed Brown   else     use_random = 1;
61c4762a1bSJed Brown   if (use_random) {
62*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rctx));
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetFromOptions(rctx));
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetInterval(rctx,0.0,PETSC_i));
65c4762a1bSJed Brown   } else {
66c4762a1bSJed Brown     sigma2 = 10.0*PETSC_i;
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown   h2 = 1.0/((n+1)*(n+1));
69c4762a1bSJed Brown   for (Ii=0; Ii<dim; Ii++) {
70c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
71c4762a1bSJed Brown     if (i>0) {
72*5f80ce2aSJacob Faibussowitsch       J = Ii-n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown     if (i<n-1) {
75*5f80ce2aSJacob Faibussowitsch       J = Ii+n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
76c4762a1bSJed Brown     }
77c4762a1bSJed Brown     if (j>0) {
78*5f80ce2aSJacob Faibussowitsch       J = Ii-1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown     if (j<n-1) {
81*5f80ce2aSJacob Faibussowitsch       J = Ii+1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
82c4762a1bSJed Brown     }
83*5f80ce2aSJacob Faibussowitsch     if (use_random) CHKERRQ(PetscRandomGetValue(rctx,&sigma2));
84c4762a1bSJed Brown     v    = 4.0 - sigma1*h2;
85*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES));
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   /* make A complex Hermitian */
88c4762a1bSJed Brown   v    = sigma2*h2;
89c4762a1bSJed Brown   Ii   = 0; J = 1;
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));
91c4762a1bSJed Brown   v    = -sigma2*h2;
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES));
93*5f80ce2aSJacob Faibussowitsch   if (use_random) CHKERRQ(PetscRandomDestroy(&rctx));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
96c4762a1bSJed Brown   m    = n = dim;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* Check whether A is symmetric */
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg));
100c4762a1bSJed Brown   if (flg) {
101c4762a1bSJed Brown     Mat Trans;
102*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Trans));
103*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatEqual(A, Trans, &isSymmetric));
104*5f80ce2aSJacob Faibussowitsch     PetscCheck(isSymmetric,PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Trans));
106c4762a1bSJed Brown   }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Convert aij matrix to MatSeqDense for LAPACK */
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg));
110c4762a1bSJed Brown   if (flg) {
111*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A_dense));
112c4762a1bSJed Brown   } else {
113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense));
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,dim,dim));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATSEQDENSE));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
121c4762a1bSJed Brown   v    = 1.0;
122c4762a1bSJed Brown   for (Ii=0; Ii<dim; Ii++) {
123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,1,&Ii,1,&Ii,&v,ADD_VALUES));
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Solve standard eigenvalue problem: A*x = lambda*x */
127c4762a1bSJed Brown   /*===================================================*/
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBLASIntCast(2*n,&lwork));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBLASIntCast(n,&bn));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&evals));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(lwork,&work));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(A_dense,&arrayA));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   if (TestZHEEV) { /* test zheev() */
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n",m));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(3*n-2,&rwork));
137c4762a1bSJed Brown     LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,rwork,&lierr);
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rwork));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     evecs_array = arrayA;
141c4762a1bSJed Brown     nevs        = m;
142c4762a1bSJed Brown     il          =1; iu=m;
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown   if (TestZHEEVX) {
145c4762a1bSJed Brown     il   = 1;
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast((0.2*m),&iu));
147*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu));
148*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(m*n+1,&evecs_array));
149*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(7*n+1,&rwork));
150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(5*n+1,&iwork));
151*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n+1,&ifail));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown     /* in the case "I", vl and vu are not referenced */
154c4762a1bSJed Brown     vl = 0.0; vu = 8.0;
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&nn));
156c4762a1bSJed Brown     LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&nn,work,&lwork,rwork,iwork,ifail,&lierr);
157*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(iwork));
158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(ifail));
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rwork));
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown   if (TestZHEGV) {
162*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute all %" PetscInt_FMT " eigensolutions...\n",m));
163*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(3*n+1,&rwork));
164*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseGetArray(B,&arrayB));
165c4762a1bSJed Brown     LAPACKsygv_(&one,"V","U",&bn,arrayA,&bn,arrayB,&bn,evals,work,&lwork,rwork,&lierr);
166c4762a1bSJed Brown     evecs_array = arrayA;
167c4762a1bSJed Brown     nevs        = m;
168c4762a1bSJed Brown     il          = 1; iu=m;
169*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseRestoreArray(B,&arrayB));
170*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rwork));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown   if (TestZHEGVX) {
173c4762a1bSJed Brown     il   = 1;
174*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast((0.2*m),&iu));
175*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute %d to %d-th eigensolutions...\n",il,iu));
176*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(m*n+1,&evecs_array));
177*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(6*n+1,&iwork));
178c4762a1bSJed Brown     ifail = iwork + 5*n;
179*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(7*n+1,&rwork));
180*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseGetArray(B,&arrayB));
181c4762a1bSJed Brown     vl    = 0.0; vu = 8.0;
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&nn));
183c4762a1bSJed 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);
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseRestoreArray(B,&arrayB));
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(iwork));
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rwork));
187c4762a1bSJed Brown   }
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(A_dense,&arrayA));
189*5f80ce2aSJacob Faibussowitsch   PetscCheck(nevs > 0,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs);
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* View evals */
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg));
193c4762a1bSJed Brown   if (flg) {
194*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %d evals: \n",nevs));
195*5f80ce2aSJacob Faibussowitsch     for (i=0; i<nevs; i++) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "  %g\n",i+il,(double)evals[i]));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Check residuals and orthogonality */
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nevs+1,&evecs));
200c4762a1bSJed Brown   for (i=0; i<nevs; i++) {
201*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreate(PETSC_COMM_SELF,&evecs[i]));
202*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetSizes(evecs[i],PETSC_DECIDE,n));
203*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetFromOptions(evecs[i]));
204*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPlaceArray(evecs[i],evecs_array+i*n));
205c4762a1bSJed Brown   }
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   tols[0] = PETSC_SQRT_MACHINE_EPSILON;  tols[1] = PETSC_SQRT_MACHINE_EPSILON;
208*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols));
209*5f80ce2aSJacob Faibussowitsch   for (i=0; i<nevs; i++) CHKERRQ(VecDestroy(&evecs[i]));
210*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(evecs));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* Free work space. */
213c4762a1bSJed Brown   if (TestZHEEVX || TestZHEGVX) {
214*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(evecs_array));
215c4762a1bSJed Brown   }
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(evals));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(work));
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A_dense));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
220*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
221c4762a1bSJed Brown   ierr = PetscFinalize();
222c4762a1bSJed Brown   return ierr;
223c4762a1bSJed Brown }
224c4762a1bSJed Brown /*------------------------------------------------
225c4762a1bSJed Brown   Check the accuracy of the eigen solution
226c4762a1bSJed Brown   ----------------------------------------------- */
227c4762a1bSJed Brown /*
228c4762a1bSJed Brown   input:
229c4762a1bSJed Brown      cklvl      - check level:
230c4762a1bSJed Brown                     1: check residual
231c4762a1bSJed Brown                     2: 1 and check B-orthogonality locally
232c4762a1bSJed Brown      A          - matrix
233c4762a1bSJed Brown      il,iu      - lower and upper index bound of eigenvalues
234c4762a1bSJed Brown      eval, evec - eigenvalues and eigenvectors stored in this process
235c4762a1bSJed Brown      tols[0]    - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
236c4762a1bSJed Brown      tols[1]    - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
237c4762a1bSJed Brown */
238c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
239c4762a1bSJed Brown {
240*5f80ce2aSJacob Faibussowitsch   PetscInt    i,j,nev;
241c4762a1bSJed Brown   Vec         vt1,vt2;  /* tmp vectors */
242c4762a1bSJed Brown   PetscReal   norm,tmp,norm_max,dot_max,rdot;
243c4762a1bSJed Brown   PetscScalar dot;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   PetscFunctionBegin;
246c4762a1bSJed Brown   nev = iu - il;
247c4762a1bSJed Brown   if (nev <= 0) PetscFunctionReturn(0);
248c4762a1bSJed Brown 
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(evec[0],&vt1));
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(evec[0],&vt2));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   switch (cklvl) {
253c4762a1bSJed Brown   case 2:
254c4762a1bSJed Brown     dot_max = 0.0;
255c4762a1bSJed Brown     for (i = il; i<iu; i++) {
256*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(evec[i], vt1));
257c4762a1bSJed Brown       for (j=il; j<iu; j++) {
258*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDot(evec[j],vt1,&dot));
259c4762a1bSJed Brown         if (j == i) {
260c4762a1bSJed Brown           rdot = PetscAbsScalar(dot - (PetscScalar)1.0);
261c4762a1bSJed Brown         } else {
262c4762a1bSJed Brown           rdot = PetscAbsScalar(dot);
263c4762a1bSJed Brown         }
264c4762a1bSJed Brown         if (rdot > dot_max) dot_max = rdot;
265c4762a1bSJed Brown         if (rdot > tols[1]) {
266*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecNorm(evec[i],NORM_INFINITY,&norm));
267*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n",i,j,(double)rdot,(double)norm));
268c4762a1bSJed Brown         }
269c4762a1bSJed Brown       }
270c4762a1bSJed Brown     }
271*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"    max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max));
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   case 1:
274c4762a1bSJed Brown     norm_max = 0.0;
275c4762a1bSJed Brown     for (i = il; i< iu; i++) {
276*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A, evec[i], vt1));
277*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(evec[i], vt2));
278c4762a1bSJed Brown       tmp  = -eval[i];
279*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(vt1,tmp,vt2));
280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(vt1, NORM_INFINITY, &norm));
281c4762a1bSJed Brown       norm = PetscAbs(norm);
282c4762a1bSJed Brown       if (norm > norm_max) norm_max = norm;
283c4762a1bSJed Brown       /* sniff, and bark if necessary */
284c4762a1bSJed Brown       if (norm > tols[0]) {
285*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  residual violation: %" PetscInt_FMT ", resi: %g\n",i, norm));
286c4762a1bSJed Brown       }
287c4762a1bSJed Brown     }
288*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"    max_resi:                    %g\n", (double)norm_max));
289c4762a1bSJed Brown     break;
290c4762a1bSJed Brown   default:
291*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%" PetscInt_FMT " is not supported \n",cklvl));
292c4762a1bSJed Brown   }
293*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vt2));
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vt1));
295c4762a1bSJed Brown   PetscFunctionReturn(0);
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown /*TEST
299c4762a1bSJed Brown 
300c4762a1bSJed Brown    build:
301c4762a1bSJed Brown       requires: complex
302c4762a1bSJed Brown 
303c4762a1bSJed Brown    test:
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    test:
306c4762a1bSJed Brown       suffix: 2
307c4762a1bSJed Brown       args: -test_zheevx
308c4762a1bSJed Brown 
309c4762a1bSJed Brown    test:
310c4762a1bSJed Brown       suffix: 3
311c4762a1bSJed Brown       args: -test_zhegv
312c4762a1bSJed Brown 
313c4762a1bSJed Brown    test:
314c4762a1bSJed Brown       suffix: 4
315c4762a1bSJed Brown       args: -test_zhegvx
316c4762a1bSJed Brown 
317c4762a1bSJed Brown TEST*/
318