Actual source code: ex4.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

  6:       This file is part of SLEPc. See the README file for conditions of use
  7:       and additional information.
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Solves a standard eigensystem Ax=kx with the matrix loaded from a file.\n"
 12:   "This example works for both real and complex numbers.\n\n"
 13:   "The command line options are:\n"
 14:   "  -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n";

 16:  #include slepceps.h

 20: int main( int argc, char **argv )
 21: {
 22:   Mat                  A;                  /* operator matrix */
 23:   EPS                  eps;                  /* eigenproblem solver context */
 24:   EPSType              type;
 25:   PetscReal            error, tol, re, im;
 26:   PetscScalar          kr, ki;
 28:   int                  nev, maxit, i, its, nconv;
 29:   char                 filename[256];
 30:   PetscViewer          viewer;
 31:   PetscTruth           flg;


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

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 37:         Load the operator matrix that defines the eigensystem, Ax=kx
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   PetscPrintf(PETSC_COMM_WORLD,"\nEigenproblem stored in file.\n\n");
 41:   PetscOptionsGetString(PETSC_NULL,"-file",filename,256,&flg);
 42:   if (!flg) {
 43:     SETERRQ(1,"Must indicate a file name with the -file option.");
 44:   }

 46: #if defined(PETSC_USE_COMPLEX)
 47:   PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
 48: #else
 49:   PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
 50: #endif
 51:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 52:   MatLoad(viewer,MATAIJ,&A);
 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 standard eigenvalue problem
 66:   */
 67:   EPSSetOperators(eps,A,PETSC_NULL);

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

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

 78:   EPSSolve(eps);
 79:   EPSGetIterationNumber(eps, &its);
 80:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

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

 92:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 93:                     Display solution and clean up
 94:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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