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