Actual source code: ex12.c
slepc-3.7.1 2016-05-27
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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[] = "Compute all eigenvalues in an interval of a symmetric-definite problem.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
25: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
27: #include <slepceps.h>
31: int main(int argc,char **argv)
32: {
33: Mat A,B; /* matrices */
34: EPS eps; /* eigenproblem solver context */
35: ST st; /* spectral transformation context */
36: KSP ksp;
37: PC pc;
38: PetscInt N,n=35,m,Istart,Iend,II,nev,i,j,k,*inertias;
39: PetscBool flag;
40: PetscReal int0,int1,*shifts;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
47: if (!flag) m=n;
48: N = n*m;
49: PetscPrintf(PETSC_COMM_WORLD,"\nSymmetric-definite problem with two intervals, N=%D (%Dx%D grid)\n\n",N,n,m);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the matrices that define the eigensystem, Ax=kBx
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: MatCreate(PETSC_COMM_WORLD,&A);
56: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
57: MatSetFromOptions(A);
58: MatSetUp(A);
60: MatCreate(PETSC_COMM_WORLD,&B);
61: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
62: MatSetFromOptions(B);
63: MatSetUp(B);
65: MatGetOwnershipRange(A,&Istart,&Iend);
66: for (II=Istart;II<Iend;II++) {
67: i = II/n; j = II-i*n;
68: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
69: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
70: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
71: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
72: MatSetValue(A,II,II,4.0,INSERT_VALUES);
73: MatSetValue(B,II,II,2.0,INSERT_VALUES);
74: }
75: if (Istart==0) {
76: MatSetValue(B,0,0,6.0,INSERT_VALUES);
77: MatSetValue(B,0,1,-1.0,INSERT_VALUES);
78: MatSetValue(B,1,0,-1.0,INSERT_VALUES);
79: MatSetValue(B,1,1,1.0,INSERT_VALUES);
80: }
82: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
84: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Create the eigensolver and set various options
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: EPSCreate(PETSC_COMM_WORLD,&eps);
92: EPSSetOperators(eps,A,B);
93: EPSSetProblemType(eps,EPS_GHEP);
95: /*
96: Set first interval and other settings for spectrum slicing
97: */
98: EPSSetWhichEigenpairs(eps,EPS_ALL);
99: EPSSetInterval(eps,1.1,1.3);
100: EPSGetST(eps,&st);
101: STSetType(st,STSINVERT);
102: STGetKSP(st,&ksp);
103: KSPGetPC(ksp,&pc);
104: KSPSetType(ksp,KSPPREONLY);
105: PCSetType(pc,PCCHOLESKY);
106: EPSSetFromOptions(eps);
108: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: Solve for first interval and display info
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: EPSSolve(eps);
113: EPSGetDimensions(eps,&nev,NULL,NULL);
114: EPSGetInterval(eps,&int0,&int1);
115: PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
116: EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
117: PetscPrintf(PETSC_COMM_WORLD," Used %D shifts (inertia):\n",k);
118: for (i=0;i<k;i++) {
119: PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
120: }
121: PetscFree(shifts);
122: PetscFree(inertias);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Solve for second interval and display info
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: EPSSetInterval(eps,1.5,1.6);
128: EPSSolve(eps);
129: EPSGetDimensions(eps,&nev,NULL,NULL);
130: EPSGetInterval(eps,&int0,&int1);
131: PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
132: EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
133: PetscPrintf(PETSC_COMM_WORLD," Used %D shifts (inertia):\n",k);
134: for (i=0;i<k;i++) {
135: PetscPrintf(PETSC_COMM_WORLD," .. %g (%D)\n",(double)shifts[i],inertias[i]);
136: }
137: PetscFree(shifts);
138: PetscFree(inertias);
140: EPSDestroy(&eps);
141: MatDestroy(&A);
142: MatDestroy(&B);
143: SlepcFinalize();
144: return ierr;
145: }