1c4762a1bSJed Brown static char help[] = "Test LAPACK routine DSYEV() or DSYEVX(). \n\ 2c4762a1bSJed Brown Reads PETSc matrix A \n\ 3c4762a1bSJed Brown then computes selected eigenvalues, and optionally, eigenvectors of \n\ 4c4762a1bSJed Brown a real generalized symmetric-definite eigenproblem \n\ 5c4762a1bSJed Brown A*x = lambda*x \n\ 6c4762a1bSJed Brown Input parameters include\n\ 7c4762a1bSJed Brown -f <input_file> : file to load\n\ 8c4762a1bSJed Brown e.g. ./ex116 -f $DATAFILESPATH/matrices/small \n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown #include <petscblaslapack.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown extern PetscErrorCode CkEigenSolutions(PetscInt,Mat,PetscInt,PetscInt,PetscReal*,Vec*,PetscReal*); 14c4762a1bSJed Brown 15c4762a1bSJed Brown int main(int argc,char **args) 16c4762a1bSJed Brown { 17c4762a1bSJed Brown Mat A,A_dense; 18c4762a1bSJed Brown Vec *evecs; 19c4762a1bSJed Brown PetscViewer fd; /* viewer */ 20c4762a1bSJed Brown char file[1][PETSC_MAX_PATH_LEN]; /* input file name */ 21c4762a1bSJed Brown PetscBool flg,TestSYEVX=PETSC_TRUE; 22c4762a1bSJed Brown PetscErrorCode ierr; 23c4762a1bSJed Brown PetscBool isSymmetric; 24c4762a1bSJed Brown PetscScalar *arrayA,*evecs_array,*work,*evals; 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown PetscInt m,n,i,cklvl=2; 27c4762a1bSJed Brown PetscBLASInt nevs,il,iu,in; 28c4762a1bSJed Brown PetscReal vl,vu,abstol=1.e-8; 29c4762a1bSJed Brown PetscBLASInt *iwork,*ifail,lwork,lierr,bn; 30c4762a1bSJed Brown PetscReal tols[2]; 31c4762a1bSJed Brown 32c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 33ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 34c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-test_syev", &flg);CHKERRQ(ierr); 37c4762a1bSJed Brown if (flg) { 38c4762a1bSJed Brown TestSYEVX = PETSC_FALSE; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Determine files from which we read the two matrices */ 42589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file[0],sizeof(file[0]),&flg);CHKERRQ(ierr); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Load matrix A */ 45c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Check whether A is symmetric */ 53c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg);CHKERRQ(ierr); 54c4762a1bSJed Brown if (flg) { 55c4762a1bSJed Brown Mat Trans; 56c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX, &Trans);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatEqual(A, Trans, &isSymmetric);CHKERRQ(ierr); 58c4762a1bSJed Brown if (!isSymmetric) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric"); 59c4762a1bSJed Brown ierr = MatDestroy(&Trans);CHKERRQ(ierr); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Solve eigenvalue problem: A_dense*x = lambda*B*x */ 63c4762a1bSJed Brown /*==================================================*/ 64c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */ 65c4762a1bSJed Brown ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);CHKERRQ(ierr); 66c4762a1bSJed Brown 67c4762a1bSJed Brown ierr = PetscBLASIntCast(8*n,&lwork);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = PetscBLASIntCast(n,&bn);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = PetscMalloc1(n,&evals);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscMalloc1(lwork,&work);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatDenseGetArray(A_dense,&arrayA);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown if (!TestSYEVX) { /* test syev() */ 74c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," LAPACKsyev: compute all %D eigensolutions...\n",m);CHKERRQ(ierr); 75c4762a1bSJed Brown LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,&lierr); 76c4762a1bSJed Brown evecs_array = arrayA; 77c4762a1bSJed Brown ierr = PetscBLASIntCast(m,&nevs);CHKERRQ(ierr); 78c4762a1bSJed Brown il = 1; 79c4762a1bSJed Brown ierr = PetscBLASIntCast(m,&iu);CHKERRQ(ierr); 80c4762a1bSJed Brown } else { /* test syevx() */ 81c4762a1bSJed Brown il = 1; 82c4762a1bSJed Brown ierr = PetscBLASIntCast(0.2*m,&iu);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = PetscBLASIntCast(n,&in);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," LAPACKsyevx: compute %d to %d-th eigensolutions...\n",il,iu);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = PetscMalloc1(m*n+1,&evecs_array);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = PetscMalloc1(6*n+1,&iwork);CHKERRQ(ierr); 87c4762a1bSJed Brown ifail = iwork + 5*n; 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* in the case "I", vl and vu are not referenced */ 90c4762a1bSJed Brown vl = 0.0; vu = 8.0; 91c4762a1bSJed Brown LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&in,work,&lwork,iwork,ifail,&lierr); 92c4762a1bSJed Brown ierr = PetscFree(iwork);CHKERRQ(ierr); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown ierr = MatDenseRestoreArray(A_dense,&arrayA);CHKERRQ(ierr); 95c4762a1bSJed Brown if (nevs <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%D, no eigensolution has found", nevs); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* View eigenvalues */ 98c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-eig_view", &flg);CHKERRQ(ierr); 99c4762a1bSJed Brown if (flg) { 100c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," %D evals: \n",nevs);CHKERRQ(ierr); 101c4762a1bSJed Brown for (i=0; i<nevs; i++) {ierr = PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i+il,(double)evals[i]);CHKERRQ(ierr);} 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Check residuals and orthogonality */ 105c4762a1bSJed Brown ierr = PetscMalloc1(nevs+1,&evecs);CHKERRQ(ierr); 106c4762a1bSJed Brown for (i=0; i<nevs; i++) { 107c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_SELF,&evecs[i]);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = VecSetSizes(evecs[i],PETSC_DECIDE,n);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = VecSetFromOptions(evecs[i]);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = VecPlaceArray(evecs[i],evecs_array+i*n);CHKERRQ(ierr); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown tols[0] = tols[1] = PETSC_SQRT_MACHINE_EPSILON; 114c4762a1bSJed Brown ierr = CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Free work space. */ 117c4762a1bSJed Brown for (i=0; i<nevs; i++) { ierr = VecDestroy(&evecs[i]);CHKERRQ(ierr);} 118c4762a1bSJed Brown ierr = PetscFree(evecs);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDestroy(&A_dense);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscFree(work);CHKERRQ(ierr); 121c4762a1bSJed Brown if (TestSYEVX) {ierr = PetscFree(evecs_array);CHKERRQ(ierr);} 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Compute SVD: A_dense = U*SIGMA*transpose(V), 124c4762a1bSJed Brown JOBU=JOBV='S': the first min(m,n) columns of U and V are returned in the arrayU and arrayV; */ 125c4762a1bSJed Brown /*==============================================================================================*/ 126c4762a1bSJed Brown { 127c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */ 128c4762a1bSJed Brown PetscScalar *arrayU,*arrayVT,*arrayErr,alpha=1.0,beta=-1.0; 129c4762a1bSJed Brown Mat Err; 130c4762a1bSJed Brown PetscBLASInt minMN,maxMN,im,in; 131c4762a1bSJed Brown PetscInt j; 132c4762a1bSJed Brown PetscReal norm; 133c4762a1bSJed Brown 134c4762a1bSJed Brown ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense);CHKERRQ(ierr); 135c4762a1bSJed Brown 136c4762a1bSJed Brown minMN = PetscMin(m,n); 137c4762a1bSJed Brown maxMN = PetscMax(m,n); 138c4762a1bSJed Brown lwork = 5*minMN + maxMN; 139c4762a1bSJed Brown ierr = PetscMalloc4(m*minMN,&arrayU,m*minMN,&arrayVT,m*minMN,&arrayErr,lwork,&work);CHKERRQ(ierr); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Create matrix Err for checking error */ 142c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Err);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatSetSizes(Err,PETSC_DECIDE,PETSC_DECIDE,m,minMN);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = MatSetType(Err,MATSEQDENSE);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(Err,(PetscScalar*)arrayErr);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Save A to arrayErr for checking accuracy later. arrayA will be destroyed by LAPACKgesvd_() */ 148c4762a1bSJed Brown ierr = MatDenseGetArray(A_dense,&arrayA);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = PetscArraycpy(arrayErr,arrayA,m*minMN);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = PetscBLASIntCast(m,&im);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscBLASIntCast(n,&in);CHKERRQ(ierr); 153c4762a1bSJed Brown /* Compute A = U*SIGMA*VT */ 154c4762a1bSJed Brown LAPACKgesvd_("S","S",&im,&in,arrayA,&im,evals,arrayU,&minMN,arrayVT,&minMN,work,&lwork,&lierr); 155c4762a1bSJed Brown ierr = MatDenseRestoreArray(A_dense,&arrayA);CHKERRQ(ierr); 156c4762a1bSJed Brown if (!lierr) { 157c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," 1st 10 of %d singular values: \n",minMN);CHKERRQ(ierr); 158c4762a1bSJed Brown for (i=0; i<10; i++) {ierr = PetscPrintf(PETSC_COMM_SELF,"%D %g\n",i,(double)evals[i]);CHKERRQ(ierr);} 159c4762a1bSJed Brown } else { 160c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"LAPACKgesvd_ fails!");CHKERRQ(ierr); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* Check Err = (U*Sigma*V^T - A) using BLASgemm() */ 164c4762a1bSJed Brown /* U = U*Sigma */ 165c4762a1bSJed Brown for (j=0; j<minMN; j++) { /* U[:,j] = sigma[j]*U[:,j] */ 166c4762a1bSJed Brown for (i=0; i<m; i++) arrayU[j*m+i] *= evals[j]; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown /* Err = U*VT - A = alpha*U*VT + beta*Err */ 169c4762a1bSJed Brown BLASgemm_("N","N",&im,&minMN,&minMN,&alpha,arrayU,&im,arrayVT,&minMN,&beta,arrayErr,&im); 170c4762a1bSJed Brown ierr = MatNorm(Err,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," || U*Sigma*VT - A || = %g\n",(double)norm);CHKERRQ(ierr); 172c4762a1bSJed Brown 173c4762a1bSJed Brown ierr = PetscFree4(arrayU,arrayVT,arrayErr,work);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = PetscFree(evals);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = MatDestroy(&A_dense);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatDestroy(&Err);CHKERRQ(ierr); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = PetscFinalize(); 181c4762a1bSJed Brown return ierr; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown /*------------------------------------------------ 184c4762a1bSJed Brown Check the accuracy of the eigen solution 185c4762a1bSJed Brown ----------------------------------------------- */ 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown input: 188c4762a1bSJed Brown cklvl - check level: 189c4762a1bSJed Brown 1: check residual 190c4762a1bSJed Brown 2: 1 and check B-orthogonality locally 191c4762a1bSJed Brown A - matrix 192c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues 193c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process 194c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] || 195c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij 196c4762a1bSJed Brown */ 197c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown PetscInt ierr,i,j,nev; 200c4762a1bSJed Brown Vec vt1,vt2; /* tmp vectors */ 201c4762a1bSJed Brown PetscReal norm,tmp,dot,norm_max,dot_max; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBegin; 204c4762a1bSJed Brown nev = iu - il; 205c4762a1bSJed Brown if (nev <= 0) PetscFunctionReturn(0); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /*ierr = VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);*/ 208c4762a1bSJed Brown ierr = VecDuplicate(evec[0],&vt1);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = VecDuplicate(evec[0],&vt2);CHKERRQ(ierr); 210c4762a1bSJed Brown 211c4762a1bSJed Brown switch (cklvl) { 212c4762a1bSJed Brown case 2: 213c4762a1bSJed Brown dot_max = 0.0; 214c4762a1bSJed Brown for (i = il; i<iu; i++) { 215c4762a1bSJed Brown ierr = VecCopy(evec[i], vt1);CHKERRQ(ierr); 216c4762a1bSJed Brown for (j=il; j<iu; j++) { 217c4762a1bSJed Brown ierr = VecDot(evec[j],vt1,&dot);CHKERRQ(ierr); 218c4762a1bSJed Brown if (j == i) { 219c4762a1bSJed Brown dot = PetscAbsScalar(dot - 1); 220c4762a1bSJed Brown } else { 221c4762a1bSJed Brown dot = PetscAbsScalar(dot); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown if (dot > dot_max) dot_max = dot; 224c4762a1bSJed Brown if (dot > tols[1]) { 225c4762a1bSJed Brown ierr = VecNorm(evec[i],NORM_INFINITY,&norm);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"|delta(%D,%D)|: %g, norm: %g\n",i,j,(double)dot,(double)norm);CHKERRQ(ierr); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 229c4762a1bSJed Brown } 230c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max);CHKERRQ(ierr); 231c4762a1bSJed Brown 232c4762a1bSJed Brown case 1: 233c4762a1bSJed Brown norm_max = 0.0; 234c4762a1bSJed Brown for (i = il; i< iu; i++) { 235c4762a1bSJed Brown ierr = MatMult(A, evec[i], vt1);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = VecCopy(evec[i], vt2);CHKERRQ(ierr); 237c4762a1bSJed Brown tmp = -eval[i]; 238c4762a1bSJed Brown ierr = VecAXPY(vt1,tmp,vt2);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = VecNorm(vt1, NORM_INFINITY, &norm);CHKERRQ(ierr); 240c4762a1bSJed Brown norm = PetscAbsScalar(norm); 241c4762a1bSJed Brown if (norm > norm_max) norm_max = norm; 242c4762a1bSJed Brown /* sniff, and bark if necessary */ 243c4762a1bSJed Brown if (norm > tols[0]) { 244c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," residual violation: %D, resi: %g\n",i, norm);CHKERRQ(ierr); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 247c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," max_resi: %g\n", (double)norm_max);CHKERRQ(ierr); 248c4762a1bSJed Brown break; 249c4762a1bSJed Brown default: 250c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%D is not supported \n",cklvl);CHKERRQ(ierr); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown ierr = VecDestroy(&vt2);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = VecDestroy(&vt1);CHKERRQ(ierr); 254c4762a1bSJed Brown PetscFunctionReturn(0); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown /*TEST 258c4762a1bSJed Brown 259c4762a1bSJed Brown build: 260c4762a1bSJed Brown requires: !complex 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 263*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 264c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 265c4762a1bSJed Brown output_file: output/ex116_1.out 266c4762a1bSJed Brown 267c4762a1bSJed Brown test: 268c4762a1bSJed Brown suffix: 2 269*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 270c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -test_syev -check_symmetry 271c4762a1bSJed Brown 272c4762a1bSJed Brown TEST*/ 273