Actual source code: ex13.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat 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[] = "Generalized Symmetric eigenproblem.\n\n"
23: "The problem is Ax = lambda Bx, with:\n"
24: " A = Laplacian operator in 2-D\n"
25: " B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
26: "The command line options are:\n"
27: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
28: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
29: " -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";
31: #include <slepceps.h>
35: int main(int argc,char **argv)
36: {
37: Mat A,B; /* matrices */
38: EPS eps; /* eigenproblem solver context */
39: ST st; /* spectral transformation context */
40: const EPSType type;
41: PetscReal tol;
42: PetscInt N,n=10,m,Istart,Iend,II,nev,maxit,i,j,its,nulldim=0;
43: PetscBool flag;
46: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
49: PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);
50: if(!flag) m=n;
51: N = n*m;
52: PetscOptionsGetInt(PETSC_NULL,"-nulldim",&nulldim,PETSC_NULL);
53: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%D (%Dx%D grid), null(B)=%D\n\n",N,n,m,nulldim);
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Compute the matrices that define the eigensystem, Ax=kBx
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
61: MatSetFromOptions(A);
62:
63: MatCreate(PETSC_COMM_WORLD,&B);
64: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
65: MatSetFromOptions(B);
67: MatGetOwnershipRange(A,&Istart,&Iend);
68: for (II=Istart;II<Iend;II++) {
69: i = II/n; j = II-i*n;
70: if(i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
71: if(i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
72: if(j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
73: if(j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
74: MatSetValue(A,II,II,4.0,INSERT_VALUES);
75: if (II>=nulldim) { MatSetValue(B,II,II,4.0,INSERT_VALUES); }
76: }
78: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
80: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Create the eigensolver and set various options
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: /*
88: Create eigensolver context
89: */
90: EPSCreate(PETSC_COMM_WORLD,&eps);
92: /*
93: Set operators. In this case, it is a generalized eigenvalue problem
94: */
95: EPSSetOperators(eps,A,B);
96: EPSSetProblemType(eps,EPS_GHEP);
98: /*
99: Select portion of spectrum
100: */
101: EPSSetTarget(eps,0.0);
102: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
104: /*
105: Set solver parameters at runtime
106: */
107: EPSSetFromOptions(eps);
109: /*
110: Use shift-and-invert to avoid solving linear systems with a singular B
111: in case nulldim>0
112: */
113: PetscTypeCompareAny((PetscObject)eps,&flag,EPSGD,EPSJD,"");
114: if (!flag) {
115: EPSGetST(eps,&st);
116: STSetType(st,STSINVERT);
117: }
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Solve the eigensystem
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: EPSSolve(eps);
124: EPSGetIterationNumber(eps,&its);
125: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
127: /*
128: Optional: Get some information from the solver and display it
129: */
130: EPSGetType(eps,&type);
131: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
132: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
133: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
134: EPSGetTolerances(eps,&tol,&maxit);
135: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Display solution and clean up
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: EPSPrintSolution(eps,PETSC_NULL);
142: EPSDestroy(&eps);
143: MatDestroy(&A);
144: MatDestroy(&B);
145: SlepcFinalize();
146: return 0;
147: }