Actual source code: ex7.c

  2: static char help[] = "Solves a generalized eigensystem Ax=kBx with matrices loaded from a file.\n\n"
  3:   "This example works for both real and complex numbers.\n\n"
  4:   "The command line options are:\n\n"
  5:   "  -f1 <filename>, where <filename> = matrix (A) file in PETSc binary form.\n"
  6:   "  -f2 <filename>, where <filename> = matrix (B) file in PETSc binary form.\n\n";

 8:  #include slepceps.h

 12: int main( int argc, char **argv )
 13: {
 14:   Mat         A,B;             /* matrices */
 15:   EPS         eps;             /* eigenproblem solver context */
 16:   EPSType     type;
 17:   PetscReal   error, tol, re, im;
 18:   PetscScalar kr, ki;
 19:   int         nev, ierr, maxit, i, its, lits, nconv;
 20:   char        filename[256];
 21:   PetscViewer viewer;
 22:   PetscTruth  flg;


 25:   SlepcInitialize(&argc,&argv,(char*)0,help);

 27:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 28:         Load the matrices that define the eigensystem, Ax=kBx
 29:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 31:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem stored in file.\n\n");
 32:   PetscOptionsGetString(PETSC_NULL,"-f1",filename,256,&flg);
 33:   if (!flg) {
 34:     SETERRQ(1,"Must indicate a file name for matrix A with the -f1 option.");
 35:   }

 37: #if defined(PETSC_USE_COMPLEX)
 38:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");
 39: #else
 40:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");
 41: #endif
 42:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,PETSC_FILE_RDONLY,&viewer);
 43:   MatLoad(viewer,MATMPIAIJ,&A);
 44:   PetscViewerDestroy(viewer);

 46:   PetscOptionsGetString(PETSC_NULL,"-f2",filename,256,&flg);
 47:   if (!flg) {
 48:     SETERRQ(1,"Must indicate a file name for matrix B with the -f2 option.");
 49:   }

 51:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,PETSC_FILE_RDONLY,&viewer);
 52:   MatLoad(viewer,MATMPIAIJ,&B);
 53:   PetscViewerDestroy(viewer);

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 56:                 Create the eigensolver and set various options
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 59:   /* 
 60:      Create eigensolver context
 61:   */
 62:   EPSCreate(PETSC_COMM_WORLD,&eps);

 64:   /* 
 65:      Set operators. In this case, it is a generalized eigenvalue problem
 66:   */
 67:   EPSSetOperators(eps,A,B);

 69:   /*
 70:      Set solver parameters at runtime
 71:   */
 72:   EPSSetFromOptions(eps);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 75:                       Solve the eigensystem
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   EPSSolve(eps);

 80:   /*
 81:      Optional: Get some information from the solver and display it
 82:   */
 83:   EPSGetIterationNumber(eps, &its);
 84:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
 85:   EPSGetNumberLinearIterations(eps, &lits);
 86:   PetscPrintf(PETSC_COMM_WORLD," Number of linear iterations of the method: %d\n",lits);
 87:   EPSGetType(eps,&type);
 88:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 89:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 90:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 91:   EPSGetTolerances(eps,&tol,&maxit);
 92:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 95:                     Display solution and clean up
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   /* 
 99:      Get number of converged eigenpairs
100:   */
101:   EPSGetConverged(eps,&nconv);
102:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

104:   if (nconv>0) {
105:     /*
106:        Display eigenvalues and relative errors
107:     */
108:     PetscPrintf(PETSC_COMM_WORLD,
109:          "           k             ||Ax-kx||/||kx||\n"
110:          "  --------------------- ------------------\n" );
111:     for( i=0; i<nconv; i++ ) {
112:       /* 
113:          Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
114:          ki (imaginary part)
115:       */
116:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);

118:       /*
119:          Compute the relative error associated to each eigenpair
120:       */
121:       EPSComputeRelativeError(eps,i,&error);

123: #if defined(PETSC_USE_COMPLEX)
124:       re = PetscRealPart(kr);
125:       im = PetscImaginaryPart(kr);
126: #else
127:       re = kr;
128:       im = ki;
129: #endif
130:       if( im != 0.0 ) {
131:         PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
132:       } else {
133:         PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",re);
134:       }
135:       PetscPrintf(PETSC_COMM_WORLD," % 12f\n",error);
136:     }
137:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
138:   }
139: 
140:   /* 
141:      Free work space
142:   */
143:   EPSDestroy(eps);
144:   MatDestroy(A);
145:   MatDestroy(B);
146:   SlepcFinalize();
147:   return 0;
148: }