Actual source code: ex7.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 generalized eigensystem Ax=kBx with matrices 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: " -f1 <filename>, where <filename> = matrix (A) file in PETSc binary form.\n"
26: " -f2 <filename>, where <filename> = matrix (B) file in PETSc binary form.\n\n";
28: #include slepceps.h
32: int main( int argc, char **argv )
33: {
34: Mat A,B; /* matrices */
35: EPS eps; /* eigenproblem solver context */
36: const EPSType type;
37: PetscReal error, tol, re, im;
38: PetscScalar kr, ki;
40: PetscInt nev, maxit, i, its, lits, nconv;
41: char filename[256];
42: PetscViewer viewer;
43: PetscTruth flg;
46: SlepcInitialize(&argc,&argv,(char*)0,help);
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Load the matrices that define the eigensystem, Ax=kBx
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem stored in file.\n\n");
53: PetscOptionsGetString(PETSC_NULL,"-f1",filename,256,&flg);
54: if (!flg) {
55: SETERRQ(1,"Must indicate a file name for matrix A with the -f1 option.");
56: }
58: #if defined(PETSC_USE_COMPLEX)
59: PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");
60: #else
61: PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");
62: #endif
63: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
64: MatLoad(viewer,MATAIJ,&A);
65: PetscViewerDestroy(viewer);
67: PetscOptionsGetString(PETSC_NULL,"-f2",filename,256,&flg);
68: if (!flg) {
69: SETERRQ(1,"Must indicate a file name for matrix B with the -f2 option.");
70: }
72: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
73: MatLoad(viewer,MATAIJ,&B);
74: PetscViewerDestroy(viewer);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Create the eigensolver and set various options
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: /*
81: Create eigensolver context
82: */
83: EPSCreate(PETSC_COMM_WORLD,&eps);
85: /*
86: Set operators. In this case, it is a generalized eigenvalue problem
87: */
88: EPSSetOperators(eps,A,B);
90: /*
91: Set solver parameters at runtime
92: */
93: EPSSetFromOptions(eps);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Solve the eigensystem
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: EPSSolve(eps);
101: /*
102: Optional: Get some information from the solver and display it
103: */
104: EPSGetIterationNumber(eps, &its);
105: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
106: EPSGetOperationCounters(eps,PETSC_NULL,PETSC_NULL,&lits);
107: PetscPrintf(PETSC_COMM_WORLD," Number of linear iterations of the method: %d\n",lits);
108: EPSGetType(eps,&type);
109: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
110: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
111: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
112: EPSGetTolerances(eps,&tol,&maxit);
113: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Display solution and clean up
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: /*
120: Get number of converged eigenpairs
121: */
122: EPSGetConverged(eps,&nconv);
123: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
125: if (nconv>0) {
126: /*
127: Display eigenvalues and relative errors
128: */
129: PetscPrintf(PETSC_COMM_WORLD,
130: " k ||Ax-kBx||/||kx||\n"
131: " --------------------- ------------------\n" );
132: for( i=0; i<nconv; i++ ) {
133: /*
134: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
135: ki (imaginary part)
136: */
137: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
139: /*
140: Compute the relative error associated to each eigenpair
141: */
142: EPSComputeRelativeError(eps,i,&error);
144: #if defined(PETSC_USE_COMPLEX)
145: re = PetscRealPart(kr);
146: im = PetscImaginaryPart(kr);
147: #else
148: re = kr;
149: im = ki;
150: #endif
151: if( im != 0.0 ) {
152: PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
153: } else {
154: PetscPrintf(PETSC_COMM_WORLD," % 6f ",re);
155: }
156: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
157: }
158: PetscPrintf(PETSC_COMM_WORLD,"\n" );
159: }
160:
161: /*
162: Free work space
163: */
164: EPSDestroy(eps);
165: MatDestroy(A);
166: MatDestroy(B);
167: SlepcFinalize();
168: return 0;
169: }