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 159371c9d4SSatish Balay int main(int argc, char **args) { 16c4762a1bSJed Brown Mat A, A_dense; 17c4762a1bSJed Brown Vec *evecs; 18c4762a1bSJed Brown PetscViewer fd; /* viewer */ 19c4762a1bSJed Brown char file[1][PETSC_MAX_PATH_LEN]; /* input file name */ 20c4762a1bSJed Brown PetscBool flg, TestSYEVX = PETSC_TRUE; 21c4762a1bSJed Brown PetscBool isSymmetric; 22c4762a1bSJed Brown PetscScalar *arrayA, *evecs_array, *work, *evals; 23c4762a1bSJed Brown PetscMPIInt size; 24c4762a1bSJed Brown PetscInt m, n, i, cklvl = 2; 25c4762a1bSJed Brown PetscBLASInt nevs, il, iu, in; 26c4762a1bSJed Brown PetscReal vl, vu, abstol = 1.e-8; 27c4762a1bSJed Brown PetscBLASInt *iwork, *ifail, lwork, lierr, bn; 28c4762a1bSJed Brown PetscReal tols[2]; 29c4762a1bSJed Brown 30327415f7SBarry Smith PetscFunctionBeginUser; 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 33be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_syev", &flg)); 369371c9d4SSatish Balay if (flg) { TestSYEVX = PETSC_FALSE; } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Determine files from which we read the two matrices */ 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file[0], sizeof(file[0]), &flg)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Load matrix A */ 429566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &fd)); 439566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 449566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ)); 459566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 469566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 479566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Check whether A is symmetric */ 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-check_symmetry", &flg)); 51c4762a1bSJed Brown if (flg) { 52c4762a1bSJed Brown Mat Trans; 539566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Trans)); 549566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Trans, &isSymmetric)); 5528b400f6SJacob Faibussowitsch PetscCheck(isSymmetric, PETSC_COMM_SELF, PETSC_ERR_USER, "A must be symmetric"); 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Trans)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Solve eigenvalue problem: A_dense*x = lambda*B*x */ 60c4762a1bSJed Brown /*==================================================*/ 61c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */ 629566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &A_dense)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(8 * n, &lwork)); 659566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &bn)); 669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &evals)); 679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lwork, &work)); 689566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(A_dense, &arrayA)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown if (!TestSYEVX) { /* test syev() */ 719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n", m)); 72c4762a1bSJed Brown LAPACKsyev_("V", "U", &bn, arrayA, &bn, evals, work, &lwork, &lierr); 73c4762a1bSJed Brown evecs_array = arrayA; 749566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &nevs)); 75c4762a1bSJed Brown il = 1; 769566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &iu)); 77c4762a1bSJed Brown } else { /* test syevx() */ 78c4762a1bSJed Brown il = 1; 799566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(0.2 * m, &iu)); 809566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &in)); 819566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LAPACKsyevx: compute %" PetscBLASInt_FMT " to %" PetscBLASInt_FMT "-th eigensolutions...\n", il, iu)); 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * n + 1, &evecs_array)); 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * n + 1, &iwork)); 84c4762a1bSJed Brown ifail = iwork + 5 * n; 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* in the case "I", vl and vu are not referenced */ 879371c9d4SSatish Balay vl = 0.0; 889371c9d4SSatish Balay vu = 8.0; 89c4762a1bSJed Brown LAPACKsyevx_("V", "I", "U", &bn, arrayA, &bn, &vl, &vu, &il, &iu, &abstol, &nevs, evals, evecs_array, &in, work, &lwork, iwork, ifail, &lierr); 909566063dSJacob Faibussowitsch PetscCall(PetscFree(iwork)); 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(A_dense, &arrayA)); 93e00437b9SBarry Smith PetscCheck(nevs > 0, PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "nev=%" PetscBLASInt_FMT ", no eigensolution has found", nevs); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* View eigenvalues */ 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-eig_view", &flg)); 97c4762a1bSJed Brown if (flg) { 989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscBLASInt_FMT " evals: \n", nevs)); 999566063dSJacob Faibussowitsch for (i = 0; i < nevs; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g\n", (PetscInt)(i + il), (double)evals[i])); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Check residuals and orthogonality */ 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevs + 1, &evecs)); 104c4762a1bSJed Brown for (i = 0; i < nevs; i++) { 1059566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &evecs[i])); 1069566063dSJacob Faibussowitsch PetscCall(VecSetSizes(evecs[i], PETSC_DECIDE, n)); 1079566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(evecs[i])); 1089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(evecs[i], evecs_array + i * n)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown tols[0] = tols[1] = PETSC_SQRT_MACHINE_EPSILON; 1129566063dSJacob Faibussowitsch PetscCall(CkEigenSolutions(cklvl, A, il - 1, iu - 1, evals, evecs, tols)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Free work space. */ 1159566063dSJacob Faibussowitsch for (i = 0; i < nevs; i++) PetscCall(VecDestroy(&evecs[i])); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(evecs)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_dense)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 1199566063dSJacob Faibussowitsch if (TestSYEVX) PetscCall(PetscFree(evecs_array)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* Compute SVD: A_dense = U*SIGMA*transpose(V), 122c4762a1bSJed Brown JOBU=JOBV='S': the first min(m,n) columns of U and V are returned in the arrayU and arrayV; */ 123c4762a1bSJed Brown /*==============================================================================================*/ 124c4762a1bSJed Brown { 125c4762a1bSJed Brown /* Convert aij matrix to MatSeqDense for LAPACK */ 126c4762a1bSJed Brown PetscScalar *arrayU, *arrayVT, *arrayErr, alpha = 1.0, beta = -1.0; 127c4762a1bSJed Brown Mat Err; 128c4762a1bSJed Brown PetscBLASInt minMN, maxMN, im, in; 129c4762a1bSJed Brown PetscInt j; 130c4762a1bSJed Brown PetscReal norm; 131c4762a1bSJed Brown 1329566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &A_dense)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown minMN = PetscMin(m, n); 135c4762a1bSJed Brown maxMN = PetscMax(m, n); 136c4762a1bSJed Brown lwork = 5 * minMN + maxMN; 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(m * minMN, &arrayU, m * minMN, &arrayVT, m * minMN, &arrayErr, lwork, &work)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Create matrix Err for checking error */ 1409566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Err)); 1419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Err, PETSC_DECIDE, PETSC_DECIDE, m, minMN)); 1429566063dSJacob Faibussowitsch PetscCall(MatSetType(Err, MATSEQDENSE)); 1439566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Err, (PetscScalar *)arrayErr)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Save A to arrayErr for checking accuracy later. arrayA will be destroyed by LAPACKgesvd_() */ 1469566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(A_dense, &arrayA)); 1479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(arrayErr, arrayA, m * minMN)); 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(m, &im)); 1509566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &in)); 151c4762a1bSJed Brown /* Compute A = U*SIGMA*VT */ 152c4762a1bSJed Brown LAPACKgesvd_("S", "S", &im, &in, arrayA, &im, evals, arrayU, &minMN, arrayVT, &minMN, work, &lwork, &lierr); 1539566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(A_dense, &arrayA)); 154c4762a1bSJed Brown if (!lierr) { 1559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " 1st 10 of %" PetscBLASInt_FMT " singular values: \n", minMN)); 1569566063dSJacob Faibussowitsch for (i = 0; i < 10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g\n", i, (double)evals[i])); 157c4762a1bSJed Brown } else { 1589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "LAPACKgesvd_ fails!")); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Check Err = (U*Sigma*V^T - A) using BLASgemm() */ 162c4762a1bSJed Brown /* U = U*Sigma */ 163c4762a1bSJed Brown for (j = 0; j < minMN; j++) { /* U[:,j] = sigma[j]*U[:,j] */ 164c4762a1bSJed Brown for (i = 0; i < m; i++) arrayU[j * m + i] *= evals[j]; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown /* Err = U*VT - A = alpha*U*VT + beta*Err */ 167c4762a1bSJed Brown BLASgemm_("N", "N", &im, &minMN, &minMN, &alpha, arrayU, &im, arrayVT, &minMN, &beta, arrayErr, &im); 1689566063dSJacob Faibussowitsch PetscCall(MatNorm(Err, NORM_FROBENIUS, &norm)); 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " || U*Sigma*VT - A || = %g\n", (double)norm)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(PetscFree4(arrayU, arrayVT, arrayErr, work)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFree(evals)); 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_dense)); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Err)); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 179b122ec5aSJacob Faibussowitsch return 0; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown /*------------------------------------------------ 182c4762a1bSJed Brown Check the accuracy of the eigen solution 183c4762a1bSJed Brown ----------------------------------------------- */ 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown input: 186c4762a1bSJed Brown cklvl - check level: 187c4762a1bSJed Brown 1: check residual 188c4762a1bSJed Brown 2: 1 and check B-orthogonality locally 189c4762a1bSJed Brown A - matrix 190c4762a1bSJed Brown il,iu - lower and upper index bound of eigenvalues 191c4762a1bSJed Brown eval, evec - eigenvalues and eigenvectors stored in this process 192c4762a1bSJed Brown tols[0] - reporting tol_res: || A * evec[i] - eval[i]*evec[i] || 193c4762a1bSJed Brown tols[1] - reporting tol_orth: evec[i]^T*evec[j] - delta_ij 194c4762a1bSJed Brown */ 1959371c9d4SSatish Balay PetscErrorCode CkEigenSolutions(PetscInt cklvl, Mat A, PetscInt il, PetscInt iu, PetscReal *eval, Vec *evec, PetscReal *tols) { 1965f80ce2aSJacob Faibussowitsch PetscInt i, j, nev; 197c4762a1bSJed Brown Vec vt1, vt2; /* tmp vectors */ 198c4762a1bSJed Brown PetscReal norm, tmp, dot, norm_max, dot_max; 199c4762a1bSJed Brown 200c4762a1bSJed Brown PetscFunctionBegin; 201c4762a1bSJed Brown nev = iu - il; 202c4762a1bSJed Brown if (nev <= 0) PetscFunctionReturn(0); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /*ierr = VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);*/ 2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0], &vt1)); 2069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(evec[0], &vt2)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown switch (cklvl) { 209c4762a1bSJed Brown case 2: 210c4762a1bSJed Brown dot_max = 0.0; 211c4762a1bSJed Brown for (i = il; i < iu; i++) { 2129566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt1)); 213c4762a1bSJed Brown for (j = il; j < iu; j++) { 2149566063dSJacob Faibussowitsch PetscCall(VecDot(evec[j], vt1, &dot)); 215c4762a1bSJed Brown if (j == i) { 216c4762a1bSJed Brown dot = PetscAbsScalar(dot - 1); 217c4762a1bSJed Brown } else { 218c4762a1bSJed Brown dot = PetscAbsScalar(dot); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown if (dot > dot_max) dot_max = dot; 221c4762a1bSJed Brown if (dot > tols[1]) { 2229566063dSJacob Faibussowitsch PetscCall(VecNorm(evec[i], NORM_INFINITY, &norm)); 2239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n", i, j, (double)dot, (double)norm)); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } 226c4762a1bSJed Brown } 2279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " max|(x_j^T*x_i) - delta_ji|: %g\n", (double)dot_max)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown case 1: 230c4762a1bSJed Brown norm_max = 0.0; 231c4762a1bSJed Brown for (i = il; i < iu; i++) { 2329566063dSJacob Faibussowitsch PetscCall(MatMult(A, evec[i], vt1)); 2339566063dSJacob Faibussowitsch PetscCall(VecCopy(evec[i], vt2)); 234c4762a1bSJed Brown tmp = -eval[i]; 2359566063dSJacob Faibussowitsch PetscCall(VecAXPY(vt1, tmp, vt2)); 2369566063dSJacob Faibussowitsch PetscCall(VecNorm(vt1, NORM_INFINITY, &norm)); 237c4762a1bSJed Brown norm = PetscAbsScalar(norm); 238c4762a1bSJed Brown if (norm > norm_max) norm_max = norm; 239c4762a1bSJed Brown /* sniff, and bark if necessary */ 240*48a46eb9SPierre Jolivet if (norm > tols[0]) PetscCall(PetscPrintf(PETSC_COMM_SELF, " residual violation: %" PetscInt_FMT ", resi: %g\n", i, (double)norm)); 241c4762a1bSJed Brown } 2429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " max_resi: %g\n", (double)norm_max)); 243c4762a1bSJed Brown break; 2449371c9d4SSatish Balay default: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: cklvl=%" PetscInt_FMT " is not supported \n", cklvl)); 245c4762a1bSJed Brown } 2469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt2)); 2479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vt1)); 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /*TEST 252c4762a1bSJed Brown 253c4762a1bSJed Brown build: 254c4762a1bSJed Brown requires: !complex 255c4762a1bSJed Brown 256c4762a1bSJed Brown test: 257dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 258c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 259c4762a1bSJed Brown output_file: output/ex116_1.out 260c4762a1bSJed Brown 261c4762a1bSJed Brown test: 262c4762a1bSJed Brown suffix: 2 263dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 264c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -test_syev -check_symmetry 265c4762a1bSJed Brown 266c4762a1bSJed Brown TEST*/ 267