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