Actual source code: ex4.c
2: static char help[] = "Solves a standard eigensystem Ax=kx with the matrix 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: " -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;
19: int nev, maxit, i, its, nconv;
20: char filename[256];
21: PetscViewer viewer;
22: PetscTruth flg;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
27: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28: Load the operator matrix that defines the eigensystem, Ax=kx
29: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
31: PetscPrintf(PETSC_COMM_WORLD,"\nEigenproblem stored in file.\n\n");
32: PetscOptionsGetString(PETSC_NULL,"-file",filename,256,&flg);
33: if (!flg) {
34: SETERRQ(1,"Must indicate a file name with the -file option.");
35: }
37: #if defined(PETSC_USE_COMPLEX)
38: PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
39: #else
40: PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
41: #endif
42: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
43: MatLoad(viewer,MATAIJ,&A);
44: PetscViewerDestroy(viewer);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Create the eigensolver and set various options
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: /*
51: Create eigensolver context
52: */
53: EPSCreate(PETSC_COMM_WORLD,&eps);
55: /*
56: Set operators. In this case, it is a standard eigenvalue problem
57: */
58: EPSSetOperators(eps,A,PETSC_NULL);
60: /*
61: Set solver parameters at runtime
62: */
63: EPSSetFromOptions(eps);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Solve the eigensystem
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: EPSSolve(eps);
70: EPSGetIterationNumber(eps, &its);
71: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
73: /*
74: Optional: Get some information from the solver and display it
75: */
76: EPSGetType(eps,&type);
77: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
78: EPSGetDimensions(eps,&nev,PETSC_NULL);
79: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
80: EPSGetTolerances(eps,&tol,&maxit);
81: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Display solution and clean up
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: /*
88: Get number of converged eigenpairs
89: */
90: EPSGetConverged(eps,&nconv);
91: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
93: if (nconv>0) {
94: /*
95: Display eigenvalues and relative errors
96: */
97: PetscPrintf(PETSC_COMM_WORLD,
98: " k ||Ax-kx||/||kx||\n"
99: " --------------------- ------------------\n" );
100: for( i=0; i<nconv; i++ ) {
101: /*
102: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
103: ki (imaginary part)
104: */
105: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
107: /*
108: Compute the relative error associated to each eigenpair
109: */
110: EPSComputeRelativeError(eps,i,&error);
112: #if defined(PETSC_USE_COMPLEX)
113: re = PetscRealPart(kr);
114: im = PetscImaginaryPart(kr);
115: #else
116: re = kr;
117: im = ki;
118: #endif
119: if( im != 0.0 ) {
120: PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
121: } else {
122: PetscPrintf(PETSC_COMM_WORLD," % 6f ",re);
123: }
124: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
125: }
126: PetscPrintf(PETSC_COMM_WORLD,"\n" );
127: }
128:
129: /*
130: Free work space
131: */
132: EPSDestroy(eps);
133: MatDestroy(A);
134: SlepcFinalize();
135: return 0;
136: }