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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 12c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP_SYS,"This example requires LAPACK routines dstebz and stien and real numbers"); 13c4762a1bSJed Brown #else 14c4762a1bSJed Brown PetscReal *work,tols[2]; 15c4762a1bSJed Brown PetscInt i,j; 16c4762a1bSJed Brown PetscBLASInt n,il=1,iu=5,*iblock,*isplit,*iwork,nevs,*ifail,cklvl=2; 17c4762a1bSJed Brown PetscMPIInt size; 18c4762a1bSJed Brown PetscBool flg; 19c4762a1bSJed Brown Vec *evecs; 20c4762a1bSJed Brown PetscScalar *evecs_array,*D,*E,*evals; 21c4762a1bSJed Brown Mat T; 22c4762a1bSJed Brown PetscReal vl=0.0,vu=4.0,tol= 1000*PETSC_MACHINE_EPSILON; 23c4762a1bSJed Brown PetscBLASInt nsplit,info; 24c4762a1bSJed Brown 25*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 265f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 272c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 28c4762a1bSJed Brown 29c4762a1bSJed Brown n = 100; 30c4762a1bSJed Brown nevs = iu - il; 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n+1,&D)); 32c4762a1bSJed Brown E = D + n; 33c4762a1bSJed Brown evals = E + n; 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5*n+1,&work)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n+1,&iwork)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n+1,&iblock)); 37c4762a1bSJed Brown isplit = iblock + n; 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Set symmetric tridiagonal matrix */ 40c4762a1bSJed Brown for (i=0; i<n; i++) { 41c4762a1bSJed Brown D[i] = 2.0; 42c4762a1bSJed Brown E[i] = 1.0; 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Solve eigenvalue problem: A*evec = eval*evec */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," LAPACKstebz_: compute %d eigenvalues...\n",nevs)); 47c4762a1bSJed Brown LAPACKstebz_("I","E",&n,&vl,&vu,&il,&iu,&tol,(PetscReal*)D,(PetscReal*)E,&nevs,&nsplit,(PetscReal*)evals,iblock,isplit,work,iwork,&info); 4828b400f6SJacob Faibussowitsch PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstebz_ fails. info %d",info); 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," LAPACKstein_: compute %d found eigenvectors...\n",nevs)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n*nevs,&evecs_array)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nevs,&ifail)); 53c4762a1bSJed Brown LAPACKstein_(&n,(PetscReal*)D,(PetscReal*)E,&nevs,(PetscReal*)evals,iblock,isplit,evecs_array,&n,work,iwork,ifail,&info); 5428b400f6SJacob Faibussowitsch PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstein_ fails. info %d",info); 55c4762a1bSJed Brown /* View evals */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg)); 57c4762a1bSJed Brown if (flg) { 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %d evals: \n",nevs)); 595f80ce2aSJacob Faibussowitsch for (i=0; i<nevs; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %g\n",i,(double)evals[i])); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Check residuals and orthogonality */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&T)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(T,MATSBAIJ)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(T)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(T)); 68c4762a1bSJed Brown for (i=0; i<n; i++) { 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(T,1,&i,1,&i,&D[i],INSERT_VALUES)); 70c4762a1bSJed Brown if (i != n-1) { 71c4762a1bSJed Brown j = i+1; 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(T,1,&i,1,&j,&E[i],INSERT_VALUES)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown } 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 77c4762a1bSJed Brown 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nevs+1,&evecs)); 79c4762a1bSJed Brown for (i=0; i<nevs; i++) { 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&evecs[i])); 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(evecs[i],PETSC_DECIDE,n)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(evecs[i])); 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(evecs[i],evecs_array+i*n)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown tols[0] = 1.e-8; tols[1] = 1.e-8; 875f80ce2aSJacob Faibussowitsch CHKERRQ(CkEigenSolutions(cklvl,T,il-1,iu-1,evals,evecs,tols)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown for (i=0; i<nevs; i++) { 905f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(evecs[i])); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* free space */ 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&T)); 96c4762a1bSJed Brown 975f80ce2aSJacob Faibussowitsch for (i=0; i<nevs; i++) CHKERRQ(VecDestroy(&evecs[i])); 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(evecs)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(D)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(work)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iwork)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iblock)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(evecs_array)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ifail)); 105*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 106*b122ec5aSJacob Faibussowitsch return 0; 107c4762a1bSJed Brown #endif 108c4762a1bSJed Brown } 109c4762a1bSJed Brown /*------------------------------------------------ 110c4762a1bSJed Brown Check the accuracy of the eigen solution 111c4762a1bSJed Brown ----------------------------------------------- */ 112c4762a1bSJed Brown /* 113c4762a1bSJed Brown input: 114c4762a1bSJed Brown cklvl - check level: 115c4762a1bSJed Brown 1: check residual 116c4762a1bSJed Brown 2: 1 and check B-orthogonality locally 117c4762a1bSJed Brown A - matrix 118c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues 119c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process 120c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] || 121c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown #undef DEBUG_CkEigenSolutions 124c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscScalar *eval,Vec *evec,PetscReal *tols) 125c4762a1bSJed Brown { 126c4762a1bSJed Brown PetscInt ierr,i,j,nev; 127c4762a1bSJed Brown Vec vt1,vt2; /* tmp vectors */ 128c4762a1bSJed Brown PetscReal norm,norm_max; 129c4762a1bSJed Brown PetscScalar dot,tmp; 130c4762a1bSJed Brown PetscReal dot_max; 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscFunctionBegin; 133c4762a1bSJed Brown nev = iu - il; 134c4762a1bSJed Brown if (nev <= 0) PetscFunctionReturn(0); 135c4762a1bSJed Brown 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(evec[0],&vt1)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(evec[0],&vt2)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown switch (cklvl) { 140c4762a1bSJed Brown case 2: 141c4762a1bSJed Brown dot_max = 0.0; 142c4762a1bSJed Brown for (i = il; i<iu; i++) { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(evec[i], vt1)); 144c4762a1bSJed Brown for (j=il; j<iu; j++) { 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(evec[j],vt1,&dot)); 146c4762a1bSJed Brown if (j == i) { 147c4762a1bSJed Brown dot = PetscAbsScalar(dot - (PetscScalar)1.0); 148c4762a1bSJed Brown } else { 149c4762a1bSJed Brown dot = PetscAbsScalar(dot); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown if (PetscAbsScalar(dot) > dot_max) dot_max = PetscAbsScalar(dot); 152c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 153c4762a1bSJed Brown if (dot > tols[1]) { 1545f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(evec[i],NORM_INFINITY,&norm)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %g, norm: %d\n",i,j,(double)dot,(double)norm)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown #endif 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 1605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown case 1: 163c4762a1bSJed Brown norm_max = 0.0; 164c4762a1bSJed Brown for (i = il; i< iu; i++) { 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, evec[i], vt1)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(evec[i], vt2)); 167c4762a1bSJed Brown tmp = -eval[i]; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(vt1,tmp,vt2)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(vt1, NORM_INFINITY, &norm)); 170c4762a1bSJed Brown norm = PetscAbsReal(norm); 171c4762a1bSJed Brown if (norm > norm_max) norm_max = norm; 172c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 173c4762a1bSJed Brown if (norm > tols[0]) { 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," residual violation: %d, resi: %g\n",i, norm)); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown #endif 177c4762a1bSJed Brown } 1785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max)); 179c4762a1bSJed Brown break; 180c4762a1bSJed Brown default: 1815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl)); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 1845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vt2)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vt1)); 186c4762a1bSJed Brown PetscFunctionReturn(0); 187c4762a1bSJed Brown } 188