Actual source code: ex12.c
slepc-main 2024-11-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
29: PetscFunctionBeginUser;
30: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
32: PetscCall(PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
34: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
35: if (!flag) m=n;
36: N = n*m;
37: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSymmetric-definite problem with two intervals, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Compute the matrices that define the eigensystem, Ax=kBx
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
44: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
45: PetscCall(MatSetFromOptions(A));
47: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
48: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
49: PetscCall(MatSetFromOptions(B));
51: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
52: for (II=Istart;II<Iend;II++) {
53: i = II/n; j = II-i*n;
54: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
55: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
56: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
57: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
58: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
59: PetscCall(MatSetValue(B,II,II,2.0,INSERT_VALUES));
60: }
61: if (Istart==0) {
62: PetscCall(MatSetValue(B,0,0,6.0,INSERT_VALUES));
63: PetscCall(MatSetValue(B,0,1,-1.0,INSERT_VALUES));
64: PetscCall(MatSetValue(B,1,0,-1.0,INSERT_VALUES));
65: PetscCall(MatSetValue(B,1,1,1.0,INSERT_VALUES));
66: }
68: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Create the eigensolver and set various options
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
78: PetscCall(EPSSetOperators(eps,A,B));
79: PetscCall(EPSSetProblemType(eps,EPS_GHEP));
81: /*
82: Set first interval and other settings for spectrum slicing
83: */
84: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
85: PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
86: PetscCall(EPSSetInterval(eps,1.1,1.3));
87: PetscCall(EPSGetST(eps,&st));
88: PetscCall(STSetType(st,STSINVERT));
89: PetscCall(EPSKrylovSchurGetKSP(eps,&ksp));
90: PetscCall(KSPGetPC(ksp,&pc));
91: PetscCall(KSPSetType(ksp,KSPPREONLY));
92: PetscCall(PCSetType(pc,PCCHOLESKY));
93: PetscCall(EPSSetFromOptions(eps));
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Solve for first interval and display info
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: PetscCall(EPSSolve(eps));
100: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
101: PetscCall(EPSGetInterval(eps,&int0,&int1));
102: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));
103: if (showinertia) {
104: PetscCall(EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias));
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k));
106: for (i=0;i<k;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
107: PetscCall(PetscFree(shifts));
108: PetscCall(PetscFree(inertias));
109: }
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Solve for second interval and display info
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: PetscCall(EPSSetInterval(eps,1.499,1.6));
115: PetscCall(EPSSolve(eps));
116: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
117: PetscCall(EPSGetInterval(eps,&int0,&int1));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));
119: if (showinertia) {
120: PetscCall(EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias));
121: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k));
122: for (i=0;i<k;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
123: PetscCall(PetscFree(shifts));
124: PetscCall(PetscFree(inertias));
125: }
127: PetscCall(EPSDestroy(&eps));
128: PetscCall(MatDestroy(&A));
129: PetscCall(MatDestroy(&B));
130: PetscCall(SlepcFinalize());
131: return 0;
132: }
134: /*TEST
136: test:
137: suffix: 1
138: args: -showinertia 0 -eps_error_relative
139: requires: !single
141: TEST*/