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 27*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 285f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 295f80ce2aSJacob Faibussowitsch PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 30c4762a1bSJed Brown 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zheevx", &flg)); 32c4762a1bSJed Brown if (flg) { 33c4762a1bSJed Brown TestZHEEV = PETSC_FALSE; 34c4762a1bSJed Brown TestZHEEVX = PETSC_TRUE; 35c4762a1bSJed Brown } 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zhegv", &flg)); 37c4762a1bSJed Brown if (flg) { 38c4762a1bSJed Brown TestZHEEV = PETSC_FALSE; 39c4762a1bSJed Brown TestZHEGV = PETSC_TRUE; 40c4762a1bSJed Brown } 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_zhegvx", &flg)); 42c4762a1bSJed Brown if (flg) { 43c4762a1bSJed Brown TestZHEEV = PETSC_FALSE; 44c4762a1bSJed Brown TestZHEGVX = PETSC_TRUE; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 49c4762a1bSJed Brown dim = n*n; 50c4762a1bSJed Brown 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQDENSE)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 56c4762a1bSJed Brown 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-norandom",&flg)); 58c4762a1bSJed Brown if (flg) use_random = 0; 59c4762a1bSJed Brown else use_random = 1; 60c4762a1bSJed Brown if (use_random) { 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rctx)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 715f80ce2aSJacob Faibussowitsch J = Ii-n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown if (i<n-1) { 745f80ce2aSJacob Faibussowitsch J = Ii+n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown if (j>0) { 775f80ce2aSJacob Faibussowitsch J = Ii-1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown if (j<n-1) { 805f80ce2aSJacob Faibussowitsch J = Ii+1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 81c4762a1bSJed Brown } 825f80ce2aSJacob Faibussowitsch if (use_random) CHKERRQ(PetscRandomGetValue(rctx,&sigma2)); 83c4762a1bSJed Brown v = 4.0 - sigma1*h2; 845f80ce2aSJacob Faibussowitsch CHKERRQ(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; 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 90c4762a1bSJed Brown v = -sigma2*h2; 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES)); 925f80ce2aSJacob Faibussowitsch if (use_random) CHKERRQ(PetscRandomDestroy(&rctx)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 95c4762a1bSJed Brown m = n = dim; 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Check whether A is symmetric */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg)); 99c4762a1bSJed Brown if (flg) { 100c4762a1bSJed Brown Mat Trans; 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Trans)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(A, Trans, &isSymmetric)); 1035f80ce2aSJacob Faibussowitsch PetscCheck(isSymmetric,PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric"); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Trans)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg)); 109c4762a1bSJed Brown if (flg) { 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A_dense)); 111c4762a1bSJed Brown } else { 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&B)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,dim,dim)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSEQDENSE)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 120c4762a1bSJed Brown v = 1.0; 121c4762a1bSJed Brown for (Ii=0; Ii<dim; Ii++) { 1225f80ce2aSJacob Faibussowitsch CHKERRQ(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 /*===================================================*/ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast(2*n,&lwork)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast(n,&bn)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&evals)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(lwork,&work)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(A_dense,&arrayA)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown if (TestZHEEV) { /* test zheev() */ 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n",m)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n-2,&rwork)); 136c4762a1bSJed Brown LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,rwork,&lierr); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast((0.2*m),&iu)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*n+1,&evecs_array)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(7*n+1,&rwork)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5*n+1,&iwork)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1545f80ce2aSJacob Faibussowitsch CHKERRQ(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); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iwork)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ifail)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rwork)); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown if (TestZHEGV) { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute all %" PetscInt_FMT " eigensolutions...\n",m)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n+1,&rwork)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(B,&arrayB)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rwork)); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown if (TestZHEGVX) { 172c4762a1bSJed Brown il = 1; 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast((0.2*m),&iu)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," LAPACKsygv: compute %d to %d-th eigensolutions...\n",il,iu)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*n+1,&evecs_array)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(6*n+1,&iwork)); 177c4762a1bSJed Brown ifail = iwork + 5*n; 1785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(7*n+1,&rwork)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(B,&arrayB)); 180c4762a1bSJed Brown vl = 0.0; vu = 8.0; 1815f80ce2aSJacob Faibussowitsch CHKERRQ(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); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(B,&arrayB)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iwork)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rwork)); 186c4762a1bSJed Brown } 1875f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg)); 192c4762a1bSJed Brown if (flg) { 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %d evals: \n",nevs)); 1945f80ce2aSJacob Faibussowitsch for (i=0; i<nevs; i++) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT " %g\n",i+il,(double)evals[i])); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Check residuals and orthogonality */ 1985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nevs+1,&evecs)); 199c4762a1bSJed Brown for (i=0; i<nevs; i++) { 2005f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&evecs[i])); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(evecs[i],PETSC_DECIDE,n)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(evecs[i])); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(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; 2075f80ce2aSJacob Faibussowitsch CHKERRQ(CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols)); 2085f80ce2aSJacob Faibussowitsch for (i=0; i<nevs; i++) CHKERRQ(VecDestroy(&evecs[i])); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(evecs)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* Free work space. */ 212c4762a1bSJed Brown if (TestZHEEVX || TestZHEGVX) { 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(evecs_array)); 214c4762a1bSJed Brown } 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(evals)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(work)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A_dense)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 220*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 221*b122ec5aSJacob 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 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(evec[0],&vt1)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 2555f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(evec[i], vt1)); 256c4762a1bSJed Brown for (j=il; j<iu; j++) { 2575f80ce2aSJacob Faibussowitsch CHKERRQ(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]) { 2655f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(evec[i],NORM_INFINITY,&norm)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 2705f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 2755f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, evec[i], vt1)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(evec[i], vt2)); 277c4762a1bSJed Brown tmp = -eval[i]; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(vt1,tmp,vt2)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(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]) { 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," residual violation: %" PetscInt_FMT ", resi: %g\n",i, norm)); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 2875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max)); 288c4762a1bSJed Brown break; 289c4762a1bSJed Brown default: 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%" PetscInt_FMT " is not supported \n",cklvl)); 291c4762a1bSJed Brown } 2925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vt2)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(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