Actual source code: ex7.c

  2: static char help[] = "Solves a generalized eigensystem Ax=kBx with matrices loaded from a file.\n"
  3:   "This example works for both real and complex numbers.\n\n"
  4:   "The command line options are:\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;
 20:   int                  nev, maxit, i, its, lits, nconv;
 21:   char                 filename[256];
 22:   PetscViewer          viewer;
 23:   PetscTruth           flg;


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

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

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

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

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

 52:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 53:   MatLoad(viewer,MATAIJ,&B);
 54:   PetscViewerDestroy(viewer);

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

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

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

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

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

 79:   EPSSolve(eps);

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

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

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

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

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

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