Actual source code: ex38.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Spectrum slicing on quadratic symmetric eigenproblem.\n\n"
12: "The command line options are:\n"
13: " -n <n> ... dimension of the matrices.\n\n";
15: #include <slepcpep.h>
17: int main(int argc,char **argv)
18: {
19: Mat M,C,K,A[3]; /* problem matrices */
20: PEP pep; /* polynomial eigenproblem solver context */
21: ST st; /* spectral transformation context */
22: KSP ksp;
23: PC pc;
24: PEPType type;
25: PetscBool show=PETSC_FALSE,terse;
26: PetscInt n=100,Istart,Iend,nev,i,*inertias,ns;
27: PetscReal mu=1,tau=10,kappa=5,inta,intb,*shifts;
30: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
32: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
33: PetscOptionsGetBool(NULL,NULL,"-show_inertias",&show,NULL);
34: PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on PEP, n=%D\n\n",n);
36: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: /* K is a tridiagonal */
41: MatCreate(PETSC_COMM_WORLD,&K);
42: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
43: MatSetFromOptions(K);
44: MatSetUp(K);
46: MatGetOwnershipRange(K,&Istart,&Iend);
47: for (i=Istart;i<Iend;i++) {
48: if (i>0) {
49: MatSetValue(K,i,i-1,-kappa,INSERT_VALUES);
50: }
51: MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES);
52: if (i<n-1) {
53: MatSetValue(K,i,i+1,-kappa,INSERT_VALUES);
54: }
55: }
57: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
60: /* C is a tridiagonal */
61: MatCreate(PETSC_COMM_WORLD,&C);
62: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
63: MatSetFromOptions(C);
64: MatSetUp(C);
66: MatGetOwnershipRange(C,&Istart,&Iend);
67: for (i=Istart;i<Iend;i++) {
68: if (i>0) {
69: MatSetValue(C,i,i-1,-tau,INSERT_VALUES);
70: }
71: MatSetValue(C,i,i,tau*3.0,INSERT_VALUES);
72: if (i<n-1) {
73: MatSetValue(C,i,i+1,-tau,INSERT_VALUES);
74: }
75: }
77: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
80: /* M is a diagonal matrix */
81: MatCreate(PETSC_COMM_WORLD,&M);
82: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
83: MatSetFromOptions(M);
84: MatSetUp(M);
85: MatGetOwnershipRange(M,&Istart,&Iend);
86: for (i=Istart;i<Iend;i++) {
87: MatSetValue(M,i,i,mu,INSERT_VALUES);
88: }
89: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create the eigensolver and solve the problem
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /*
97: Create eigensolver context
98: */
99: PEPCreate(PETSC_COMM_WORLD,&pep);
101: /*
102: Set operators and set problem type
103: */
104: A[0] = K; A[1] = C; A[2] = M;
105: PEPSetOperators(pep,3,A);
106: PEPSetProblemType(pep,PEP_HERMITIAN);
108: /*
109: Set interval for spectrum slicing
110: */
111: inta = -11.3;
112: intb = -9.5;
113: PEPSetInterval(pep,inta,intb);
114: PEPSetWhichEigenpairs(pep,PEP_ALL);
116: /*
117: Spectrum slicing requires STOAR
118: */
119: PEPSetType(pep,PEPSTOAR);
121: /*
122: Set shift-and-invert with Cholesky; select MUMPS if available
123: */
124: PEPGetST(pep,&st);
125: STSetType(st,STSINVERT);
127: STGetKSP(st,&ksp);
128: KSPSetType(ksp,KSPPREONLY);
129: KSPGetPC(ksp,&pc);
130: PCSetType(pc,PCCHOLESKY);
132: #if defined(PETSC_HAVE_MUMPS)
133: #if defined(PETSC_USE_COMPLEX)
134: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Spectrum slicing with MUMPS is not available for complex scalars");
135: #endif
136: PEPSTOARSetDetectZeros(pep,PETSC_TRUE); /* enforce zero detection */
137: PCFactorSetMatSolverType(pc,MATSOLVERMUMPS);
138: /*
139: Add several MUMPS options (currently there is no better way of setting this in program):
140: '-mat_mumps_icntl_13 1': turn off ScaLAPACK for matrix inertia
141: '-mat_mumps_icntl_24 1': detect null pivots in factorization (for the case that a shift is equal to an eigenvalue)
142: '-mat_mumps_cntl_3 <tol>': a tolerance used for null pivot detection (must be larger than machine epsilon)
144: Note: depending on the interval, it may be necessary also to increase the workspace:
145: '-mat_mumps_icntl_14 <percentage>': increase workspace with a percentage (50, 100 or more)
146: */
147: PetscOptionsInsertString(NULL,"-mat_mumps_icntl_13 1 -mat_mumps_icntl_24 1 -mat_mumps_cntl_3 1e-12");
148: #endif
150: /*
151: Set solver parameters at runtime
152: */
153: PEPSetFromOptions(pep);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: Solve the eigensystem
157: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158: PEPSetUp(pep);
159: if (show) {
160: PEPSTOARGetInertias(pep,&ns,&shifts,&inertias);
161: PetscPrintf(PETSC_COMM_WORLD,"Subintervals (after setup):\n");
162: for (i=0;i<ns;i++) { PetscPrintf(PETSC_COMM_WORLD,"Shift %g Inertia %D \n",shifts[i],inertias[i]); }
163: PetscPrintf(PETSC_COMM_WORLD,"\n");
164: PetscFree(shifts);
165: PetscFree(inertias);
166: }
167: PEPSolve(pep);
168: if (show) {
169: PEPSTOARGetInertias(pep,&ns,&shifts,&inertias);
170: PetscPrintf(PETSC_COMM_WORLD,"All shifts (after solve):\n");
171: for (i=0;i<ns;i++) { PetscPrintf(PETSC_COMM_WORLD,"Shift %g Inertia %D \n",shifts[i],inertias[i]); }
172: PetscPrintf(PETSC_COMM_WORLD,"\n");
173: PetscFree(shifts);
174: PetscFree(inertias);
175: }
177: /*
178: Show eigenvalues in interval and print solution
179: */
180: PEPGetType(pep,&type);
181: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
182: PEPGetDimensions(pep,&nev,NULL,NULL);
183: PEPGetInterval(pep,&inta,&intb);
184: PetscPrintf(PETSC_COMM_WORLD," %D eigenvalues found in [%g, %g]\n",nev,(double)inta,(double)intb);
186: /*
187: Show detailed info unless -terse option is given by user
188: */
189: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
190: if (terse) {
191: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
192: } else {
193: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
194: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
195: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
196: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
197: }
199: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200: Clean up
201: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: PEPDestroy(&pep);
203: MatDestroy(&M);
204: MatDestroy(&C);
205: MatDestroy(&K);
206: SlepcFinalize();
207: return ierr;
208: }