Actual source code: ex12.c
slepc-3.20.1 2023-11-27
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,(char*)0,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));
46: PetscCall(MatSetUp(A));
48: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
49: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
50: PetscCall(MatSetFromOptions(B));
51: PetscCall(MatSetUp(B));
53: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
54: for (II=Istart;II<Iend;II++) {
55: i = II/n; j = II-i*n;
56: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
57: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
58: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
59: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
60: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
61: PetscCall(MatSetValue(B,II,II,2.0,INSERT_VALUES));
62: }
63: if (Istart==0) {
64: PetscCall(MatSetValue(B,0,0,6.0,INSERT_VALUES));
65: PetscCall(MatSetValue(B,0,1,-1.0,INSERT_VALUES));
66: PetscCall(MatSetValue(B,1,0,-1.0,INSERT_VALUES));
67: PetscCall(MatSetValue(B,1,1,1.0,INSERT_VALUES));
68: }
70: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
72: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
73: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Create the eigensolver and set various options
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
80: PetscCall(EPSSetOperators(eps,A,B));
81: PetscCall(EPSSetProblemType(eps,EPS_GHEP));
83: /*
84: Set first interval and other settings for spectrum slicing
85: */
86: PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
87: PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
88: PetscCall(EPSSetInterval(eps,1.1,1.3));
89: PetscCall(EPSGetST(eps,&st));
90: PetscCall(STSetType(st,STSINVERT));
91: PetscCall(EPSKrylovSchurGetKSP(eps,&ksp));
92: PetscCall(KSPGetPC(ksp,&pc));
93: PetscCall(KSPSetType(ksp,KSPPREONLY));
94: PetscCall(PCSetType(pc,PCCHOLESKY));
95: PetscCall(EPSSetFromOptions(eps));
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Solve for first interval and display info
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: PetscCall(EPSSolve(eps));
102: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
103: PetscCall(EPSGetInterval(eps,&int0,&int1));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));
105: if (showinertia) {
106: PetscCall(EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias));
107: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k));
108: for (i=0;i<k;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
109: PetscCall(PetscFree(shifts));
110: PetscCall(PetscFree(inertias));
111: }
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Solve for second interval and display info
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: PetscCall(EPSSetInterval(eps,1.499,1.6));
117: PetscCall(EPSSolve(eps));
118: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
119: PetscCall(EPSGetInterval(eps,&int0,&int1));
120: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1));
121: if (showinertia) {
122: PetscCall(EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias));
123: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k));
124: for (i=0;i<k;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]));
125: PetscCall(PetscFree(shifts));
126: PetscCall(PetscFree(inertias));
127: }
129: PetscCall(EPSDestroy(&eps));
130: PetscCall(MatDestroy(&A));
131: PetscCall(MatDestroy(&B));
132: PetscCall(SlepcFinalize());
133: return 0;
134: }
136: /*TEST
138: test:
139: suffix: 1
140: args: -showinertia 0 -eps_error_relative
141: requires: !single
143: TEST*/