1*c4762a1bSJed Brown static char help[] = "Test LAPACK routine DSTEBZ() and DTEIN(). \n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown #include <petscblaslapack.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown extern PetscErrorCode CkEigenSolutions(PetscInt,Mat,PetscInt,PetscInt,PetscScalar*,Vec*,PetscReal*); 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown int main(int argc,char **args) 9*c4762a1bSJed Brown { 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) || defined(PETSC_MISSING_LAPACK_STEBZ) || defined(PETSC_MISSING_LAPACK_STEIN) 12*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 13*c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP_SYS,"This example requires LAPACK routines dstebz and stien and real numbers"); 14*c4762a1bSJed Brown #else 15*c4762a1bSJed Brown PetscReal *work,tols[2]; 16*c4762a1bSJed Brown PetscInt i,j; 17*c4762a1bSJed Brown PetscBLASInt n,il=1,iu=5,*iblock,*isplit,*iwork,nevs,*ifail,cklvl=2; 18*c4762a1bSJed Brown PetscMPIInt size; 19*c4762a1bSJed Brown PetscBool flg; 20*c4762a1bSJed Brown Vec *evecs; 21*c4762a1bSJed Brown PetscScalar *evecs_array,*D,*E,*evals; 22*c4762a1bSJed Brown Mat T; 23*c4762a1bSJed Brown PetscReal vl=0.0,vu=4.0,tol= 1000*PETSC_MACHINE_EPSILON; 24*c4762a1bSJed Brown PetscBLASInt nsplit,info; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 27*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 28*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown n = 100; 31*c4762a1bSJed Brown nevs = iu - il; 32*c4762a1bSJed Brown ierr = PetscMalloc1(3*n+1,&D);CHKERRQ(ierr); 33*c4762a1bSJed Brown E = D + n; 34*c4762a1bSJed Brown evals = E + n; 35*c4762a1bSJed Brown ierr = PetscMalloc1(5*n+1,&work);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = PetscMalloc1(3*n+1,&iwork);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = PetscMalloc1(3*n+1,&iblock);CHKERRQ(ierr); 38*c4762a1bSJed Brown isplit = iblock + n; 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown /* Set symmetric tridiagonal matrix */ 41*c4762a1bSJed Brown for (i=0; i<n; i++) { 42*c4762a1bSJed Brown D[i] = 2.0; 43*c4762a1bSJed Brown E[i] = 1.0; 44*c4762a1bSJed Brown } 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* Solve eigenvalue problem: A*evec = eval*evec */ 47*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," LAPACKstebz_: compute %d eigenvalues...\n",nevs);CHKERRQ(ierr); 48*c4762a1bSJed Brown LAPACKstebz_("I","E",&n,&vl,&vu,&il,&iu,&tol,(PetscReal*)D,(PetscReal*)E,&nevs,&nsplit,(PetscReal*)evals,iblock,isplit,work,iwork,&info); 49*c4762a1bSJed Brown if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstebz_ fails. info %d",info); 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," LAPACKstein_: compute %d found eigenvectors...\n",nevs);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = PetscMalloc1(n*nevs,&evecs_array);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = PetscMalloc1(nevs,&ifail);CHKERRQ(ierr); 54*c4762a1bSJed Brown LAPACKstein_(&n,(PetscReal*)D,(PetscReal*)E,&nevs,(PetscReal*)evals,iblock,isplit,evecs_array,&n,work,iwork,ifail,&info); 55*c4762a1bSJed Brown if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"LAPACKstein_ fails. info %d",info); 56*c4762a1bSJed Brown /* View evals */ 57*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-eig_view", &flg);CHKERRQ(ierr); 58*c4762a1bSJed Brown if (flg) { 59*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," %d evals: \n",nevs);CHKERRQ(ierr); 60*c4762a1bSJed Brown for (i=0; i<nevs; i++) {ierr = PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i,(double)evals[i]);CHKERRQ(ierr);} 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown /* Check residuals and orthogonality */ 64*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatSetType(T,MATSBAIJ);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatSetFromOptions(T);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatSetUp(T);CHKERRQ(ierr); 69*c4762a1bSJed Brown for (i=0; i<n; i++) { 70*c4762a1bSJed Brown ierr = MatSetValues(T,1,&i,1,&i,&D[i],INSERT_VALUES);CHKERRQ(ierr); 71*c4762a1bSJed Brown if (i != n-1) { 72*c4762a1bSJed Brown j = i+1; 73*c4762a1bSJed Brown ierr = MatSetValues(T,1,&i,1,&j,&E[i],INSERT_VALUES);CHKERRQ(ierr); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown ierr = PetscMalloc1(nevs+1,&evecs);CHKERRQ(ierr); 80*c4762a1bSJed Brown for (i=0; i<nevs; i++) { 81*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_SELF,&evecs[i]);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecSetSizes(evecs[i],PETSC_DECIDE,n);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = VecSetFromOptions(evecs[i]);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = VecPlaceArray(evecs[i],evecs_array+i*n);CHKERRQ(ierr); 85*c4762a1bSJed Brown } 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown tols[0] = 1.e-8; tols[1] = 1.e-8; 88*c4762a1bSJed Brown ierr = CkEigenSolutions(cklvl,T,il-1,iu-1,evals,evecs,tols);CHKERRQ(ierr); 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown for (i=0; i<nevs; i++) { 91*c4762a1bSJed Brown ierr = VecResetArray(evecs[i]);CHKERRQ(ierr); 92*c4762a1bSJed Brown } 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown /* free space */ 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown ierr = MatDestroy(&T);CHKERRQ(ierr); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown for (i=0; i<nevs; i++) { ierr = VecDestroy(&evecs[i]);CHKERRQ(ierr);} 99*c4762a1bSJed Brown ierr = PetscFree(evecs);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = PetscFree(D);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = PetscFree(work);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscFree(iwork);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = PetscFree(iblock);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscFree(evecs_array);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = PetscFree(ifail);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = PetscFinalize(); 107*c4762a1bSJed Brown return ierr; 108*c4762a1bSJed Brown #endif 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown /*------------------------------------------------ 111*c4762a1bSJed Brown Check the accuracy of the eigen solution 112*c4762a1bSJed Brown ----------------------------------------------- */ 113*c4762a1bSJed Brown /* 114*c4762a1bSJed Brown input: 115*c4762a1bSJed Brown cklvl - check level: 116*c4762a1bSJed Brown 1: check residual 117*c4762a1bSJed Brown 2: 1 and check B-orthogonality locally 118*c4762a1bSJed Brown A - matrix 119*c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues 120*c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process 121*c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] || 122*c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij 123*c4762a1bSJed Brown */ 124*c4762a1bSJed Brown #undef DEBUG_CkEigenSolutions 125*c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscScalar *eval,Vec *evec,PetscReal *tols) 126*c4762a1bSJed Brown { 127*c4762a1bSJed Brown PetscInt ierr,i,j,nev; 128*c4762a1bSJed Brown Vec vt1,vt2; /* tmp vectors */ 129*c4762a1bSJed Brown PetscReal norm,norm_max; 130*c4762a1bSJed Brown PetscScalar dot,tmp; 131*c4762a1bSJed Brown PetscReal dot_max; 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown PetscFunctionBegin; 134*c4762a1bSJed Brown nev = iu - il; 135*c4762a1bSJed Brown if (nev <= 0) PetscFunctionReturn(0); 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown ierr = VecDuplicate(evec[0],&vt1);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = VecDuplicate(evec[0],&vt2);CHKERRQ(ierr); 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown switch (cklvl) { 141*c4762a1bSJed Brown case 2: 142*c4762a1bSJed Brown dot_max = 0.0; 143*c4762a1bSJed Brown for (i = il; i<iu; i++) { 144*c4762a1bSJed Brown ierr = VecCopy(evec[i], vt1);CHKERRQ(ierr); 145*c4762a1bSJed Brown for (j=il; j<iu; j++) { 146*c4762a1bSJed Brown ierr = VecDot(evec[j],vt1,&dot);CHKERRQ(ierr); 147*c4762a1bSJed Brown if (j == i) { 148*c4762a1bSJed Brown dot = PetscAbsScalar(dot - (PetscScalar)1.0); 149*c4762a1bSJed Brown } else { 150*c4762a1bSJed Brown dot = PetscAbsScalar(dot); 151*c4762a1bSJed Brown } 152*c4762a1bSJed Brown if (PetscAbsScalar(dot) > dot_max) dot_max = PetscAbsScalar(dot); 153*c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 154*c4762a1bSJed Brown if (dot > tols[1]) { 155*c4762a1bSJed Brown ierr = VecNorm(evec[i],NORM_INFINITY,&norm);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"|delta(%d,%d)|: %g, norm: %d\n",i,j,(double)dot,(double)norm);CHKERRQ(ierr); 157*c4762a1bSJed Brown } 158*c4762a1bSJed Brown #endif 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max);CHKERRQ(ierr); 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown case 1: 164*c4762a1bSJed Brown norm_max = 0.0; 165*c4762a1bSJed Brown for (i = il; i< iu; i++) { 166*c4762a1bSJed Brown ierr = MatMult(A, evec[i], vt1);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = VecCopy(evec[i], vt2);CHKERRQ(ierr); 168*c4762a1bSJed Brown tmp = -eval[i]; 169*c4762a1bSJed Brown ierr = VecAXPY(vt1,tmp,vt2);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = VecNorm(vt1, NORM_INFINITY, &norm);CHKERRQ(ierr); 171*c4762a1bSJed Brown norm = PetscAbsReal(norm); 172*c4762a1bSJed Brown if (norm > norm_max) norm_max = norm; 173*c4762a1bSJed Brown #if defined(DEBUG_CkEigenSolutions) 174*c4762a1bSJed Brown if (norm > tols[0]) { 175*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," residual violation: %d, resi: %g\n",i, norm);CHKERRQ(ierr); 176*c4762a1bSJed Brown } 177*c4762a1bSJed Brown #endif 178*c4762a1bSJed Brown } 179*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max);CHKERRQ(ierr); 180*c4762a1bSJed Brown break; 181*c4762a1bSJed Brown default: 182*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%d is not supported \n",cklvl);CHKERRQ(ierr); 183*c4762a1bSJed Brown } 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown ierr = VecDestroy(&vt2);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = VecDestroy(&vt1);CHKERRQ(ierr); 187*c4762a1bSJed Brown PetscFunctionReturn(0); 188*c4762a1bSJed Brown } 189