xref: /petsc/src/mat/tests/ex116.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscBool      isSymmetric;
23c4762a1bSJed Brown   PetscScalar    *arrayA,*evecs_array,*work,*evals;
24c4762a1bSJed Brown   PetscMPIInt    size;
25c4762a1bSJed Brown   PetscInt       m,n,i,cklvl=2;
26c4762a1bSJed Brown   PetscBLASInt   nevs,il,iu,in;
27c4762a1bSJed Brown   PetscReal      vl,vu,abstol=1.e-8;
28c4762a1bSJed Brown   PetscBLASInt   *iwork,*ifail,lwork,lierr,bn;
29c4762a1bSJed Brown   PetscReal      tols[2];
30c4762a1bSJed Brown 
31*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
325f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
332c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
34c4762a1bSJed Brown 
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_syev", &flg));
36c4762a1bSJed Brown   if (flg) {
37c4762a1bSJed Brown     TestSYEVX = PETSC_FALSE;
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Determine files from which we read the two matrices */
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file[0],sizeof(file[0]),&flg));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Load matrix A */
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&m,&n));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Check whether A is symmetric */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-check_symmetry", &flg));
53c4762a1bSJed Brown   if (flg) {
54c4762a1bSJed Brown     Mat Trans;
555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Trans));
565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatEqual(A, Trans, &isSymmetric));
5728b400f6SJacob Faibussowitsch     PetscCheck(isSymmetric,PETSC_COMM_SELF,PETSC_ERR_USER,"A must be symmetric");
585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Trans));
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Solve eigenvalue problem: A_dense*x = lambda*B*x */
62c4762a1bSJed Brown   /*==================================================*/
63c4762a1bSJed Brown   /* Convert aij matrix to MatSeqDense for LAPACK */
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense));
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBLASIntCast(8*n,&lwork));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBLASIntCast(n,&bn));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&evals));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(lwork,&work));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArray(A_dense,&arrayA));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   if (!TestSYEVX) { /* test syev() */
735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF," LAPACKsyev: compute all %" PetscInt_FMT " eigensolutions...\n",m));
74c4762a1bSJed Brown     LAPACKsyev_("V","U",&bn,arrayA,&bn,evals,work,&lwork,&lierr);
75c4762a1bSJed Brown     evecs_array = arrayA;
765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(m,&nevs));
77c4762a1bSJed Brown     il          = 1;
785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(m,&iu));
79c4762a1bSJed Brown   } else { /* test syevx()  */
80c4762a1bSJed Brown     il   = 1;
815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(0.2*m,&iu));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&in));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF," LAPACKsyevx: compute %" PetscBLASInt_FMT " to %" PetscBLASInt_FMT "-th eigensolutions...\n",il,iu));
845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(m*n+1,&evecs_array));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(6*n+1,&iwork));
86c4762a1bSJed Brown     ifail = iwork + 5*n;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown     /* in the case "I", vl and vu are not referenced */
89c4762a1bSJed Brown     vl = 0.0; vu = 8.0;
90c4762a1bSJed Brown     LAPACKsyevx_("V","I","U",&bn,arrayA,&bn,&vl,&vu,&il,&iu,&abstol,&nevs,evals,evecs_array,&in,work,&lwork,iwork,ifail,&lierr);
915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(iwork));
92c4762a1bSJed Brown   }
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(A_dense,&arrayA));
942c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nevs <= 0,PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED, "nev=%" PetscBLASInt_FMT ", no eigensolution has found", nevs);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* View eigenvalues */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-eig_view", &flg));
98c4762a1bSJed Brown   if (flg) {
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %" PetscBLASInt_FMT " evals: \n",nevs));
1005f80ce2aSJacob Faibussowitsch     for (i=0; i<nevs; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT "  %g\n",(PetscInt)(i+il),(double)evals[i]));
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Check residuals and orthogonality */
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nevs+1,&evecs));
105c4762a1bSJed Brown   for (i=0; i<nevs; i++) {
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreate(PETSC_COMM_SELF,&evecs[i]));
1075f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetSizes(evecs[i],PETSC_DECIDE,n));
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetFromOptions(evecs[i]));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPlaceArray(evecs[i],evecs_array+i*n));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   tols[0] = tols[1] = PETSC_SQRT_MACHINE_EPSILON;
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(CkEigenSolutions(cklvl,A,il-1,iu-1,evals,evecs,tols));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Free work space. */
1165f80ce2aSJacob Faibussowitsch   for (i=0; i<nevs; i++) CHKERRQ(VecDestroy(&evecs[i]));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(evecs));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A_dense));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(work));
1205f80ce2aSJacob Faibussowitsch   if (TestSYEVX) CHKERRQ(PetscFree(evecs_array));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* Compute SVD: A_dense = U*SIGMA*transpose(V),
123c4762a1bSJed Brown      JOBU=JOBV='S':  the first min(m,n) columns of U and V are returned in the arrayU and arrayV; */
124c4762a1bSJed Brown   /*==============================================================================================*/
125c4762a1bSJed Brown   {
126c4762a1bSJed Brown     /* Convert aij matrix to MatSeqDense for LAPACK */
127c4762a1bSJed Brown     PetscScalar  *arrayU,*arrayVT,*arrayErr,alpha=1.0,beta=-1.0;
128c4762a1bSJed Brown     Mat          Err;
129c4762a1bSJed Brown     PetscBLASInt minMN,maxMN,im,in;
130c4762a1bSJed Brown     PetscInt     j;
131c4762a1bSJed Brown     PetscReal    norm;
132c4762a1bSJed Brown 
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&A_dense));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown     minMN = PetscMin(m,n);
136c4762a1bSJed Brown     maxMN = PetscMax(m,n);
137c4762a1bSJed Brown     lwork = 5*minMN + maxMN;
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc4(m*minMN,&arrayU,m*minMN,&arrayVT,m*minMN,&arrayErr,lwork,&work));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     /* Create matrix Err for checking error */
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Err));
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(Err,PETSC_DECIDE,PETSC_DECIDE,m,minMN));
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(Err,MATSEQDENSE));
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSeqDenseSetPreallocation(Err,(PetscScalar*)arrayErr));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown     /* Save A to arrayErr for checking accuracy later. arrayA will be destroyed by LAPACKgesvd_() */
1475f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseGetArray(A_dense,&arrayA));
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(arrayErr,arrayA,m*minMN));
149c4762a1bSJed Brown 
1505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(m,&im));
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&in));
152c4762a1bSJed Brown     /* Compute A = U*SIGMA*VT */
153c4762a1bSJed Brown     LAPACKgesvd_("S","S",&im,&in,arrayA,&im,evals,arrayU,&minMN,arrayVT,&minMN,work,&lwork,&lierr);
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseRestoreArray(A_dense,&arrayA));
155c4762a1bSJed Brown     if (!lierr) {
1565f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF," 1st 10 of %" PetscBLASInt_FMT " singular values: \n",minMN));
1575f80ce2aSJacob Faibussowitsch       for (i=0; i<10; i++) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT "  %g\n",i,(double)evals[i]));
158c4762a1bSJed Brown     } else {
1595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"LAPACKgesvd_ fails!"));
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     /* Check Err = (U*Sigma*V^T - A) using BLASgemm() */
163c4762a1bSJed Brown     /* U = U*Sigma */
164c4762a1bSJed Brown     for (j=0; j<minMN; j++) { /* U[:,j] = sigma[j]*U[:,j] */
165c4762a1bSJed Brown       for (i=0; i<m; i++) arrayU[j*m+i] *= evals[j];
166c4762a1bSJed Brown     }
167c4762a1bSJed Brown     /* Err = U*VT - A = alpha*U*VT + beta*Err */
168c4762a1bSJed Brown     BLASgemm_("N","N",&im,&minMN,&minMN,&alpha,arrayU,&im,arrayVT,&minMN,&beta,arrayErr,&im);
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(Err,NORM_FROBENIUS,&norm));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF," || U*Sigma*VT - A || = %g\n",(double)norm));
171c4762a1bSJed Brown 
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree4(arrayU,arrayVT,arrayErr,work));
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(evals));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A_dense));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Err));
176c4762a1bSJed Brown   }
177c4762a1bSJed Brown 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
179*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
180*b122ec5aSJacob Faibussowitsch   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown /*------------------------------------------------
183c4762a1bSJed Brown   Check the accuracy of the eigen solution
184c4762a1bSJed Brown   ----------------------------------------------- */
185c4762a1bSJed Brown /*
186c4762a1bSJed Brown   input:
187c4762a1bSJed Brown      cklvl      - check level:
188c4762a1bSJed Brown                     1: check residual
189c4762a1bSJed Brown                     2: 1 and check B-orthogonality locally
190c4762a1bSJed Brown      A          - matrix
191c4762a1bSJed Brown      il,iu      - lower and upper index bound of eigenvalues
192c4762a1bSJed Brown      eval, evec - eigenvalues and eigenvectors stored in this process
193c4762a1bSJed Brown      tols[0]    - reporting tol_res: || A * evec[i] - eval[i]*evec[i] ||
194c4762a1bSJed Brown      tols[1]    - reporting tol_orth: evec[i]^T*evec[j] - delta_ij
195c4762a1bSJed Brown */
196c4762a1bSJed Brown PetscErrorCode CkEigenSolutions(PetscInt cklvl,Mat A,PetscInt il,PetscInt iu,PetscReal *eval,Vec *evec,PetscReal *tols)
197c4762a1bSJed Brown {
1985f80ce2aSJacob Faibussowitsch   PetscInt  i,j,nev;
199c4762a1bSJed Brown   Vec       vt1,vt2;    /* tmp vectors */
200c4762a1bSJed Brown   PetscReal norm,tmp,dot,norm_max,dot_max;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   PetscFunctionBegin;
203c4762a1bSJed Brown   nev = iu - il;
204c4762a1bSJed Brown   if (nev <= 0) PetscFunctionReturn(0);
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /*ierr = VecView(evec[0],PETSC_VIEWER_STDOUT_WORLD);*/
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(evec[0],&vt1));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(evec[0],&vt2));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   switch (cklvl) {
211c4762a1bSJed Brown   case 2:
212c4762a1bSJed Brown     dot_max = 0.0;
213c4762a1bSJed Brown     for (i = il; i<iu; i++) {
2145f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(evec[i], vt1));
215c4762a1bSJed Brown       for (j=il; j<iu; j++) {
2165f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDot(evec[j],vt1,&dot));
217c4762a1bSJed Brown         if (j == i) {
218c4762a1bSJed Brown           dot = PetscAbsScalar(dot - 1);
219c4762a1bSJed Brown         } else {
220c4762a1bSJed Brown           dot = PetscAbsScalar(dot);
221c4762a1bSJed Brown         }
222c4762a1bSJed Brown         if (dot > dot_max) dot_max = dot;
223c4762a1bSJed Brown         if (dot > tols[1]) {
2245f80ce2aSJacob Faibussowitsch           CHKERRQ(VecNorm(evec[i],NORM_INFINITY,&norm));
2255f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"|delta(%" PetscInt_FMT ",%" PetscInt_FMT ")|: %g, norm: %g\n",i,j,(double)dot,(double)norm));
226c4762a1bSJed Brown         }
227c4762a1bSJed Brown       }
228c4762a1bSJed Brown     }
2295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"    max|(x_j^T*x_i) - delta_ji|: %g\n",(double)dot_max));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   case 1:
232c4762a1bSJed Brown     norm_max = 0.0;
233c4762a1bSJed Brown     for (i = il; i< iu; i++) {
2345f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A, evec[i], vt1));
2355f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(evec[i], vt2));
236c4762a1bSJed Brown       tmp  = -eval[i];
2375f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(vt1,tmp,vt2));
2385f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(vt1, NORM_INFINITY, &norm));
239c4762a1bSJed Brown       norm = PetscAbsScalar(norm);
240c4762a1bSJed Brown       if (norm > norm_max) norm_max = norm;
241c4762a1bSJed Brown       /* sniff, and bark if necessary */
242c4762a1bSJed Brown       if (norm > tols[0]) {
2435f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  residual violation: %" PetscInt_FMT ", resi: %g\n",i, (double)norm));
244c4762a1bSJed Brown       }
245c4762a1bSJed Brown     }
2465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"    max_resi:                    %g\n", (double)norm_max));
247c4762a1bSJed Brown     break;
248c4762a1bSJed Brown   default:
2495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: cklvl=%" PetscInt_FMT " is not supported \n",cklvl));
250c4762a1bSJed Brown   }
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vt2));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vt1));
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown /*TEST
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    build:
259c4762a1bSJed Brown       requires: !complex
260c4762a1bSJed Brown 
261c4762a1bSJed Brown    test:
262dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
263c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small
264c4762a1bSJed Brown       output_file: output/ex116_1.out
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    test:
267c4762a1bSJed Brown       suffix: 2
268dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
269c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -test_syev -check_symmetry
270c4762a1bSJed Brown 
271c4762a1bSJed Brown TEST*/
272