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