Actual source code: ex4.c

  2: static char help[] = "Solves an standard eigensystem Ax=kx with the matrix 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:   "  -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n";

 7:  #include slepceps.h

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


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

 26:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 27:         Load the operator matrix that defines the eigensystem, Ax=kx
 28:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 46:                 Create the eigensolver and set various options
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 49:   /* 
 50:      Create eigensolver context
 51:   */
 52:   EPSCreate(PETSC_COMM_WORLD,&eps);

 54:   /* 
 55:      Set operators. In this case, it is a standard eigenvalue problem
 56:   */
 57:   EPSSetOperators(eps,A,PETSC_NULL);

 59:   /*
 60:      Set solver parameters at runtime
 61:   */
 62:   EPSSetFromOptions(eps);

 64:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 65:                       Solve the eigensystem
 66:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 68:   EPSSolve(eps);
 69:   EPSGetIterationNumber(eps, &its);
 70:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 72:   /*
 73:      Optional: Get some information from the solver and display it
 74:   */
 75:   EPSGetType(eps,&type);
 76:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 77:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 78:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 79:   EPSGetTolerances(eps,&tol,&maxit);
 80:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 83:                     Display solution and clean up
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   /* 
 87:      Get number of converged eigenpairs
 88:   */
 89:   EPSGetConverged(eps,&nconv);
 90:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

 92:   if (nconv>0) {
 93:     /*
 94:        Display eigenvalues and relative errors
 95:     */
 96:     PetscPrintf(PETSC_COMM_WORLD,
 97:          "           k             ||Ax-kx||/||kx||\n"
 98:          "  --------------------- ------------------\n" );
 99:     for( i=0; i<nconv; i++ ) {
100:       /* 
101:          Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
102:          ki (imaginary part)
103:       */
104:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);

106:       /*
107:          Compute the relative error associated to each eigenpair
108:       */
109:       EPSComputeRelativeError(eps,i,&error);

111: #if defined(PETSC_USE_COMPLEX)
112:       re = PetscRealPart(kr);
113:       im = PetscImaginaryPart(kr);
114: #else
115:       re = kr;
116:       im = ki;
117: #endif
118:       if( im != 0.0 ) {
119:         PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
120:       } else {
121:         PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",re);
122:       }
123:       PetscPrintf(PETSC_COMM_WORLD," % 12f\n",error);
124:     }
125:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
126:   }
127: 
128:   /* 
129:      Free work space
130:   */
131:   EPSDestroy(eps);
132:   MatDestroy(A);
133:   SlepcFinalize();
134:   return 0;
135: }