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