1c4762a1bSJed Brown static char help[] = "Test LAPACK routine DSTEBZ() and DTEIN(). \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petscblaslapack.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown extern PetscErrorCode CkEigenSolutions(PetscInt,Mat,PetscInt,PetscInt,PetscScalar*,Vec*,PetscReal*); 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **args) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) || defined(PETSC_MISSING_LAPACK_STEBZ) || defined(PETSC_MISSING_LAPACK_STEIN) 11*327415f7SBarry Smith PetscFunctionBeginUser; 129566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 13c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP_SYS,"This example requires LAPACK routines dstebz and stien and real numbers"); 14c4762a1bSJed Brown #else 15c4762a1bSJed Brown PetscReal *work,tols[2]; 16c4762a1bSJed Brown PetscInt i,j; 17c4762a1bSJed Brown PetscBLASInt n,il=1,iu=5,*iblock,*isplit,*iwork,nevs,*ifail,cklvl=2; 18c4762a1bSJed Brown PetscMPIInt size; 19c4762a1bSJed Brown PetscBool flg; 20c4762a1bSJed Brown Vec *evecs; 21c4762a1bSJed Brown PetscScalar *evecs_array,*D,*E,*evals; 22c4762a1bSJed Brown Mat T; 23c4762a1bSJed Brown PetscReal vl=0.0,vu=4.0,tol= 1000*PETSC_MACHINE_EPSILON; 24c4762a1bSJed Brown PetscBLASInt nsplit,info; 25c4762a1bSJed Brown 26*327415f7SBarry Smith PetscFunctionBeginUser; 279566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 29be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 30c4762a1bSJed Brown 31c4762a1bSJed Brown n = 100; 32c4762a1bSJed Brown nevs = iu - il; 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*n+1,&D)); 34c4762a1bSJed Brown E = D + n; 35c4762a1bSJed Brown evals = E + n; 369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5*n+1,&work)); 379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*n+1,&iwork)); 389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*n+1,&iblock)); 39c4762a1bSJed Brown isplit = iblock + n; 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Set symmetric tridiagonal matrix */ 42c4762a1bSJed Brown for (i=0; i<n; i++) { 43c4762a1bSJed Brown D[i] = 2.0; 44c4762a1bSJed Brown E[i] = 1.0; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* Solve eigenvalue problem: A*evec = eval*evec */ 489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," LAPACKstebz_: compute %d eigenvalues...\n",nevs)); 49c4762a1bSJed Brown LAPACKstebz_("I","E",&n,&vl,&vu,&il,&iu,&tol,(PetscReal*)D,(PetscReal*)E,&nevs,&nsplit,(PetscReal*)evals,iblock,isplit,work,iwork,&info); 5028b400f6SJacob Faibussowitsch PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstebz_ fails. info %d",info); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," LAPACKstein_: compute %d found eigenvectors...\n",nevs)); 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n*nevs,&evecs_array)); 549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevs,&ifail)); 55c4762a1bSJed Brown LAPACKstein_(&n,(PetscReal*)D,(PetscReal*)E,&nevs,(PetscReal*)evals,iblock,isplit,evecs_array,&n,work,iwork,ifail,&info); 5628b400f6SJacob Faibussowitsch PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstein_ fails. info %d",info); 57c4762a1bSJed Brown /* View evals */ 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg)); 59c4762a1bSJed Brown if (flg) { 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," %d evals: \n",nevs)); 619566063dSJacob Faibussowitsch for (i=0; i<nevs; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %g\n",i,(double)evals[i])); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Check residuals and orthogonality */ 659566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&T)); 669566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n)); 679566063dSJacob Faibussowitsch PetscCall(MatSetType(T,MATSBAIJ)); 689566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(T)); 699566063dSJacob Faibussowitsch PetscCall(MatSetUp(T)); 70c4762a1bSJed Brown for (i=0; i<n; i++) { 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(T,1,&i,1,&i,&D[i],INSERT_VALUES)); 72c4762a1bSJed Brown if (i != n-1) { 73c4762a1bSJed Brown j = i+1; 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(T,1,&i,1,&j,&E[i],INSERT_VALUES)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown } 779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevs+1,&evecs)); 81c4762a1bSJed Brown for (i=0; i<nevs; i++) { 829566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&evecs[i])); 839566063dSJacob Faibussowitsch PetscCall(VecSetSizes(evecs[i],PETSC_DECIDE,n)); 849566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(evecs[i])); 859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(evecs[i],evecs_array+i*n)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown tols[0] = 1.e-8; tols[1] = 1.e-8; 899566063dSJacob Faibussowitsch PetscCall(CkEigenSolutions(cklvl,T,il-1,iu-1,evals,evecs,tols)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown for (i=0; i<nevs; i++) { 929566063dSJacob Faibussowitsch PetscCall(VecResetArray(evecs[i])); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* free space */ 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch for (i=0; i<nevs; i++) PetscCall(VecDestroy(&evecs[i])); 1009566063dSJacob Faibussowitsch PetscCall(PetscFree(evecs)); 1019566063dSJacob Faibussowitsch PetscCall(PetscFree(D)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 1039566063dSJacob Faibussowitsch PetscCall(PetscFree(iwork)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(iblock)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(evecs_array)); 1069566063dSJacob Faibussowitsch PetscCall(PetscFree(ifail)); 1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 108b122ec5aSJacob Faibussowitsch return 0; 109c4762a1bSJed Brown #endif 110c4762a1bSJed Brown } 111c4762a1bSJed Brown /*------------------------------------------------ 112c4762a1bSJed Brown Check the accuracy of the eigen solution 113c4762a1bSJed Brown ----------------------------------------------- */ 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown input: 116c4762a1bSJed Brown cklvl - check level: 117c4762a1bSJed Brown 1: check residual 118c4762a1bSJed Brown 2: 1 and check B-orthogonality locally 119c4762a1bSJed Brown A - matrix 120c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues 121c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process 122c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] || 123c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij 124c4762a1bSJed Brown */ 125c4762a1bSJed Brown #undef DEBUG_CkEigenSolutions 126c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscScalar *eval,Vec *evec,PetscReal *tols) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown PetscInt ierr,i,j,nev; 129c4762a1bSJed Brown Vec vt1,vt2; /* tmp vectors */ 130c4762a1bSJed Brown PetscReal norm,norm_max; 131c4762a1bSJed Brown PetscScalar dot,tmp; 132c4762a1bSJed Brown PetscReal dot_max; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBegin; 135c4762a1bSJed Brown nev = iu - il; 136c4762a1bSJed Brown if (nev <= 0) PetscFunctionReturn(0); 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0],&vt1)); 1399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0],&vt2)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown switch (cklvl) { 142c4762a1bSJed Brown case 2: 143c4762a1bSJed Brown dot_max = 0.0; 144c4762a1bSJed Brown for (i = il; i<iu; i++) { 1459566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt1)); 146c4762a1bSJed Brown for (j=il; j<iu; j++) { 1479566063dSJacob Faibussowitsch PetscCall(VecDot(evec[j],vt1,&dot)); 148c4762a1bSJed Brown if (j == i) { 149c4762a1bSJed Brown dot = PetscAbsScalar(dot - (PetscScalar)1.0); 150c4762a1bSJed Brown } else { 151c4762a1bSJed Brown dot = PetscAbsScalar(dot); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown if (PetscAbsScalar(dot) > dot_max) dot_max = PetscAbsScalar(dot); 154c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 155c4762a1bSJed Brown if (dot > tols[1]) { 1569566063dSJacob Faibussowitsch PetscCall(VecNorm(evec[i],NORM_INFINITY,&norm)); 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %g, norm: %d\n",i,j,(double)dot,(double)norm)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown #endif 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 1629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown case 1: 165c4762a1bSJed Brown norm_max = 0.0; 166c4762a1bSJed Brown for (i = il; i< iu; i++) { 1679566063dSJacob Faibussowitsch PetscCall(MatMult(A, evec[i], vt1)); 1689566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt2)); 169c4762a1bSJed Brown tmp = -eval[i]; 1709566063dSJacob Faibussowitsch PetscCall(VecAXPY(vt1,tmp,vt2)); 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(vt1, NORM_INFINITY, &norm)); 172c4762a1bSJed Brown norm = PetscAbsReal(norm); 173c4762a1bSJed Brown if (norm > norm_max) norm_max = norm; 174c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 175c4762a1bSJed Brown if (norm > tols[0]) { 1769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," residual violation: %d, resi: %g\n",i, norm)); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown #endif 179c4762a1bSJed Brown } 1809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max)); 181c4762a1bSJed Brown break; 182c4762a1bSJed Brown default: 1839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl)); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt2)); 1879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt1)); 188c4762a1bSJed Brown PetscFunctionReturn(0); 189c4762a1bSJed Brown } 190