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