Actual source code: ex4.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Solves a standard eigensystem Ax=kx with the matrix loaded from a file.\n"
23: "This example works for both real and complex numbers.\n\n"
24: "The command line options are:\n"
25: " -file <filename>, where <filename> = matrix file in PETSc binary form.\n\n";
27: #include slepceps.h
31: int main( int argc, char **argv )
32: {
33: Mat A; /* operator matrix */
34: EPS eps; /* eigenproblem solver context */
35: const EPSType type;
36: PetscReal error, tol, re, im;
37: PetscScalar kr, ki;
39: PetscInt nev, maxit, i, its, nconv;
40: char filename[256];
41: PetscViewer viewer;
42: PetscTruth flg;
45: SlepcInitialize(&argc,&argv,(char*)0,help);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Load the operator matrix that defines the eigensystem, Ax=kx
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: PetscPrintf(PETSC_COMM_WORLD,"\nEigenproblem stored in file.\n\n");
52: PetscOptionsGetString(PETSC_NULL,"-file",filename,256,&flg);
53: if (!flg) {
54: SETERRQ(1,"Must indicate a file name with the -file option.");
55: }
57: #if defined(PETSC_USE_COMPLEX)
58: PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n");
59: #else
60: PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n");
61: #endif
62: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
63: MatLoad(viewer,MATAIJ,&A);
64: PetscViewerDestroy(viewer);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create the eigensolver and set various options
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: /*
71: Create eigensolver context
72: */
73: EPSCreate(PETSC_COMM_WORLD,&eps);
75: /*
76: Set operators. In this case, it is a standard eigenvalue problem
77: */
78: EPSSetOperators(eps,A,PETSC_NULL);
80: /*
81: Set solver parameters at runtime
82: */
83: EPSSetFromOptions(eps);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Solve the eigensystem
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: EPSSolve(eps);
90: EPSGetIterationNumber(eps, &its);
91: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
93: /*
94: Optional: Get some information from the solver and display it
95: */
96: EPSGetType(eps,&type);
97: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
98: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
99: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
100: EPSGetTolerances(eps,&tol,&maxit);
101: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
103: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: Display solution and clean up
105: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: /*
108: Get number of converged eigenpairs
109: */
110: EPSGetConverged(eps,&nconv);
111: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
113: if (nconv>0) {
114: /*
115: Display eigenvalues and relative errors
116: */
117: PetscPrintf(PETSC_COMM_WORLD,
118: " k ||Ax-kx||/||kx||\n"
119: " --------------------- ------------------\n" );
120: for( i=0; i<nconv; i++ ) {
121: /*
122: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
123: ki (imaginary part)
124: */
125: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
127: /*
128: Compute the relative error associated to each eigenpair
129: */
130: EPSComputeRelativeError(eps,i,&error);
132: #if defined(PETSC_USE_COMPLEX)
133: re = PetscRealPart(kr);
134: im = PetscImaginaryPart(kr);
135: #else
136: re = kr;
137: im = ki;
138: #endif
139: if( im != 0.0 ) {
140: PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
141: } else {
142: PetscPrintf(PETSC_COMM_WORLD," % 6f ",re);
143: }
144: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
145: }
146: PetscPrintf(PETSC_COMM_WORLD,"\n" );
147: }
148:
149: /*
150: Free work space
151: */
152: EPSDestroy(eps);
153: MatDestroy(A);
154: SlepcFinalize();
155: return 0;
156: }