Actual source code: ex11.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[] = "Computes the smallest nonzero eigenvalue of the Laplacian of a graph.\n\n"
23: "This example illustrates EPSSetDeflationSpace(). The example graph corresponds to a "
24: "2-D regular mesh. The command line options are:\n"
25: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
26: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
28: #include slepceps.h
32: int main( int argc, char **argv )
33: {
34: Mat A; /* operator matrix */
35: Vec x;
36: EPS eps; /* eigenproblem solver context */
37: const EPSType type;
38: PetscReal error, tol, re, im;
39: PetscScalar kr, ki;
41: PetscInt N, n=10, m, i, j, II, Istart, Iend, nev, maxit, its, nconv;
42: PetscScalar w;
43: PetscTruth flag;
45: SlepcInitialize(&argc,&argv,(char*)0,help);
47: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
48: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
49: if( flag==PETSC_FALSE ) m=n;
50: N = n*m;
51: PetscPrintf(PETSC_COMM_WORLD,"\nFiedler vector of a 2-D regular mesh, N=%d (%dx%d grid)\n\n",N,n,m);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Compute the operator matrix that defines the eigensystem, Ax=kx
55: In this example, A = L(G), where L is the Laplacian of graph G, i.e.
56: Lii = degree of node i, Lij = -1 if edge (i,j) exists in G
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
61: MatSetFromOptions(A);
62:
63: MatGetOwnershipRange(A,&Istart,&Iend);
64: for( II=Istart; II<Iend; II++ ) {
65: i = II/n; j = II-i*n;
66: w = 0.0;
67: if(i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); w=w+1.0; }
68: if(i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); w=w+1.0; }
69: if(j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); w=w+1.0; }
70: if(j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); w=w+1.0; }
71: MatSetValue(A,II,II,w,INSERT_VALUES);
72: }
74: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
75: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Create the eigensolver and set various options
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: /*
82: Create eigensolver context
83: */
84: EPSCreate(PETSC_COMM_WORLD,&eps);
86: /*
87: Set operators. In this case, it is a standard eigenvalue problem
88: */
89: EPSSetOperators(eps,A,PETSC_NULL);
90: EPSSetProblemType(eps,EPS_HEP);
91:
92: /*
93: Select portion of spectrum
94: */
95: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
97: /*
98: Set solver parameters at runtime
99: */
100: EPSSetFromOptions(eps);
102: /*
103: Attach deflation space: in this case, the matrix has a constant
104: nullspace, [1 1 ... 1]^T is the eigenvector of the zero eigenvalue
105: */
106: MatGetVecs(A,&x,PETSC_NULL);
107: VecSet(x,1.0);
108: EPSSetDeflationSpace(eps,1,&x);
109: VecDestroy(x);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Solve the eigensystem
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: EPSSolve(eps);
116: EPSGetIterationNumber(eps, &its);
117: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
119: /*
120: Optional: Get some information from the solver and display it
121: */
122: EPSGetType(eps,&type);
123: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
124: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
125: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
126: EPSGetTolerances(eps,&tol,&maxit);
127: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Display solution and clean up
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: /*
134: Get number of converged approximate eigenpairs
135: */
136: EPSGetConverged(eps,&nconv);
137: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
138:
140: if (nconv>0) {
141: /*
142: Display eigenvalues and relative errors
143: */
144: PetscPrintf(PETSC_COMM_WORLD,
145: " k ||Ax-kx||/||kx||\n"
146: " ----------------- ------------------\n" );
148: for( i=0; i<nconv; i++ ) {
149: /*
150: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
151: ki (imaginary part)
152: */
153: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
154: /*
155: Compute the relative error associated to each eigenpair
156: */
157: EPSComputeRelativeError(eps,i,&error);
159: #ifdef PETSC_USE_COMPLEX
160: re = PetscRealPart(kr);
161: im = PetscImaginaryPart(kr);
162: #else
163: re = kr;
164: im = ki;
165: #endif
166: if (im!=0.0) {
167: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
168: } else {
169: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);
170: }
171: }
172: PetscPrintf(PETSC_COMM_WORLD,"\n" );
173: }
174:
175: /*
176: Free work space
177: */
178: EPSDestroy(eps);
179: MatDestroy(A);
180: SlepcFinalize();
181: return 0;
182: }