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\n"
5: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n\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: ST st;
18: EPSType type;
19: PetscReal error, tol, re, im;
20: PetscScalar kr, ki;
21: int N, n=10, m, nev, ierr, maxit, i, j, I, J, its, nconv, Istart, Iend;
22: PetscScalar v, w, one = 1.0;
23: PetscTruth flag;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
28: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
29: if( flag==PETSC_FALSE ) m=n;
30: N = n*m;
31: PetscPrintf(PETSC_COMM_WORLD,"\nFiedler vector of a 2-D regular mesh, N=%d (%dx%d grid)\n\n",N,n,m);
33: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: Compute the operator matrix that defines the eigensystem, Ax=kx
35: In this example, A = L(G), where L is the Laplacian of graph G, i.e.
36: Lii = degree of node i, Lij = -1 if edge (i,j) exists in G
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
40: MatSetFromOptions(A);
41:
42: MatGetOwnershipRange(A,&Istart,&Iend);
43: for( I=Istart; I<Iend; I++ ) {
44: v = -1.0; i = I/n; j = I-i*n;
45: w = 0.0;
46: if(i>0) { J=I-n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); w=w+1.0; }
47: if(i<m-1) { J=I+n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); w=w+1.0; }
48: if(j>0) { J=I-1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); w=w+1.0; }
49: if(j<n-1) { J=I+1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES); w=w+1.0; }
50: MatSetValues(A,1,&I,1,&I,&w,INSERT_VALUES);
51: }
53: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create the eigensolver and set various options
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: /*
61: Create eigensolver context
62: */
63: EPSCreate(PETSC_COMM_WORLD,&eps);
65: /*
66: Set operators. In this case, it is a standard eigenvalue problem
67: */
68: EPSSetOperators(eps,A,PETSC_NULL);
69: EPSSetProblemType(eps,EPS_HEP);
70:
71: /*
72: Specify some options: use shift-and-invert to compute eigenpairs
73: close to the origin
74: */
75: EPSGetST(eps,&st);
76: STSetType(st,STSINV);
77: STSetShift(st,0);
79: /*
80: Set solver parameters at runtime
81: */
82: EPSSetFromOptions(eps);
84: /*
85: Attach deflation space: in this case, the matrix has a constant
86: nullspace, [1 1 ... 1]^T is the eigenvector of the zero eigenvalue
87: */
88: MatGetVecs(A,&x,PETSC_NULL);
89: VecSet(&one,x);
90: EPSAttachDeflationSpace(eps,1,&x,PETSC_FALSE);
91: VecDestroy(x);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Solve the eigensystem
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: EPSSolve(eps);
98: EPSGetIterationNumber(eps, &its);
99: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
101: /*
102: Optional: Get some information from the solver and display it
103: */
104: EPSGetType(eps,&type);
105: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
106: EPSGetDimensions(eps,&nev,PETSC_NULL);
107: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
108: EPSGetTolerances(eps,&tol,&maxit);
109: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Display solution and clean up
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: /*
116: Get number of converged approximate eigenpairs
117: */
118: EPSGetConverged(eps,&nconv);
119: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
120:
122: if (nconv>0) {
123: /*
124: Display eigenvalues and relative errors
125: */
126: PetscPrintf(PETSC_COMM_WORLD,
127: " k ||Ax-kx||/||kx||\n"
128: " ----------------- ------------------\n" );
130: for( i=0; i<nconv; i++ ) {
131: /*
132: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
133: ki (imaginary part)
134: */
135: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
136: /*
137: Compute the relative error associated to each eigenpair
138: */
139: EPSComputeRelativeError(eps,i,&error);
141: #ifdef PETSC_USE_COMPLEX
142: re = PetscRealPart(kr);
143: im = PetscImaginaryPart(kr);
144: #else
145: re = kr;
146: im = ki;
147: #endif
148: if (im!=0.0) {
149: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
150: } else {
151: PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,error);
152: }
153: }
154: PetscPrintf(PETSC_COMM_WORLD,"\n" );
155: }
156:
157: /*
158: Free work space
159: */
160: EPSDestroy(eps);
161: MatDestroy(A);
162: SlepcFinalize();
163: return 0;
164: }