Actual source code: ex11.c
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
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: EPS eps; /* eigenproblem solver context */
35: Mat A; /* operator matrix */
36: Vec x;
37: EPSType type;
38: PetscInt N,n=10,m,i,j,II,Istart,Iend,nev;
39: PetscScalar w;
40: PetscBool flag,terse;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,"-n",&n,NULL);
46: PetscOptionsGetInt(NULL,"-m",&m,&flag);
47: if (!flag) m=n;
48: N = n*m;
49: PetscPrintf(PETSC_COMM_WORLD,"\nFiedler vector of a 2-D regular mesh, N=%D (%Dx%D grid)\n\n",N,n,m);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the operator matrix that defines the eigensystem, Ax=kx
53: In this example, A = L(G), where L is the Laplacian of graph G, i.e.
54: Lii = degree of node i, Lij = -1 if edge (i,j) exists in G
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
59: MatSetFromOptions(A);
60: MatSetUp(A);
62: MatGetOwnershipRange(A,&Istart,&Iend);
63: for (II=Istart;II<Iend;II++) {
64: i = II/n; j = II-i*n;
65: w = 0.0;
66: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); w=w+1.0; }
67: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); w=w+1.0; }
68: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); w=w+1.0; }
69: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); w=w+1.0; }
70: MatSetValue(A,II,II,w,INSERT_VALUES);
71: }
73: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Create the eigensolver and set various options
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: /*
81: Create eigensolver context
82: */
83: EPSCreate(PETSC_COMM_WORLD,&eps);
85: /*
86: Set operators. In this case, it is a standard eigenvalue problem
87: */
88: EPSSetOperators(eps,A,NULL);
89: EPSSetProblemType(eps,EPS_HEP);
91: /*
92: Select portion of spectrum
93: */
94: EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
96: /*
97: Set solver parameters at runtime
98: */
99: EPSSetFromOptions(eps);
101: /*
102: Attach deflation space: in this case, the matrix has a constant
103: nullspace, [1 1 ... 1]^T is the eigenvector of the zero eigenvalue
104: */
105: MatCreateVecs(A,&x,NULL);
106: VecSet(x,1.0);
107: EPSSetDeflationSpace(eps,1,&x);
108: VecDestroy(&x);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Solve the eigensystem
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: EPSSolve(eps);
116: /*
117: Optional: Get some information from the solver and display it
118: */
119: EPSGetType(eps,&type);
120: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
121: EPSGetDimensions(eps,&nev,NULL,NULL);
122: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Display solution and clean up
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: /* show detailed info unless -terse option is given by user */
129: PetscOptionsHasName(NULL,"-terse",&terse);
130: if (terse) {
131: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
132: } else {
133: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
134: EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
135: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
136: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
137: }
138: EPSDestroy(&eps);
139: MatDestroy(&A);
140: SlepcFinalize();
141: return 0;
142: }