Actual source code: ex25.c
slepc-3.7.2 2016-07-19
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[] = "Spectrum slicing on generalized symmetric eigenproblem.\n\n"
23: "The problem is similar to ex13.c.\n\n"
24: "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";
28: #include <slepceps.h>
32: int main(int argc,char **argv)
33: {
34: Mat A,B; /* matrices */
35: EPS eps; /* eigenproblem solver context */
36: ST st; /* spectral transformation context */
37: KSP ksp;
38: PC pc;
39: EPSType type;
40: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j,*inertias,ns;
41: PetscReal inta,intb,*shifts;
42: PetscBool flag,show=PETSC_FALSE,terse;
45: SlepcInitialize(&argc,&argv,(char*)0,help);
47: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
48: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
49: PetscOptionsGetBool(NULL,NULL,"-show_inertias",&show,NULL);
50: if (!flag) m=n;
51: N = n*m;
52: PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on GHEP, N=%D (%Dx%D grid)\n\n",N,n,m);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Compute the matrices that define the eigensystem, Ax=kBx
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: MatCreate(PETSC_COMM_WORLD,&A);
59: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
60: MatSetFromOptions(A);
61: MatSetUp(A);
63: MatCreate(PETSC_COMM_WORLD,&B);
64: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
65: MatSetFromOptions(B);
66: MatSetUp(B);
68: MatGetOwnershipRange(A,&Istart,&Iend);
69: for (II=Istart;II<Iend;II++) {
70: i = II/n; j = II-i*n;
71: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
72: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
73: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
74: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
75: MatSetValue(A,II,II,4.0,INSERT_VALUES);
76: MatSetValue(B,II,II,4.0,INSERT_VALUES);
77: }
79: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
80: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
81: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create the eigensolver and set various options
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: /*
89: Create eigensolver context
90: */
91: EPSCreate(PETSC_COMM_WORLD,&eps);
93: /*
94: Set operators and set problem type
95: */
96: EPSSetOperators(eps,A,B);
97: EPSSetProblemType(eps,EPS_GHEP);
99: /*
100: Set interval for spectrum slicing
101: */
102: inta = 0.1;
103: intb = 0.2;
104: EPSSetInterval(eps,inta,intb);
105: EPSSetWhichEigenpairs(eps,EPS_ALL);
107: /*
108: Spectrum slicing requires Krylov-Schur
109: */
110: EPSSetType(eps,EPSKRYLOVSCHUR);
112: /*
113: Set shift-and-invert with Cholesky; select MUMPS if available
114: */
116: EPSGetST(eps,&st);
117: STSetType(st,STSINVERT);
118:
119: STGetKSP(st,&ksp);
120: KSPSetType(ksp,KSPPREONLY);
121: KSPGetPC(ksp,&pc);
122: PCSetType(pc,PCCHOLESKY);
123:
124: #if defined(PETSC_HAVE_MUMPS)
125: #if defined(PETSC_USE_COMPLEX)
126: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Spectrum slicing with MUMPS is not available for complex scalars");
127: #endif
128: EPSKrylovSchurSetDetectZeros(eps,PETSC_TRUE); /* enforce zero detection */
129: PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);
130: /*
131: Add several MUMPS options (currently there is no better way of setting this in program):
132: '-mat_mumps_icntl_13 1': turn off ScaLAPACK for matrix inertia
133: '-mat_mumps_icntl_24 1': detect null pivots in factorization (for the case that a shift is equal to an eigenvalue)
134: '-mat_mumps_cntl_3 <tol>': a tolerance used for null pivot detection (must be larger than machine epsilon)
136: Note: depending on the interval, it may be necessary also to increase the workspace:
137: '-mat_mumps_icntl_14 <percentage>': increase workspace with a percentage (50, 100 or more)
138: */
139: PetscOptionsInsertString(NULL,"-mat_mumps_icntl_13 1 -mat_mumps_icntl_24 1 -mat_mumps_cntl_3 1e-12");
140: #endif
142: /*
143: Set solver parameters at runtime
144: */
145: EPSSetFromOptions(eps);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Solve the eigensystem
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: EPSSetUp(eps);
151: if (show) {
152: EPSKrylovSchurGetInertias(eps,&ns,&shifts,&inertias);
153: PetscPrintf(PETSC_COMM_WORLD,"Subintervals (after setup):\n");
154: for (i=0;i<ns;i++) { PetscPrintf(PETSC_COMM_WORLD,"Shift %g Inertia %D \n",shifts[i],inertias[i]); }
155: PetscPrintf(PETSC_COMM_WORLD,"\n");
156: PetscFree(shifts);
157: PetscFree(inertias);
158: }
159: EPSSolve(eps);
160: if (show) {
161: EPSKrylovSchurGetInertias(eps,&ns,&shifts,&inertias);
162: PetscPrintf(PETSC_COMM_WORLD,"All shifts (after solve):\n");
163: for (i=0;i<ns;i++) { PetscPrintf(PETSC_COMM_WORLD,"Shift %g Inertia %D \n",shifts[i],inertias[i]); }
164: PetscPrintf(PETSC_COMM_WORLD,"\n");
165: PetscFree(shifts);
166: PetscFree(inertias);
167: }
169: /*
170: Show eigenvalues in interval and print solution
171: */
172: EPSGetType(eps,&type);
173: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
174: EPSGetDimensions(eps,&nev,NULL,NULL);
175: EPSGetInterval(eps,&inta,&intb);
176: PetscPrintf(PETSC_COMM_WORLD," %D eigenvalues found in [%g, %g]\n",nev,(double)inta,(double)intb);
178: /*
179: Show detailed info unless -terse option is given by user
180: */
181: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
182: if (terse) {
183: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
184: } else {
185: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
186: EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
187: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
188: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
189: }
191: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192: Clean up
193: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
194: EPSDestroy(&eps);
195: MatDestroy(&A);
196: MatDestroy(&B);
197: SlepcFinalize();
198: return ierr;
199: }