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