Actual source code: ex13.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[] = "Generalized Symmetric eigenproblem.\n\n"
12: "The problem is Ax = lambda Bx, with:\n"
13: " A = Laplacian operator in 2-D\n"
14: " B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
15: "The command line options are:\n"
16: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
17: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
18: " -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";
20: #include slepceps.h
24: int main( int argc, char **argv )
25: {
26: Mat A, B; /* matrices */
27: EPS eps; /* eigenproblem solver context */
28: EPSType type;
29: PetscReal error, tol, re, im;
30: PetscScalar kr, ki;
32: PetscInt N, n=10, m, Istart, Iend, II, J;
33: int nev, maxit, i, j, its, nconv, nulldim=0;
34: PetscScalar v;
35: PetscTruth flag;
37: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
40: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
41: if( flag==PETSC_FALSE ) m=n;
42: N = n*m;
43: PetscOptionsGetInt(PETSC_NULL,"-nulldim",&nulldim,PETSC_NULL);
44: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%d (%dx%d grid), null(B)=%d\n\n",N,n,m,nulldim);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrices that define the eigensystem, Ax=kBx
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: MatCreate(PETSC_COMM_WORLD,&A);
51: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
52: MatSetFromOptions(A);
53:
54: MatCreate(PETSC_COMM_WORLD,&B);
55: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
56: MatSetFromOptions(B);
58: MatGetOwnershipRange(A,&Istart,&Iend);
59: for( II=Istart; II<Iend; II++ ) {
60: v = -1.0; i = II/n; j = II-i*n;
61: if(i>0) { J=II-n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
62: if(i<m-1) { J=II+n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
63: if(j>0) { J=II-1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
64: if(j<n-1) { J=II+1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES); }
65: v=4.0; MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES);
66: if (II>=nulldim) { v=4.0; MatSetValues(B,1,&II,1,&II,&v,INSERT_VALUES); }
67: }
69: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
71: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create the eigensolver and set various options
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: /*
79: Create eigensolver context
80: */
81: EPSCreate(PETSC_COMM_WORLD,&eps);
83: /*
84: Set operators. In this case, it is a generalized eigenvalue problem
85: */
86: EPSSetOperators(eps,A,B);
87: EPSSetProblemType(eps,EPS_GHEP);
89: /*
90: Set solver parameters at runtime
91: */
92: EPSSetFromOptions(eps);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Solve the eigensystem
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: EPSSolve(eps);
99: EPSGetIterationNumber(eps, &its);
100: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
102: /*
103: Optional: Get some information from the solver and display it
104: */
105: EPSGetType(eps,&type);
106: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
107: EPSGetDimensions(eps,&nev,PETSC_NULL);
108: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
109: EPSGetTolerances(eps,&tol,&maxit);
110: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Display solution and clean up
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: /*
117: Get number of converged approximate eigenpairs
118: */
119: EPSGetConverged(eps,&nconv);
120: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
121:
123: if (nconv>0) {
124: /*
125: Display eigenvalues and relative errors
126: */
127: PetscPrintf(PETSC_COMM_WORLD,
128: " k ||Ax-kBx||/||kx||\n"
129: " ----------------- ------------------\n" );
131: for( i=0; i<nconv; i++ ) {
132: /*
133: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
134: ki (imaginary part)
135: */
136: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
137: /*
138: Compute the relative error associated to each eigenpair
139: */
140: EPSComputeRelativeError(eps,i,&error);
142: #ifdef PETSC_USE_COMPLEX
143: re = PetscRealPart(kr);
144: im = PetscImaginaryPart(kr);
145: #else
146: re = kr;
147: im = ki;
148: #endif
149: if (im!=0.0) {
150: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
151: } else {
152: PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,error);
153: }
154: }
155: PetscPrintf(PETSC_COMM_WORLD,"\n" );
156: }
157:
158: /*
159: Free work space
160: */
161: EPSDestroy(eps);
162: MatDestroy(A);
163: SlepcFinalize();
164: return 0;
165: }