static char help[] = "Test LAPACK routine ZHEEV, ZHEEVX, ZHEGV and ZHEGVX. \n\ ZHEEV computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. \n\n"; #include #include extern PetscErrorCode CkEigenSolutions(PetscInt,Mat,PetscInt,PetscInt,PetscReal*,Vec*,PetscReal*); int main(int argc,char **args) { Mat A,A_dense,B; Vec *evecs; PetscBool flg,TestZHEEV=PETSC_TRUE,TestZHEEVX=PETSC_FALSE,TestZHEGV=PETSC_FALSE,TestZHEGVX=PETSC_FALSE; PetscBool isSymmetric; PetscScalar *arrayA,*arrayB,*evecs_array=NULL,*work; PetscReal *evals,*rwork; PetscMPIInt size; PetscInt m,i,j,cklvl=2; PetscReal vl,vu,abstol=1.e-8; PetscBLASInt nn,nevs,il,iu,*iwork,*ifail,lwork,lierr,bn,one=1; PetscReal tols[2]; PetscScalar v,sigma2; PetscRandom rctx; PetscReal h2,sigma1 = 100.0; PetscInt dim,Ii,J,n = 6,use_random; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); PetscCall(PetscOptionsHasName(NULL,NULL, "-test_zheevx", &flg)); if (flg) { TestZHEEV = PETSC_FALSE; TestZHEEVX = PETSC_TRUE; } PetscCall(PetscOptionsHasName(NULL,NULL, "-test_zhegv", &flg)); if (flg) { TestZHEEV = PETSC_FALSE; TestZHEGV = PETSC_TRUE; } PetscCall(PetscOptionsHasName(NULL,NULL, "-test_zhegvx", &flg)); if (flg) { TestZHEEV = PETSC_FALSE; TestZHEGVX = PETSC_TRUE; } PetscCall(PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); dim = n*n; PetscCall(MatCreate(PETSC_COMM_SELF,&A)); PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim)); PetscCall(MatSetType(A,MATSEQDENSE)); PetscCall(MatSetFromOptions(A)); PetscCall(MatSetUp(A)); PetscCall(PetscOptionsHasName(NULL,NULL,"-norandom",&flg)); if (flg) use_random = 0; else use_random = 1; if (use_random) { PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rctx)); PetscCall(PetscRandomSetFromOptions(rctx)); PetscCall(PetscRandomSetInterval(rctx,0.0,PETSC_i)); } else { sigma2 = 10.0*PETSC_i; } h2 = 1.0/((n+1)*(n+1)); for (Ii=0; Ii0) { J = Ii-n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); } if (i0) { J = Ii-1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); } if (j 0,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%d, no eigensolution has found", nevs); /* View evals */ PetscCall(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg)); if (flg) { PetscCall(PetscPrintf(PETSC_COMM_WORLD," %d evals: \n",nevs)); for (i=0; i dot_max) dot_max = rdot; if (rdot > tols[1]) { PetscCall(VecNorm(evec[i],NORM_INFINITY,&norm)); PetscCall(PetscPrintf(PETSC_COMM_SELF,"|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n",i,j,(double)rdot,(double)norm)); } } } PetscCall(PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max)); case 1: norm_max = 0.0; for (i = il; i< iu; i++) { PetscCall(MatMult(A, evec[i], vt1)); PetscCall(VecCopy(evec[i], vt2)); tmp = -eval[i]; PetscCall(VecAXPY(vt1,tmp,vt2)); PetscCall(VecNorm(vt1, NORM_INFINITY, &norm)); norm = PetscAbs(norm); if (norm > norm_max) norm_max = norm; /* sniff, and bark if necessary */ if (norm > tols[0]) { PetscCall(PetscPrintf(PETSC_COMM_WORLD," residual violation: %" PetscInt_FMT ", resi: %g\n",i, norm)); } } PetscCall(PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max)); break; default: PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%" PetscInt_FMT " is not supported \n",cklvl)); } PetscCall(VecDestroy(&vt2)); PetscCall(VecDestroy(&vt1)); PetscFunctionReturn(0); } /*TEST build: requires: complex test: test: suffix: 2 args: -test_zheevx test: suffix: 3 args: -test_zhegv test: suffix: 4 args: -test_zhegvx TEST*/